-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathconstraints.hh
422 lines (321 loc) · 9.84 KB
/
constraints.hh
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
/*
constraints.hh - This file is part of MUSIC -
a code to generate multi-scale initial conditions
for cosmological simulations
Copyright (C) 2010 Oliver Hahn
*/
#ifndef __CONSTRAINTS_HH
#define __CONSTRAINTS_HH
#include <vector>
#include <complex>
#include <gsl/gsl_linalg.h>
#include "general.hh"
#include "config_file.hh"
#include "transfer_function.hh"
#include "cosmology.hh"
//! matrix class serving as a gsl wrapper
class matrix
{
protected:
gsl_matrix * m_;
//double *data_;
size_t M_, N_;
public:
matrix( size_t M, size_t N )
: M_(M), N_(N)
{
m_ = gsl_matrix_alloc(M_,N_);
}
matrix( size_t N )
: M_(N), N_(N)
{
m_ = gsl_matrix_alloc(M_,N_);
}
matrix( const matrix& o )
{
M_ = o.M_;
N_ = o.N_;
m_ = gsl_matrix_alloc(M_,N_);
gsl_matrix_memcpy(m_, o.m_ );
}
~matrix()
{
gsl_matrix_free( m_ );
}
double& operator()( size_t i, size_t j )
{ return *gsl_matrix_ptr( m_, i, j ); }
const double& operator()( size_t i, size_t j ) const
{ return *gsl_matrix_const_ptr( m_, i, j ); }
matrix& operator=( const matrix& o )
{
gsl_matrix_free( m_ );
M_ = o.M_;
N_ = o.N_;
m_ = gsl_matrix_alloc(M_,N_);
gsl_matrix_memcpy(m_, o.m_ );
return *this;
}
matrix& invert()
{
if( M_!=N_ )
throw std::runtime_error("Attempt to invert a non-square matrix!");
int s;
gsl_matrix* im = gsl_matrix_alloc(M_,N_);
gsl_permutation * p = gsl_permutation_alloc (M_);
gsl_linalg_LU_decomp( m_, p, &s );
gsl_linalg_LU_invert( m_, p, im );
gsl_matrix_memcpy(m_, im);
gsl_permutation_free(p);
gsl_matrix_free(im);
return *this;
}
};
//! class to impose constraints on the white noise field (van de Weygaert & Bertschinger 1996)
class constraint_set
{
public:
enum constr_type{ halo, peak };
protected:
struct constraint{
constr_type type;
double x,y,z;
double gx,gy,gz;
double Rg, Rg2;
double gRg, gRg2;
double sigma;
};
config_file *pcf_;
std::vector<constraint> cset_;
transfer_function *ptf_;
CosmoCalc *pccalc_;
Cosmology *pcosmo_;
double dplus0_;
unsigned constr_level_;
inline std::complex<double> eval_constr( size_t icon, double kx, double ky, double kz )
{
double re, im, kdotx, k2;
kdotx = cset_[icon].gx*kx+cset_[icon].gy*ky+cset_[icon].gz*kz;
k2 = kx*kx+ky*ky+kz*kz;
re = im = exp(-k2*cset_[icon].gRg2/2.0);
re *= cos( kdotx );
im *= sin( kdotx );
return std::complex<double>(re,im);
}
#if defined(FFTW3) && defined(SINGLE_PRECISION)
//! apply constraints to the white noise
void wnoise_constr_corr( double dx, size_t nx, size_t ny, size_t nz, std::vector<double>& g0, matrix& cinv, fftwf_complex* cw );
//! measure sigma for each constraint in the unconstrained noise
void wnoise_constr_corr( double dx, fftwf_complex* cw, size_t nx, size_t ny, size_t nz, std::vector<double>& g0 );
#else
//! apply constraints to the white noise
void wnoise_constr_corr( double dx, size_t nx, size_t ny, size_t nz, std::vector<double>& g0, matrix& cinv, fftw_complex* cw );
//! measure sigma for each constraint in the unconstrained noise
void wnoise_constr_corr( double dx, fftw_complex* cw, size_t nx, size_t ny, size_t nz, std::vector<double>& g0 );
#endif
//! compute the covariance between the constraints
void icov_constr( double dx, size_t nx, size_t ny, size_t nz, matrix& cij );
public:
//! constructor
constraint_set( config_file& cf, transfer_function *ptf );
//! destructor
~constraint_set()
{
delete pccalc_;
delete pcosmo_;
}
template< typename rng >
void apply( unsigned ilevel, int x0[], int lx[], rng* wnoise )
{
if( cset_.size() == 0 || constr_level_ != ilevel )
return;
unsigned nlvl = 1<<ilevel;
double boxlength = pcf_->getValue<double>("setup","boxlength");
//... compute constraint coordinates for grid
for( size_t i=0; i<cset_.size(); ++i )
{
cset_[i].gx = cset_[i].x * (double)nlvl;
cset_[i].gy = cset_[i].y * (double)nlvl;
cset_[i].gz = cset_[i].z * (double)nlvl;
cset_[i].gRg = cset_[i].Rg/boxlength * (double)nlvl;
cset_[i].gRg2 = cset_[i].gRg*cset_[i].gRg;
if(cset_[i].gRg > 0.5*lx[0])
LOGWARN("Constraint %d appears to be too large scale",i);
}
std::vector<double> g0;
// unsigned levelmax = pcf_->getValue<unsigned>("setup","levelmax");
unsigned levelmin = pcf_->getValue<unsigned>("setup","levelmin_TF");
bool bperiodic = ilevel==levelmin;
double dx = pcf_->getValue<double>("setup","boxlength")/(1<<ilevel);
LOGINFO("Computing constrained realization...");
if( bperiodic )
{
//... we are operating on the periodic coarse grid
size_t nx = lx[0], ny = lx[1], nz = lx[2], nzp = nz+2;
fftw_real * w = new fftw_real[nx*ny*nzp];
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_complex * cw = reinterpret_cast<fftwf_complex*> (w);
fftwf_plan p = fftwf_plan_dft_r2c_3d( nx, ny, nz, w, cw, FFTW_ESTIMATE),
ip = fftwf_plan_dft_c2r_3d( nx, ny, nz, cw, w, FFTW_ESTIMATE);
#else
fftw_complex * cw = reinterpret_cast<fftw_complex*> (w);
fftw_plan p = fftw_plan_dft_r2c_3d( nx, ny, nz, w, cw, FFTW_ESTIMATE),
ip = fftw_plan_dft_c2r_3d( nx, ny, nz, cw, w, FFTW_ESTIMATE);
#endif
#else
fftw_complex * cw = reinterpret_cast<fftw_complex*> (w);
rfftwnd_plan p = rfftw3d_create_plan( nx, ny, nz, FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE),
ip = rfftw3d_create_plan( nx, ny, nz, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE|FFTW_IN_PLACE);
#endif
double fftnorm = 1.0/sqrt(nx*ny*nz);
#pragma omp parallel for
for( int i=0; i<(int)nx; i++ )
for( int j=0; j<(int)ny; j++ )
for( int k=0; k<(int)nz; k++ )
{
size_t q = ((size_t)i*ny+(size_t)j)*nzp+(size_t)k;
w[q] = (*wnoise)((x0[0]+i)%nx,(x0[1]+j)%ny,(x0[2]+k)%nz)*fftnorm;
}
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_execute( p );
#else
fftw_execute( p );
#endif
#else
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_real_to_complex( omp_get_max_threads(), p, w, NULL );
#else
rfftwnd_one_real_to_complex( p, w, NULL );
#endif
#endif
wnoise_constr_corr( dx, cw, nx, ny, nz, g0 );
matrix c(2,2);
icov_constr( dx, nx, ny, nz, c );
wnoise_constr_corr( dx, nx, ny, nz, g0, c, cw );
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_execute( ip );
#else
fftw_execute( ip );
#endif
#else
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_complex_to_real( omp_get_max_threads(), ip, cw, NULL );
#else
rfftwnd_one_complex_to_real( ip, cw, NULL );
#endif
#endif
#pragma omp parallel for
for( int i=0; i<(int)nx; i++ )
for( int j=0; j<(int)ny; j++ )
for( int k=0; k<(int)nz; k++ )
{
size_t q = ((size_t)i*ny+(size_t)j)*nzp+(size_t)k;
(*wnoise)((x0[0]+i),(x0[1]+j),(x0[2]+k)) = w[q]*fftnorm;
}
LOGINFO("Applied constraints to level %d.",ilevel);
delete[] w;
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_destroy_plan(p);
#else
fftw_destroy_plan(p);
#endif
#else
fftwnd_destroy_plan(p);
#endif
}else{
//... we are operating on a refinement grid, not necessarily the finest
size_t nx = lx[0], ny = lx[1], nz = lx[2], nzp = nz+2;
fftw_real * w = new fftw_real[nx*ny*nzp];
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_complex * cw = reinterpret_cast<fftwf_complex*> (w);
fftwf_plan p = fftwf_plan_dft_r2c_3d( nx, ny, nz, w, cw, FFTW_ESTIMATE),
ip = fftwf_plan_dft_c2r_3d( nx, ny, nz, cw, w, FFTW_ESTIMATE);
#else
fftw_complex * cw = reinterpret_cast<fftw_complex*> (w);
fftw_plan p = fftw_plan_dft_r2c_3d( nx, ny, nz, w, cw, FFTW_ESTIMATE),
ip = fftw_plan_dft_c2r_3d( nx, ny, nz, cw, w, FFTW_ESTIMATE);
#endif
#else
fftw_complex * cw = reinterpret_cast<fftw_complex*> (w);
rfftwnd_plan p = rfftw3d_create_plan( nx, ny, nz, FFTW_REAL_TO_COMPLEX, FFTW_ESTIMATE|FFTW_IN_PLACE),
ip = rfftw3d_create_plan( nx, ny, nz, FFTW_COMPLEX_TO_REAL, FFTW_ESTIMATE|FFTW_IN_PLACE);
#endif
double fftnorm = 1.0/sqrt(nx*ny*nz);
int il = nx/4, ir = 3*nx/4, jl=ny/4, jr = 3*ny/4, kl = nz/4, kr = 3*nz/4;
#pragma omp parallel for
for( int i=0; i<(int)nx; i++ )
for( int j=0; j<(int)ny; j++ )
for( int k=0; k<(int)nz; k++ )
{
size_t q = ((size_t)i*ny+(size_t)j)*nzp+(size_t)k;
if( i>=il && i<ir && j>=jl && j<jr && k>=kl && k<kr )
w[q] = (*wnoise)((x0[0]+i),(x0[1]+j),(x0[2]+k))*fftnorm;
else
w[q] = 0.0;
}
int nlvl05 = 1<<(ilevel-1);
int xs = nlvl05-x0[0], ys = nlvl05-x0[1], zs = nlvl05-x0[2];
for( size_t i=0; i<cset_.size(); ++i )
{
cset_[i].gx -= xs;
cset_[i].gy -= ys;
cset_[i].gz -= zs;
}
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_execute( p );
#else
fftw_execute( p );
#endif
#else
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_real_to_complex( omp_get_max_threads(), p, w, NULL );
#else
rfftwnd_one_real_to_complex( p, w, NULL );
#endif
#endif
wnoise_constr_corr( dx, cw, nx, ny, nz, g0 );
matrix c(2,2);
icov_constr( dx, nx, ny, nz, c );
wnoise_constr_corr( dx, nx, ny, nz, g0, c, cw );
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_execute( ip );
#else
fftw_execute( ip );
#endif
#else
#ifndef SINGLETHREAD_FFTW
rfftwnd_threads_one_complex_to_real( omp_get_max_threads(), ip, cw, NULL );
#else
rfftwnd_one_complex_to_real( ip, cw, NULL );
#endif
#endif
#pragma omp parallel for
for( int i=0; i<(int)nx; i++ )
for( int j=0; j<(int)ny; j++ )
for( int k=0; k<(int)nz; k++ )
{
size_t q = ((size_t)i*ny+(size_t)j)*nzp+(size_t)k;
if( i>=il && i<ir && j>=jl && j<jr && k>=kl && k<kr )
(*wnoise)((x0[0]+i),(x0[1]+j),(x0[2]+k)) = w[q]*fftnorm;
}
LOGINFO("Applied constraints to level %d.",ilevel);
delete[] w;
#ifdef FFTW3
#ifdef SINGLE_PRECISION
fftwf_destroy_plan(p);
#else
fftw_destroy_plan(p);
#endif
#else
fftwnd_destroy_plan(p);
#endif
}
}
};
#endif // __CONSTRAINTS_HH