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

WIP fixes/cleanup for MUSIC RNG

This commit is contained in:
Oliver Hahn 2023-02-12 22:38:50 -08:00
parent 244f11a821
commit ced852069e
2 changed files with 292 additions and 165 deletions

View file

@ -7,6 +7,157 @@
#include "random.hh"
#include "random_music_wnoise_generator.hh"
template <typename T>
void rapid_proto_ngenic_rng(size_t res, long baseseed, music_wnoise_generator<T> &R)
{
LOGUSER("Invoking the N-GenIC random number generator");
unsigned *seedtable = new unsigned[res * res];
gsl_rng *random_generator = gsl_rng_alloc(gsl_rng_ranlxd1);
gsl_rng_set(random_generator, baseseed);
for (size_t i = 0; i < res / 2; i++)
{
size_t j;
for (j = 0; j < i; j++)
seedtable[i * res + j] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i + 1; j++)
seedtable[j * res + i] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i; j++)
seedtable[(res - 1 - i) * res + j] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i + 1; j++)
seedtable[(res - 1 - j) * res + i] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i; j++)
seedtable[i * res + (res - 1 - j)] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i + 1; j++)
seedtable[j * res + (res - 1 - i)] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i; j++)
seedtable[(res - 1 - i) * res + (res - 1 - j)] = 0x7fffffff * gsl_rng_uniform(random_generator);
for (j = 0; j < i + 1; j++)
seedtable[(res - 1 - j) * res + (res - 1 - i)] = 0x7fffffff * gsl_rng_uniform(random_generator);
}
fftw_real *rnoise = new fftw_real[res * res * (res + 2)];
fftw_complex *knoise = reinterpret_cast<fftw_complex *>(rnoise);
double fnorm = 1. / sqrt(res * res * res);
// #warning need to check for race conditions below
//#pragma omp parallel for
for (size_t i = 0; i < res; i++)
{
int ii = (int)res - (int)i;
if (ii == (int)res)
ii = 0;
for (size_t j = 0; j < res; j++)
{
gsl_rng_set(random_generator, seedtable[i * res + j]);
for (size_t k = 0; k < res / 2; k++)
{
double phase = gsl_rng_uniform(random_generator) * 2 * M_PI;
double ampl;
do
ampl = gsl_rng_uniform(random_generator);
while (ampl == 0);
if (i == res / 2 || j == res / 2 || k == res / 2)
continue;
if (i == 0 && j == 0 && k == 0)
continue;
T rp = -sqrt(-log(ampl)) * cos(phase) * fnorm;
T ip = -sqrt(-log(ampl)) * sin(phase) * fnorm;
if (k > 0)
{
RE(knoise[(i * res + j) * (res / 2 + 1) + k]) = rp;
IM(knoise[(i * res + j) * (res / 2 + 1) + k]) = ip;
}
else /* k=0 plane needs special treatment */
{
if (i == 0)
{
if (j >= res / 2)
continue;
else
{
int jj = (int)res - (int)j; /* note: j!=0 surely holds at this point */
RE(knoise[(i * res + j) * (res / 2 + 1) + k]) = rp;
IM(knoise[(i * res + j) * (res / 2 + 1) + k]) = ip;
RE(knoise[(i * res + jj) * (res / 2 + 1) + k]) = rp;
IM(knoise[(i * res + jj) * (res / 2 + 1) + k]) = -ip;
}
}
else
{
if (i >= res / 2)
continue;
else
{
int ii = (int)res - (int)i;
if (ii == (int)res)
ii = 0;
int jj = (int)res - (int)j;
if (jj == (int)res)
jj = 0;
RE(knoise[(i * res + j) * (res / 2 + 1) + k]) = rp;
IM(knoise[(i * res + j) * (res / 2 + 1) + k]) = ip;
if (ii >= 0 && ii < (int)res)
{
RE(knoise[(ii * res + jj) * (res / 2 + 1) + k]) = rp;
IM(knoise[(ii * res + jj) * (res / 2 + 1) + k]) = -ip;
}
}
}
}
}
}
}
delete[] seedtable;
//... perform FT to real space
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_plan plan = fftwf_plan_dft_c2r_3d(res, res, res, knoise, rnoise, FFTW_ESTIMATE);
fftwf_execute(plan);
fftwf_destroy_plan(plan);
#else
fftw_plan plan = fftw_plan_dft_c2r_3d(res, res, res, knoise, rnoise, FFTW_ESTIMATE);
fftw_execute(plan);
fftw_destroy_plan(plan);
#endif
#else
rfftwnd_plan plan = rfftw3d_create_plan(res, res, res, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE | FFTW_IN_PLACE);
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_complex_to_real(omp_get_max_threads(), plan, knoise, NULL);
#else
rfftwnd_one_complex_to_real(plan, knoise, NULL);
#endif
rfftwnd_destroy_plan(plan);
#endif
// copy to array that holds the random numbers
#pragma omp parallel for
for (int i = 0; i < (int)res; ++i)
for (size_t j = 0; j < res; ++j)
for (size_t k = 0; k < res; ++k)
R(i, j, k) = rnoise[((size_t)i * res + j) * res + k];
delete[] rnoise;
}
template <typename T>
music_wnoise_generator<T>::music_wnoise_generator(unsigned res, unsigned cubesize, long baseseed, int *x0, int *lx)
: res_(res), cubesize_(cubesize), ncubes_(1), baseseed_(baseseed)
@ -42,22 +193,8 @@ music_wnoise_generator<T>::music_wnoise_generator(unsigned res, unsigned cubesiz
rnums_.push_back(new Meshvar<T>(res, 0, 0, 0));
cubemap_[0] = 0; // create dummy map index
register_cube(0, 0, 0);
// rapid_proto_ngenic_rng( res_, baseseed_, *this );
}
/*
if( musicnoise )
mean = fill_all();
else
{
rnums_.push_back( new Meshvar<T>( res, 0, 0, 0 ) );
cubemap_[0] = 0; // create dummy map index
register_cube(0,0,0);
rapid_proto_ngenic_rng( res_, baseseed_, *this );
}
*/
}
if (zeromean)
{
@ -847,8 +984,7 @@ double music_wnoise_generator<T>::fill_all(void)
register_cube(ii, jj, kk);
}
#pragma omp parallel for reduction(+ \
: sum)
#pragma omp parallel for reduction(+ : sum)
for (int i = 0; i < (int)ncubes_; ++i)
for (int j = 0; j < (int)ncubes_; ++j)
for (int k = 0; k < (int)ncubes_; ++k)
@ -862,9 +998,8 @@ double music_wnoise_generator<T>::fill_all(void)
sum += fill_cube(ii, jj, kk);
}
//... subtract mean
#pragma omp parallel for reduction(+ \
: sum)
//... subtract mean
#pragma omp parallel for reduction(+ : sum)
for (int i = 0; i < (int)ncubes_; ++i)
for (int j = 0; j < (int)ncubes_; ++j)
for (int k = 0; k < (int)ncubes_; ++k)

