1
0
Fork 0
mirror of https://github.com/cosmo-sims/monofonIC.git synced 2024-09-19 17:03:45 +02:00

fixed inconsistency between NGENIC random numbers generated with MPI and without

This commit is contained in:
Oliver Hahn 2019-11-15 22:32:09 +01:00
parent 07b430f25c
commit 89ec1775f3
2 changed files with 28 additions and 9 deletions

View file

@ -73,6 +73,8 @@ public:
const grid_fft_t *get_grid(size_t ilevel) const { return this; }
bool is_distributed( void ) const { return bdistributed; }
void Setup();
//! return the (local) size of dimension i

View file

@ -82,7 +82,11 @@ public:
for (size_t j = 0; j < nres_; ++j)
{
ptrdiff_t jj = (j>0)? nres_ - j : 0;
gsl_rng_set( pRandomGenerator_, SeedTable_[i * nres_ + j]);
if( g.is_distributed() )
gsl_rng_set( pRandomGenerator_, SeedTable_[j * nres_ + i]);
else
gsl_rng_set( pRandomGenerator_, SeedTable_[i * nres_ + j]);
for (size_t k = 0; k < g.size(2); ++k)
{
double phase = gsl_rng_uniform(pRandomGenerator_) * 2 * M_PI;
@ -101,15 +105,28 @@ public:
if (k > 0) {
if (i_in_range) g.kelem(ip,j,k) = zrand;
} else{ /* k=0 plane needs special treatment */
if (i == 0) {
if (j < nres_ / 2 && i_in_range)
{
g.kelem(ip,j,k) = zrand;
g.kelem(ip,jj,k) = std::conj(zrand);
if( g.is_distributed() ){
if (j == 0) {
if (i < nres_ / 2 && i_in_range)
{
if(i_in_range) g.kelem(ip,jj,k) = zrand;
if(ii_in_range) g.kelem(iip,j,k) = std::conj(zrand);
}
} else if (j < nres_ / 2) {
if(i_in_range) g.kelem(ip,j,k) = zrand;
if(ii_in_range) g.kelem(iip,jj,k) = std::conj(zrand);
}
}else{
if (i == 0) {
if (j < nres_ / 2 && i_in_range)
{
g.kelem(ip,j,k) = zrand;
g.kelem(ip,jj,k) = std::conj(zrand);
}
} else if (i < nres_ / 2) {
if(i_in_range) g.kelem(ip,j,k) = zrand;
if (ii_in_range) g.kelem(iip,jj,k) = std::conj(zrand);
}
} else if (i < nres_ / 2) {
if(i_in_range) g.kelem(ip,j,k) = zrand;
if (ii_in_range) g.kelem(iip,jj,k) = std::conj(zrand);
}
}
}