View file

@ -1,205 +1,197 @@
#ifndef __RANDOM_MUSIC_WNOISE_GENERATOR_HH
#define __RANDOM_MUSIC_WNOISE_GENERATOR_HH
#pragma once
#include <vector>
#include <map>
#include "general.hh"
#include "mesh.hh"
#define DEF_RAN_CUBE_SIZE 32
#define DEF_RAN_CUBE_SIZE 32
/*!
* @brief encapsulates all things random number generator related
*/
template< typename T >
template <typename T>
class music_wnoise_generator
{
public:
unsigned
res_, //!< resolution of the full mesh
cubesize_, //!< size of one independent random number cube
ncubes_; //!< number of random number cubes to cover the full mesh
long baseseed_; //!< base seed from which cube seeds are computed
unsigned
res_, //!< resolution of the full mesh
cubesize_, //!< size of one independent random number cube
ncubes_; //!< number of random number cubes to cover the full mesh
long baseseed_; //!< base seed from which cube seeds are computed
protected:
//! vector of 3D meshes (the random number cubes) with random numbers
std::vector< Meshvar<T>* > rnums_;
std::vector<Meshvar<T> *> rnums_;
//! map of 3D indices to cube index
std::map<size_t,size_t> cubemap_;
typedef std::map<size_t,size_t>::iterator cubemap_iterator;
std::map<size_t, size_t> cubemap_;
typedef std::map<size_t, size_t>::iterator cubemap_iterator;
protected:
//! register a cube with the hash map
void register_cube( int i, int j, int k);
void register_cube(int i, int j, int k);
//! fills a subcube with random numbers
double fill_cube( int i, int j, int k);
double fill_cube(int i, int j, int k);
//! subtract a constant from an entire cube
void subtract_from_cube( int i, int j, int k, double val );
void subtract_from_cube(int i, int j, int k, double val);
//! copy random numbers from a cube to a full grid array
template< class C >
void copy_cube( int i, int j, int k, C& dat )
template <class C>
void copy_cube(int i, int j, int k, C &dat)
{
int offi, offj, offk;
offi = i*cubesize_;
offj = j*cubesize_;
offk = k*cubesize_;
i = (i+ncubes_)%ncubes_;
j = (j+ncubes_)%ncubes_;
k = (k+ncubes_)%ncubes_;
size_t icube = (i*ncubes_+j)*ncubes_+k;
cubemap_iterator it = cubemap_.find( icube );
if( it == cubemap_.end() )
{
LOGERR("attempting to copy data from non-existing RND cube %d,%d,%d",i,j,k);
throw std::runtime_error("attempting to copy data from non-existing RND cube");
}
offi = i * cubesize_;
offj = j * cubesize_;
offk = k * cubesize_;
i = (i + ncubes_) % ncubes_;
j = (j + ncubes_) % ncubes_;
k = (k + ncubes_) % ncubes_;
size_t icube = (i * ncubes_ + j) * ncubes_ + k;
cubemap_iterator it = cubemap_.find(icube);
if (it == cubemap_.end())
{
LOGERR("attempting to copy data from non-existing RND cube %d,%d,%d", i, j, k);
throw std::runtime_error("attempting to copy data from non-existing RND cube");
}
size_t cubeidx = it->second;
for( int ii=0; ii<(int)cubesize_; ++ii )
for( int jj=0; jj<(int)cubesize_; ++jj )
for( int kk=0; kk<(int)cubesize_; ++kk )
dat(offi+ii,offj+jj,offk+kk) = (*rnums_[cubeidx])(ii,jj,kk);
for (int ii = 0; ii < (int)cubesize_; ++ii)
for (int jj = 0; jj < (int)cubesize_; ++jj)
for (int kk = 0; kk < (int)cubesize_; ++kk)
dat(offi + ii, offj + jj, offk + kk) = (*rnums_[cubeidx])(ii, jj, kk);
}
//! free the memory associated with a subcube
void free_cube( int i, int j, int k );
void free_cube(int i, int j, int k);
//! initialize member variables and allocate memory
void initialize( void );
void initialize(void);
//! fill a cubic subvolume of the full grid with random numbers
double fill_subvolume( int *i0, int *n );
double fill_subvolume(int *i0, int *n);
//! fill an entire grid with random numbers
double fill_all( void );
double fill_all(void);
//! fill an external array instead of the internal field
template< class C >
double fill_all( C& dat )
template <class C>
double fill_all(C &dat)
{
double sum = 0.0;
for( int i=0; i<(int)ncubes_; ++i )
for( int j=0; j<(int)ncubes_; ++j )
for( int k=0; k<(int)ncubes_; ++k )
{
int ii(i),jj(j),kk(k);
register_cube(ii,jj,kk);
}
#pragma omp parallel for reduction(+:sum)
for( int i=0; i<(int)ncubes_; ++i )
for( int j=0; j<(int)ncubes_; ++j )
for( int k=0; k<(int)ncubes_; ++k )
{
int ii(i),jj(j),kk(k);
ii = (ii+ncubes_)%ncubes_;
jj = (jj+ncubes_)%ncubes_;
kk = (kk+ncubes_)%ncubes_;
sum+=fill_cube(ii, jj, kk);
copy_cube(ii,jj,kk,dat);
free_cube(ii, jj, kk);
}
return sum/(ncubes_*ncubes_*ncubes_);
for (int i = 0; i < (int)ncubes_; ++i)
for (int j = 0; j < (int)ncubes_; ++j)
for (int k = 0; k < (int)ncubes_; ++k)
{
int ii(i), jj(j), kk(k);
register_cube(ii, jj, kk);
}
#pragma omp parallel for reduction(+ \
: sum)
for (int i = 0; i < (int)ncubes_; ++i)
for (int j = 0; j < (int)ncubes_; ++j)
for (int k = 0; k < (int)ncubes_; ++k)
{
int ii(i), jj(j), kk(k);
ii = (ii + ncubes_) % ncubes_;
jj = (jj + ncubes_) % ncubes_;
kk = (kk + ncubes_) % ncubes_;
sum += fill_cube(ii, jj, kk);
copy_cube(ii, jj, kk, dat);
free_cube(ii, jj, kk);
}
return sum / (ncubes_ * ncubes_ * ncubes_);
}
//! write the number of allocated random number cubes to stdout
void print_allocated( void );
void print_allocated(void);
public:
//! constructor
music_wnoise_generator( unsigned res, unsigned cubesize, long baseseed, int *x0, int *lx );
music_wnoise_generator(unsigned res, unsigned cubesize, long baseseed, int *x0, int *lx);
//! constructor for constrained fine field
music_wnoise_generator( music_wnoise_generator<T>& rc, unsigned cubesize, long baseseed, int *x0_=NULL, int *lx_=NULL, bool zeromean=true );
music_wnoise_generator(music_wnoise_generator<T> &rc, unsigned cubesize, long baseseed, int *x0_ = NULL, int *lx_ = NULL, bool zeromean = true);
//! constructor
music_wnoise_generator( unsigned res, unsigned cubesize, long baseseed, bool zeromean=true );
music_wnoise_generator(unsigned res, unsigned cubesize, long baseseed, bool zeromean = true);
//! constructor to read white noise from file
music_wnoise_generator( unsigned res, std::string randfname, bool rndsign );
music_wnoise_generator(unsigned res, std::string randfname, bool rndsign);
//! copy constructor for averaged field (not copying) hence explicit!
explicit music_wnoise_generator( /*const*/ music_wnoise_generator <T>& rc );
explicit music_wnoise_generator(/*const*/ music_wnoise_generator<T> &rc);
//! destructor
~music_wnoise_generator()
{
for( unsigned i=0; i<rnums_.size(); ++i )
if( rnums_[i] != NULL )
delete rnums_[i];
for (unsigned i = 0; i < rnums_.size(); ++i)
if (rnums_[i] != NULL)
delete rnums_[i];
rnums_.clear();
}
//! access a random number, this allocates a cube and fills it with consistent random numbers
inline T& operator()( int i, int j, int k, bool fillrand=true )
inline T &operator()(int i, int j, int k, bool fillrand = true)
{
int ic, jc, kc, is, js, ks;
if( ncubes_ == 0 )
if (ncubes_ == 0)
throw std::runtime_error("music_wnoise_generator: internal error, not properly initialized");
//... determine cube
ic = (int)((double)i/cubesize_ + ncubes_) % ncubes_;
jc = (int)((double)j/cubesize_ + ncubes_) % ncubes_;
kc = (int)((double)k/cubesize_ + ncubes_) % ncubes_;
size_t icube = ((size_t)ic*ncubes_+(size_t)jc)*ncubes_+(size_t)kc;
cubemap_iterator it = cubemap_.find( icube );
if( it == cubemap_.end() )
{
LOGERR("Attempting to copy data from non-existing RND cube %d,%d,%d @ %d,%d,%d",ic,jc,kc,i,j,k);
throw std::runtime_error("attempting to copy data from non-existing RND cube");
}
ic = (int)((double)i / cubesize_ + ncubes_) % ncubes_;
jc = (int)((double)j / cubesize_ + ncubes_) % ncubes_;
kc = (int)((double)k / cubesize_ + ncubes_) % ncubes_;
size_t icube = ((size_t)ic * ncubes_ + (size_t)jc) * ncubes_ + (size_t)kc;
cubemap_iterator it = cubemap_.find(icube);
if (it == cubemap_.end())
{
LOGERR("Attempting to copy data from non-existing RND cube %d,%d,%d @ %d,%d,%d", ic, jc, kc, i, j, k);
throw std::runtime_error("attempting to copy data from non-existing RND cube");
}
size_t cubeidx = it->second;
if( rnums_[ cubeidx ] == NULL )
{
LOGERR("Attempting to access data from non-allocated RND cube %d,%d,%d",ic,jc,kc);
throw std::runtime_error("attempting to access data from non-allocated RND cube");
}
if (rnums_[cubeidx] == NULL)
{
LOGERR("Attempting to access data from non-allocated RND cube %d,%d,%d", ic, jc, kc);
throw std::runtime_error("attempting to access data from non-allocated RND cube");
}
//... determine cell in cube
is = (i - ic * cubesize_ + cubesize_) % cubesize_;
js = (j - jc * cubesize_ + cubesize_) % cubesize_;
ks = (k - kc * cubesize_ + cubesize_) % cubesize_;
return (*rnums_[ cubeidx ])(is,js,ks);
return (*rnums_[cubeidx])(is, js, ks);
}
//! free all cubes
void free_all_mem( void )
void free_all_mem(void)
{
for( unsigned i=0; i<rnums_.size(); ++i )
if( rnums_[i] != NULL )
{
delete rnums_[i];
rnums_[i] = NULL;
}
}
for (unsigned i = 0; i < rnums_.size(); ++i)
if (rnums_[i] != NULL)
{
delete rnums_[i];
rnums_[i] = NULL;
}
}
};
#endif //__MUSIC_WNOISE_GENERATOR_HH