1 | #include "machdefs.h"
|
---|
2 | #include <iostream>
|
---|
3 | #include <stdlib.h>
|
---|
4 | #include <stdio.h>
|
---|
5 | #include <string.h>
|
---|
6 | #include <math.h>
|
---|
7 | #include <unistd.h>
|
---|
8 |
|
---|
9 | #include "tarray.h"
|
---|
10 | #include "pexceptions.h"
|
---|
11 | #include "perandom.h"
|
---|
12 | #include "srandgen.h"
|
---|
13 |
|
---|
14 | #include "fabtcolread.h"
|
---|
15 | #include "fabtwriter.h"
|
---|
16 | #include "fioarr.h"
|
---|
17 | #include "ntuple.h"
|
---|
18 |
|
---|
19 | #include "arrctcast.h"
|
---|
20 |
|
---|
21 | #include "constcosmo.h"
|
---|
22 | #include "geneutils.h"
|
---|
23 | #include "schechter.h"
|
---|
24 |
|
---|
25 | #include "genefluct3d.h"
|
---|
26 |
|
---|
27 | #if defined(GEN3D_FLOAT)
|
---|
28 | #define GEN3D_FFTW_INIT_THREADS fftwf_init_threads
|
---|
29 | #define GEN3D_FFTW_CLEANUP_THREADS fftwf_cleanup_threads
|
---|
30 | #define GEN3D_FFTW_PLAN_WITH_NTHREADS fftwf_plan_with_nthreads
|
---|
31 | #define GEN3D_FFTW_PLAN_DFT_R2C_3D fftwf_plan_dft_r2c_3d
|
---|
32 | #define GEN3D_FFTW_PLAN_DFT_C2R_3D fftwf_plan_dft_c2r_3d
|
---|
33 | #define GEN3D_FFTW_DESTROY_PLAN fftwf_destroy_plan
|
---|
34 | #define GEN3D_FFTW_EXECUTE fftwf_execute
|
---|
35 | #else
|
---|
36 | #define GEN3D_FFTW_INIT_THREADS fftw_init_threads
|
---|
37 | #define GEN3D_FFTW_CLEANUP_THREADS fftw_cleanup_threads
|
---|
38 | #define GEN3D_FFTW_PLAN_WITH_NTHREADS fftw_plan_with_nthreads
|
---|
39 | #define GEN3D_FFTW_PLAN_DFT_R2C_3D fftw_plan_dft_r2c_3d
|
---|
40 | #define GEN3D_FFTW_PLAN_DFT_C2R_3D fftw_plan_dft_c2r_3d
|
---|
41 | #define GEN3D_FFTW_DESTROY_PLAN fftw_destroy_plan
|
---|
42 | #define GEN3D_FFTW_EXECUTE fftw_execute
|
---|
43 | #endif
|
---|
44 |
|
---|
45 | #define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag()))
|
---|
46 |
|
---|
47 | namespace SOPHYA {
|
---|
48 |
|
---|
49 | //-------------------------------------------------------
|
---|
50 | GeneFluct3D::GeneFluct3D(long nx,long ny,long nz,double dx,double dy,double dz
|
---|
51 | ,unsigned short nthread,int lp)
|
---|
52 | {
|
---|
53 | init_default();
|
---|
54 |
|
---|
55 | lp_ = lp;
|
---|
56 | nthread_ = nthread;
|
---|
57 |
|
---|
58 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
59 | setalloc();
|
---|
60 | setpointers(false);
|
---|
61 | init_fftw();
|
---|
62 | }
|
---|
63 |
|
---|
64 | GeneFluct3D::GeneFluct3D(unsigned short nthread)
|
---|
65 | {
|
---|
66 | init_default();
|
---|
67 | setsize(2,2,2,1.,1.,1.);
|
---|
68 | nthread_ = nthread;
|
---|
69 | setalloc();
|
---|
70 | setpointers(false);
|
---|
71 | init_fftw();
|
---|
72 | }
|
---|
73 |
|
---|
74 | GeneFluct3D::~GeneFluct3D(void)
|
---|
75 | {
|
---|
76 | delete_fftw();
|
---|
77 | }
|
---|
78 |
|
---|
79 | //-------------------------------------------------------
|
---|
80 | void GeneFluct3D::init_default(void)
|
---|
81 | {
|
---|
82 | Nx_ = Ny_ = Nz_ = 0;
|
---|
83 | is_set_fft_plan = false;
|
---|
84 | nthread_ = 0;
|
---|
85 | lp_ = 0;
|
---|
86 | array_allocated_ = false; array_type = 0;
|
---|
87 | cosmo_ = NULL;
|
---|
88 | growth_ = NULL;
|
---|
89 | good_dzinc_ = 0.01;
|
---|
90 | compute_pk_redsh_ref_ = -999.;
|
---|
91 | redsh_ref_ = -999.;
|
---|
92 | kredsh_ref_ = 0.;
|
---|
93 | dred_ref_ = -999.;
|
---|
94 | h_ref_ = -999.;
|
---|
95 | loscom_ref_ = -999.;
|
---|
96 | dtrc_ref_ = dlum_ref_ = dang_ref_ = -999.;
|
---|
97 | nu_ref_ = dnu_ref_ = -999.;
|
---|
98 | loscom_min_ = loscom_max_ = -999.;
|
---|
99 | loscom2zred_min_ = loscom2zred_max_ = 0.;
|
---|
100 | xobs_[0] = xobs_[1] = xobs_[2] = 0.;
|
---|
101 | zred_.resize(0);
|
---|
102 | loscom_.resize(0);
|
---|
103 | loscom2zred_.resize(0);
|
---|
104 | }
|
---|
105 |
|
---|
106 | void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz)
|
---|
107 | {
|
---|
108 | if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<<nx<<","<<ny<<","<<nz
|
---|
109 | <<" D="<<dx<<","<<dy<<","<<dz<<endl;
|
---|
110 | if(nx<=0 || dx<=0.) {
|
---|
111 | const char *bla = "GeneFluct3D::setsize_Error: bad value(s) for nn/dx";
|
---|
112 | cout<<bla<<endl; throw ParmError(bla);
|
---|
113 | }
|
---|
114 |
|
---|
115 | // Les tailles des tableaux
|
---|
116 | Nx_ = nx;
|
---|
117 | Ny_ = ny; if(Ny_ <= 0) Ny_ = Nx_;
|
---|
118 | Nz_ = nz; if(Nz_ <= 0) Nz_ = Nx_;
|
---|
119 | N_.resize(0); N_.push_back(Nx_); N_.push_back(Ny_); N_.push_back(Nz_);
|
---|
120 | NRtot_ = (int_8)Nx_*(int_8)Ny_*(int_8)Nz_; // nombre de pixels dans le survey
|
---|
121 | NCz_ = Nz_/2 +1;
|
---|
122 | NFcoef_ = (int_8)Nx_*(int_8)Ny_*(int_8)NCz_; // nombre de coeff de Fourier dans le survey
|
---|
123 | NTz_ = 2*NCz_;
|
---|
124 |
|
---|
125 | // le pas dans l'espace (Mpc)
|
---|
126 | Dx_ = dx;
|
---|
127 | Dy_ = dy; if(Dy_ <= 0.) Dy_ = Dx_;
|
---|
128 | Dz_ = dz; if(Dz_ <= 0.) Dz_ = Dx_;
|
---|
129 | D_.resize(0); D_.push_back(Dx_); D_.push_back(Dy_); D_.push_back(Dz_);
|
---|
130 | dVol_ = Dx_*Dy_*Dz_;
|
---|
131 | Vol_ = (Nx_*Dx_)*(Ny_*Dy_)*(Nz_*Dz_);
|
---|
132 |
|
---|
133 | // Le pas dans l'espace de Fourier (Mpc^-1)
|
---|
134 | Dkx_ = 2.*M_PI/(Nx_*Dx_);
|
---|
135 | Dky_ = 2.*M_PI/(Ny_*Dy_);
|
---|
136 | Dkz_ = 2.*M_PI/(Nz_*Dz_);
|
---|
137 | Dk_.resize(0); Dk_.push_back(Dkx_); Dk_.push_back(Dky_); Dk_.push_back(Dkz_);
|
---|
138 | Dk3_ = Dkx_*Dky_*Dkz_;
|
---|
139 |
|
---|
140 | // La frequence de Nyquist en k (Mpc^-1)
|
---|
141 | Knyqx_ = M_PI/Dx_;
|
---|
142 | Knyqy_ = M_PI/Dy_;
|
---|
143 | Knyqz_ = M_PI/Dz_;
|
---|
144 | Knyq_.resize(0); Knyq_.push_back(Knyqx_); Knyq_.push_back(Knyqy_); Knyq_.push_back(Knyqz_);
|
---|
145 | }
|
---|
146 |
|
---|
147 | void GeneFluct3D::setalloc(void)
|
---|
148 | {
|
---|
149 | #if defined(GEN3D_FLOAT)
|
---|
150 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc FLOAT ---"<<endl;
|
---|
151 | #else
|
---|
152 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc DOUBLE ---"<<endl;
|
---|
153 | #endif
|
---|
154 | // Dimensionnement du tableau complex<r_8>
|
---|
155 | // ATTENTION: TArray adresse en memoire a l'envers du C
|
---|
156 | // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
|
---|
157 | sa_size_t SzK_[3] = {NCz_,Ny_,Nx_}; // a l'envers
|
---|
158 | try {
|
---|
159 | T_.ReSize(3,SzK_);
|
---|
160 | array_allocated_ = true; array_type=0;
|
---|
161 | if(lp_>1) cout<<" allocating: "<<T_.Size()*sizeof(complex<GEN3D_TYPE>)/1.e6<<" Mo"<<endl;
|
---|
162 | } catch (...) {
|
---|
163 | cout<<"GeneFluct3D::setalloc_Error: Problem allocating T_"<<endl;
|
---|
164 | }
|
---|
165 | }
|
---|
166 |
|
---|
167 | void GeneFluct3D::setpointers(bool from_real)
|
---|
168 | {
|
---|
169 | if(lp_>1) cout<<"--- GeneFluct3D::setpointers ---"<<endl;
|
---|
170 | if(from_real) T_ = ArrCastR2C(R_);
|
---|
171 | else R_ = ArrCastC2R(T_);
|
---|
172 | // On remplit les pointeurs
|
---|
173 | fdata_ = (GEN3D_FFTW_COMPLEX *) (&T_(0,0,0));
|
---|
174 | data_ = (GEN3D_TYPE *) (&R_(0,0,0));
|
---|
175 | }
|
---|
176 |
|
---|
177 | void GeneFluct3D::init_fftw(void)
|
---|
178 | {
|
---|
179 | if( is_set_fft_plan ) delete_fftw();
|
---|
180 |
|
---|
181 | // --- Initialisation de fftw3 (attention data est sur-ecrit a l'init)
|
---|
182 | if(lp_>1) cout<<"--- GeneFluct3D::init_fftw ---"<<endl;
|
---|
183 | #ifdef WITH_FFTW_THREAD
|
---|
184 | if(nthread_>0) {
|
---|
185 | cout<<"...Computing with "<<nthread_<<" threads"<<endl;
|
---|
186 | GEN3D_FFTW_INIT_THREADS();
|
---|
187 | GEN3D_FFTW_PLAN_WITH_NTHREADS(nthread_);
|
---|
188 | }
|
---|
189 | #endif
|
---|
190 | if(lp_>1) cout<<"...forward plan"<<endl;
|
---|
191 | pf_ = GEN3D_FFTW_PLAN_DFT_R2C_3D(Nx_,Ny_,Nz_,data_,fdata_,FFTW_ESTIMATE);
|
---|
192 | if(lp_>1) cout<<"...backward plan"<<endl;
|
---|
193 | pb_ = GEN3D_FFTW_PLAN_DFT_C2R_3D(Nx_,Ny_,Nz_,fdata_,data_,FFTW_ESTIMATE);
|
---|
194 | is_set_fft_plan = true;
|
---|
195 | }
|
---|
196 |
|
---|
197 | void GeneFluct3D::delete_fftw(void)
|
---|
198 | {
|
---|
199 | if( !is_set_fft_plan ) return;
|
---|
200 | GEN3D_FFTW_DESTROY_PLAN(pf_);
|
---|
201 | GEN3D_FFTW_DESTROY_PLAN(pb_);
|
---|
202 | #ifdef WITH_FFTW_THREAD
|
---|
203 | if(nthread_>0) GEN3D_FFTW_CLEANUP_THREADS();
|
---|
204 | #endif
|
---|
205 | is_set_fft_plan = false;
|
---|
206 | }
|
---|
207 |
|
---|
208 | void GeneFluct3D::check_array_alloc(void)
|
---|
209 | // Pour tester si le tableau T_ est alloue
|
---|
210 | {
|
---|
211 | if(array_allocated_) return;
|
---|
212 | char bla[90];
|
---|
213 | sprintf(bla,"GeneFluct3D::check_array_alloc_Error: array is not allocated");
|
---|
214 | cout<<bla<<endl; throw ParmError(bla);
|
---|
215 | }
|
---|
216 |
|
---|
217 | //-------------------------------------------------------
|
---|
218 | void GeneFluct3D::SetObservator(double redshref,double kredshref)
|
---|
219 | // L'observateur est au redshift z=0
|
---|
220 | // est situe sur la "perpendiculaire" a la face x,y
|
---|
221 | // issue du centre de cette face
|
---|
222 | // Il faut positionner le cube sur l'axe des z cad des redshifts:
|
---|
223 | // redshref = redshift de reference
|
---|
224 | // Si redshref<0 alors redshref=0
|
---|
225 | // kredshref = indice (en double) correspondant a ce redshift
|
---|
226 | // Si kredshref<0 alors kredshref=nz/2 (milieu du cube)
|
---|
227 | // Exemple: redshref=1.5 kredshref=250.75
|
---|
228 | // -> Le pixel i=nx/2 j=ny/2 k=250.75 est au redshift 1.5
|
---|
229 | {
|
---|
230 | if(redshref<0.) redshref = 0.;
|
---|
231 | if(kredshref<0.) {
|
---|
232 | if(Nz_<=0) {
|
---|
233 | const char *bla = "GeneFluct3D::SetObservator_Error: for kredsh_ref<0 define cube geometry first";
|
---|
234 | cout<<bla<<endl; throw ParmError(bla);
|
---|
235 | }
|
---|
236 | kredshref = Nz_/2.;
|
---|
237 | }
|
---|
238 | redsh_ref_ = redshref;
|
---|
239 | kredsh_ref_ = kredshref;
|
---|
240 | if(lp_>0)
|
---|
241 | cout<<"--- GeneFluct3D::SetObservator zref="<<redsh_ref_<<" kref="<<kredsh_ref_<<endl;
|
---|
242 | }
|
---|
243 |
|
---|
244 | void GeneFluct3D::SetCosmology(CosmoCalc& cosmo)
|
---|
245 | {
|
---|
246 | cosmo_ = &cosmo;
|
---|
247 | if(lp_>1) cosmo_->Print();
|
---|
248 | }
|
---|
249 |
|
---|
250 | void GeneFluct3D::SetGrowthFactor(GrowthFactor& growth)
|
---|
251 | {
|
---|
252 | growth_ = &growth;
|
---|
253 | }
|
---|
254 |
|
---|
255 | long GeneFluct3D::LosComRedshift(double zinc,long npoints)
|
---|
256 | // Given a position of the cube relative to the observer
|
---|
257 | // and a cosmology
|
---|
258 | // (SetObservator() and SetCosmology() should have been called !)
|
---|
259 | // This routine filled:
|
---|
260 | // the vector "zred_" of scanned redshift (by zinc increments)
|
---|
261 | // the vector "loscom_" of corresponding los comoving distance
|
---|
262 | // -- Input:
|
---|
263 | // zinc : redshift increment for computation
|
---|
264 | // npoints : number of points required for inverting loscom -> zred
|
---|
265 | //
|
---|
266 | {
|
---|
267 | if(zinc<=0.) zinc = good_dzinc_;
|
---|
268 | if(lp_>0) cout<<"--- LosComRedshift: zinc="<<zinc<<" , npoints="<<npoints<<endl;
|
---|
269 |
|
---|
270 | if(cosmo_ == NULL || growth_==NULL || redsh_ref_<0.) {
|
---|
271 | const char *bla = "GeneFluct3D::LosComRedshift_Error: set Observator, Cosmology and Growth first";
|
---|
272 | cout<<bla<<endl; throw ParmError(bla);
|
---|
273 | }
|
---|
274 |
|
---|
275 | // La distance angulaire/luminosite/Dnu au pixel de reference
|
---|
276 | dred_ref_ = Dz_/(cosmo_->Dhubble()/cosmo_->E(redsh_ref_));
|
---|
277 | loscom_ref_ = cosmo_->Dloscom(redsh_ref_);
|
---|
278 | dtrc_ref_ = cosmo_->Dtrcom(redsh_ref_);
|
---|
279 | dlum_ref_ = cosmo_->Dlum(redsh_ref_);
|
---|
280 | dang_ref_ = cosmo_->Dang(redsh_ref_);
|
---|
281 | h_ref_ = cosmo_->H(redsh_ref_);
|
---|
282 | growth_ref_ = (*growth_)(redsh_ref_);
|
---|
283 | dsdz_growth_ref_ = growth_->DsDz(redsh_ref_,zinc);
|
---|
284 | nu_ref_ = Fr_HyperFin_Par/(1.+redsh_ref_); // GHz
|
---|
285 | dnu_ref_ = Fr_HyperFin_Par *dred_ref_/pow(1.+redsh_ref_,2.); // GHz
|
---|
286 | if(lp_>0) {
|
---|
287 | cout<<"...reference pixel redshref="<<redsh_ref_
|
---|
288 | <<", dredref="<<dred_ref_
|
---|
289 | <<", nuref="<<nu_ref_ <<" GHz"
|
---|
290 | <<", dnuref="<<dnu_ref_ <<" GHz"<<endl
|
---|
291 | <<" H="<<h_ref_<<" km/s/Mpc"
|
---|
292 | <<", D="<<growth_ref_
|
---|
293 | <<", dD/dz="<<dsdz_growth_ref_<<endl
|
---|
294 | <<" dlosc="<<loscom_ref_<<" Mpc com"
|
---|
295 | <<", dtrc="<<dtrc_ref_<<" Mpc com"
|
---|
296 | <<", dlum="<<dlum_ref_<<" Mpc"
|
---|
297 | <<", dang="<<dang_ref_<<" Mpc"<<endl;
|
---|
298 | }
|
---|
299 |
|
---|
300 | // On calcule les coordonnees de l'observateur dans le repere du cube
|
---|
301 | // cad dans le repere ou l'origine est au centre du pixel i=j=l=0.
|
---|
302 | // L'observateur est sur un axe centre sur le milieu de la face Oxy
|
---|
303 | xobs_[0] = Nx_/2.*Dx_;
|
---|
304 | xobs_[1] = Ny_/2.*Dy_;
|
---|
305 | xobs_[2] = kredsh_ref_*Dz_ - loscom_ref_;
|
---|
306 |
|
---|
307 | // L'observateur est-il dans le cube?
|
---|
308 | bool obs_in_cube = false;
|
---|
309 | if(xobs_[2]>=0. && xobs_[2]<=Nz_*Dz_) obs_in_cube = true;
|
---|
310 |
|
---|
311 | // Find MINIMUM los com distance to the observer:
|
---|
312 | // c'est le centre de la face a k=0
|
---|
313 | // (ou zero si l'observateur est dans le cube)
|
---|
314 | loscom_min_ = 0.;
|
---|
315 | if(!obs_in_cube) loscom_min_ = -xobs_[2];
|
---|
316 |
|
---|
317 | // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED
|
---|
318 | if(loscom_min_<=1.e-50)
|
---|
319 | for(int i=0;i<10;i++)
|
---|
320 | cout<<"ATTENTION TOUTES LES PARTIES DU CODE NE MARCHENT PAS POUR UN OBSERVATEUR DANS LE CUBE"<<endl;
|
---|
321 | // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED
|
---|
322 |
|
---|
323 |
|
---|
324 | // Find MAXIMUM los com distance to the observer:
|
---|
325 | // ou que soit positionne l'observateur, la distance
|
---|
326 | // maximal est sur un des coins du cube
|
---|
327 | loscom_max_ = 0.;
|
---|
328 | for(long i=0;i<=1;i++) {
|
---|
329 | double dx2 = DXcom(i*(Nx_-1)); dx2 *= dx2;
|
---|
330 | for(long j=0;j<=1;j++) {
|
---|
331 | double dy2 = DYcom(j*(Ny_-1)); dy2 *= dy2;
|
---|
332 | for(long k=0;k<=1;k++) {
|
---|
333 | double dz2 = DZcom(k*(Nz_-1)); dz2 *= dz2;
|
---|
334 | dz2 = sqrt(dx2+dy2+dz2);
|
---|
335 | if(dz2>loscom_max_) loscom_max_ = dz2;
|
---|
336 | }
|
---|
337 | }
|
---|
338 | }
|
---|
339 | if(lp_>0) {
|
---|
340 | cout<<"...zref="<<redsh_ref_<<" kzref="<<kredsh_ref_<<" losref="<<loscom_ref_<<" Mpc\n"
|
---|
341 | <<" xobs="<<xobs_[0]<<" , "<<xobs_[1]<<" , "<<xobs_[2]<<" Mpc "
|
---|
342 | <<" in_cube="<<obs_in_cube
|
---|
343 | <<" loscom_min="<<loscom_min_<<" loscom_max="<<loscom_max_<<" Mpc (com)"<<endl;
|
---|
344 | }
|
---|
345 |
|
---|
346 | // Fill the corresponding vectors for loscom and zred
|
---|
347 | // Be shure to have one dlc<loscom_min and one dlc>loscom_max
|
---|
348 | double zmin = 0., dlcmin=0.;
|
---|
349 | while(1) {
|
---|
350 | if(lp_>0)
|
---|
351 | cout<<"...Filling zred starting at zmin="<<zmin<<" with zinc="<<zinc
|
---|
352 | <<", loscom_min-max=["<<loscom_min_<<","<<loscom_max_<<"]"<<endl;
|
---|
353 | zred_.resize(0); loscom_.resize(0);
|
---|
354 | for(double z=zmin; ; z+=zinc) {
|
---|
355 | double dlc = cosmo_->Dloscom(z);
|
---|
356 | if(dlc<loscom_min_) {
|
---|
357 | dlcmin = dlc;
|
---|
358 | zmin = z;
|
---|
359 | zred_.resize(0); loscom_.resize(0);
|
---|
360 | }
|
---|
361 | zred_.push_back(z);
|
---|
362 | loscom_.push_back(dlc);
|
---|
363 | z += zinc;
|
---|
364 | if(dlc>loscom_max_) {
|
---|
365 | if(lp_>0)
|
---|
366 | cout<<" Min: z="<<zmin<<" dlc="<<dlcmin<<", Max: z="<<z<<" dlc="<<dlc<<endl;
|
---|
367 | break; // on sort apres avoir stoque un dlc>dlcmax
|
---|
368 | }
|
---|
369 | }
|
---|
370 | if(zred_.size()>=10) break;
|
---|
371 | zinc /= 10.;
|
---|
372 | cout<<" not enough points ("<<zred_.size()<<") for zref, retry with zinc="<<zinc<<endl;
|
---|
373 | }
|
---|
374 |
|
---|
375 | if(lp_>0) {
|
---|
376 | long n = zred_.size();
|
---|
377 | cout<<"...zred/loscom tables[zinc="<<zinc<<"]: n="<<n;
|
---|
378 | if(n>0) cout<<" z="<<zred_[0]<<" -> d="<<loscom_[0];
|
---|
379 | if(n>1) cout<<" , z="<<zred_[n-1]<<" -> d="<<loscom_[n-1];
|
---|
380 | cout<<endl;
|
---|
381 | }
|
---|
382 |
|
---|
383 | // Compute the parameters and tables needed for inversion loscom->zred
|
---|
384 | if(npoints<3) npoints = zred_.size();
|
---|
385 | InverseFunc invfun(zred_,loscom_);
|
---|
386 | invfun.ComputeParab(npoints,loscom2zred_);
|
---|
387 | loscom2zred_min_ = invfun.YMin();
|
---|
388 | loscom2zred_max_ = invfun.YMax();
|
---|
389 |
|
---|
390 | if(lp_>0) {
|
---|
391 | long n = loscom2zred_.size();
|
---|
392 | cout<<"...loscom -> zred[npoints="<<npoints<<"]: n="<<n
|
---|
393 | <<" los_min="<<loscom2zred_min_
|
---|
394 | <<" los_max="<<loscom2zred_max_
|
---|
395 | <<" -> zred=[";
|
---|
396 | if(n>0) cout<<loscom2zred_[0];
|
---|
397 | cout<<",";
|
---|
398 | if(n>1) cout<<loscom2zred_[n-1];
|
---|
399 | cout<<"]"<<endl;
|
---|
400 | if(lp_>1 && n>0)
|
---|
401 | for(int i=0;i<n;i++)
|
---|
402 | if(i<2 || abs(i-n/2)<2 || i>=n-2)
|
---|
403 | cout<<" i="<<i
|
---|
404 | <<" d="<<loscom2zred_min_+i*(loscom2zred_max_-loscom2zred_min_)/(n-1.)
|
---|
405 | <<" Mpc z="<<loscom2zred_[i]<<endl;
|
---|
406 | }
|
---|
407 |
|
---|
408 |
|
---|
409 | // Compute the table for D'(z)/D(z) where D'(z)=dD/dEta (Eta is conformal time)
|
---|
410 | zredmin_dpsd_ = zred_[0];
|
---|
411 | zredmax_dpsd_ = zred_[zred_.size()-1];
|
---|
412 | long nptd = long(sqrt(Nx_*Nx_ + Ny_*Ny_ +Nz_*Nz_));
|
---|
413 | if(nptd<10) nptd = 10;
|
---|
414 | good_dzinc_ = (zredmax_dpsd_ - zredmin_dpsd_)/nptd;
|
---|
415 | if(lp_>0) cout<<"...good_dzinc changed to "<<good_dzinc_<<endl;
|
---|
416 | dpsdfrzred_.resize(0);
|
---|
417 | double dz = (zredmax_dpsd_ - zredmin_dpsd_)/(nptd-1);
|
---|
418 | int nmod = nptd/5; if(nmod==0) nmod = 1;
|
---|
419 | if(lp_>0) cout<<"...Compute table D'/D on "<<nptd<<" pts for z["
|
---|
420 | <<zredmin_dpsd_<<","<<zredmax_dpsd_<<"]"<<endl;
|
---|
421 | for(int i=0;i<nptd;i++) {
|
---|
422 | double z = zredmin_dpsd_ + i*dz;
|
---|
423 | double om = cosmo_->Omatter(z);
|
---|
424 | double h = cosmo_->H(z);
|
---|
425 | double dsdz = growth_->DsDz(z,good_dzinc_);
|
---|
426 | double d = (*growth_)(z);
|
---|
427 | double v = -h * (1.+z) * dsdz / d;
|
---|
428 | dpsdfrzred_.push_back(v);
|
---|
429 | if(lp_ && (i%nmod==0 || i==nptd-1))
|
---|
430 | cout<<" z="<<z<<" H="<<h<<" Beta="<<-(1.+z)*dsdz/d
|
---|
431 | <<" (Om^0.6="<<pow(om,0.6)
|
---|
432 | <<") -H*(1+z)*D'/D="<<v<<" km/s/Mpc"<<endl;
|
---|
433 | }
|
---|
434 |
|
---|
435 | return zred_.size();
|
---|
436 | }
|
---|
437 |
|
---|
438 | //-------------------------------------------------------
|
---|
439 | void GeneFluct3D::WriteFits(string cfname,int bitpix)
|
---|
440 | {
|
---|
441 | cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
|
---|
442 | try {
|
---|
443 | FitsImg3DWriter fwrt(cfname.c_str(),bitpix,5);
|
---|
444 | fwrt.WriteKey("NX",Nx_," axe transverse 1"); // NAXIS3
|
---|
445 | fwrt.WriteKey("NY",Ny_," axe transverse 2"); // NAXIS2
|
---|
446 | fwrt.WriteKey("NZ",Nz_," axe longitudinal (redshift)"); // NAXIS1
|
---|
447 | fwrt.WriteKey("DX",Dx_," Mpc");
|
---|
448 | fwrt.WriteKey("DY",Dy_," Mpc");
|
---|
449 | fwrt.WriteKey("DZ",Dz_," Mpc");
|
---|
450 | fwrt.WriteKey("DKX",Dkx_," Mpc^-1");
|
---|
451 | fwrt.WriteKey("DKY",Dky_," Mpc^-1");
|
---|
452 | fwrt.WriteKey("DKZ",Dkz_," Mpc^-1");
|
---|
453 | fwrt.WriteKey("ZREFPK",compute_pk_redsh_ref_," Pk computed redshift");
|
---|
454 | fwrt.WriteKey("ZREF",redsh_ref_," reference redshift");
|
---|
455 | fwrt.WriteKey("KZREF",kredsh_ref_," reference redshift on z axe");
|
---|
456 | fwrt.WriteKey("HREF",h_ref_," ref H(z)");
|
---|
457 | fwrt.WriteKey("Growth",Dref()," growth at reference redshift");
|
---|
458 | fwrt.WriteKey("GrowthPk",DrefPk()," growth at Pk computed redshift");
|
---|
459 | fwrt.Write(R_);
|
---|
460 | } catch (PThrowable & exc) {
|
---|
461 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
462 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
463 | return;
|
---|
464 | } catch (...) {
|
---|
465 | cout<<" some other exception was caught !"<<endl;
|
---|
466 | return;
|
---|
467 | }
|
---|
468 | }
|
---|
469 |
|
---|
470 | void GeneFluct3D::ReadFits(string cfname)
|
---|
471 | {
|
---|
472 | cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
|
---|
473 | try {
|
---|
474 | FitsImg3DRead fimg(cfname.c_str(),0,5);
|
---|
475 | fimg.Read(R_);
|
---|
476 | long nx = fimg.ReadKeyL("NX"); // NAXIS3
|
---|
477 | long ny = fimg.ReadKeyL("NY"); // NAXIS2
|
---|
478 | long nz = fimg.ReadKeyL("NZ"); // NAXIS1
|
---|
479 | double dx = fimg.ReadKey("DX");
|
---|
480 | double dy = fimg.ReadKey("DY");
|
---|
481 | double dz = fimg.ReadKey("DZ");
|
---|
482 | double pkzref = fimg.ReadKey("ZREFPK");
|
---|
483 | double zref = fimg.ReadKey("ZREF");
|
---|
484 | double kzref = fimg.ReadKey("KZREF");
|
---|
485 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
486 | setpointers(true);
|
---|
487 | init_fftw();
|
---|
488 | SetObservator(zref,kzref);
|
---|
489 | array_allocated_ = true;
|
---|
490 | compute_pk_redsh_ref_ = pkzref;
|
---|
491 | } catch (PThrowable & exc) {
|
---|
492 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
493 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
494 | return;
|
---|
495 | } catch (...) {
|
---|
496 | cout<<" some other exception was caught !"<<endl;
|
---|
497 | return;
|
---|
498 | }
|
---|
499 | }
|
---|
500 |
|
---|
501 | void GeneFluct3D::WritePPF(string cfname,bool write_real)
|
---|
502 | // On ecrit soit le TArray<r_8> ou le TArray<complex <r_8> >
|
---|
503 | {
|
---|
504 | cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
|
---|
505 | try {
|
---|
506 | R_.Info()["NX"] = (int_8)Nx_;
|
---|
507 | R_.Info()["NY"] = (int_8)Ny_;
|
---|
508 | R_.Info()["NZ"] = (int_8)Nz_;
|
---|
509 | R_.Info()["DX"] = (r_8)Dx_;
|
---|
510 | R_.Info()["DY"] = (r_8)Dy_;
|
---|
511 | R_.Info()["DZ"] = (r_8)Dz_;
|
---|
512 | R_.Info()["ZREFPK"] = (r_8)compute_pk_redsh_ref_;
|
---|
513 | R_.Info()["ZREF"] = (r_8)redsh_ref_;
|
---|
514 | R_.Info()["KZREF"] = (r_8)kredsh_ref_;
|
---|
515 | R_.Info()["HREF"] = (r_8)h_ref_;
|
---|
516 | R_.Info()["Growth"] = (r_8)Dref();
|
---|
517 | R_.Info()["GrowthPk"] = (r_8)DrefPk();
|
---|
518 | POutPersist pos(cfname.c_str());
|
---|
519 | if(write_real) pos << PPFNameTag("rgen") << R_;
|
---|
520 | else pos << PPFNameTag("pkgen") << T_;
|
---|
521 | } catch (PThrowable & exc) {
|
---|
522 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
523 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
524 | return;
|
---|
525 | } catch (...) {
|
---|
526 | cout<<" some other exception was caught !"<<endl;
|
---|
527 | return;
|
---|
528 | }
|
---|
529 | }
|
---|
530 |
|
---|
531 | void GeneFluct3D::ReadPPF(string cfname)
|
---|
532 | {
|
---|
533 | cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
|
---|
534 | try {
|
---|
535 | bool from_real = true;
|
---|
536 | PInPersist pis(cfname.c_str());
|
---|
537 | string name_tag_k = "pkgen";
|
---|
538 | bool found_tag_k = pis.GotoNameTag("pkgen");
|
---|
539 | if(found_tag_k) {
|
---|
540 | cout<<" ...reading spectrum into TArray<complex <r_8> >"<<endl;
|
---|
541 | pis >> PPFNameTag("pkgen") >> T_;
|
---|
542 | from_real = false;
|
---|
543 | } else {
|
---|
544 | cout<<" ...reading space into TArray<r_8>"<<endl;
|
---|
545 | pis >> PPFNameTag("rgen") >> R_;
|
---|
546 | }
|
---|
547 | setpointers(from_real); // a mettre ici pour relire les DVInfo
|
---|
548 | int_8 nx = R_.Info()["NX"];
|
---|
549 | int_8 ny = R_.Info()["NY"];
|
---|
550 | int_8 nz = R_.Info()["NZ"];
|
---|
551 | r_8 dx = R_.Info()["DX"];
|
---|
552 | r_8 dy = R_.Info()["DY"];
|
---|
553 | r_8 dz = R_.Info()["DZ"];
|
---|
554 | r_8 pkzref = R_.Info()["ZREFPK"];
|
---|
555 | r_8 zref = R_.Info()["ZREF"];
|
---|
556 | r_8 kzref = R_.Info()["KZREF"];
|
---|
557 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
558 | init_fftw();
|
---|
559 | SetObservator(zref,kzref);
|
---|
560 | array_allocated_ = true;
|
---|
561 | compute_pk_redsh_ref_ = pkzref;
|
---|
562 | } catch (PThrowable & exc) {
|
---|
563 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
564 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
565 | return;
|
---|
566 | } catch (...) {
|
---|
567 | cout<<" some other exception was caught !"<<endl;
|
---|
568 | return;
|
---|
569 | }
|
---|
570 | }
|
---|
571 |
|
---|
572 | void GeneFluct3D::WriteSlicePPF(string cfname)
|
---|
573 | // On ecrit 3 tranches du cube selon chaque axe
|
---|
574 | {
|
---|
575 | cout<<"--- GeneFluct3D::WriteSlicePPF: Writing Cube Slices "<<cfname<<endl;
|
---|
576 | try {
|
---|
577 |
|
---|
578 | POutPersist pos(cfname.c_str());
|
---|
579 | TMatrix<r_4> S;
|
---|
580 | char str[16];
|
---|
581 | long i,j,l;
|
---|
582 |
|
---|
583 | // Tranches en Z
|
---|
584 | for(int s=0;s<3;s++) {
|
---|
585 | S.ReSize(Nx_,Ny_);
|
---|
586 | if(s==0) l=0; else if(s==1) l=(Nz_+1)/2; else l=Nz_-1;
|
---|
587 | sprintf(str,"z%ld",l);
|
---|
588 | for(i=0;i<Nx_;i++) for(j=0;j<Ny_;j++) S(i,j)=data_[IndexR(i,j,l)];
|
---|
589 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
590 | }
|
---|
591 |
|
---|
592 | // Tranches en Y
|
---|
593 | for(int s=0;s<3;s++) {
|
---|
594 | S.ReSize(Nz_,Nx_);
|
---|
595 | if(s==0) j=0; else if(s==1) j=(Ny_+1)/2; else j=Ny_-1;
|
---|
596 | sprintf(str,"y%ld",j);
|
---|
597 | for(i=0;i<Nx_;i++) for(l=0;l<Nz_;l++) S(l,i)=data_[IndexR(i,j,l)];
|
---|
598 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
599 | }
|
---|
600 |
|
---|
601 | // Tranches en X
|
---|
602 | for(int s=0;s<3;s++) {
|
---|
603 | S.ReSize(Nz_,Ny_);
|
---|
604 | if(s==0) i=0; else if(s==1) i=(Nx_+1)/2; else i=Nx_-1;
|
---|
605 | sprintf(str,"x%ld",i);
|
---|
606 | for(j=0;j<Ny_;j++) for(l=0;l<Nz_;l++) S(l,j)=data_[IndexR(i,j,l)];
|
---|
607 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
608 | }
|
---|
609 |
|
---|
610 | } catch (PThrowable & exc) {
|
---|
611 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
612 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
613 | return;
|
---|
614 | } catch (...) {
|
---|
615 | cout<<" some other exception was caught !"<<endl;
|
---|
616 | return;
|
---|
617 | }
|
---|
618 | }
|
---|
619 |
|
---|
620 | //-------------------------------------------------------
|
---|
621 | void GeneFluct3D::NTupleCheck(POutPersist &pos,string ntname,unsigned long nent)
|
---|
622 | // Remplit le NTuple "ntname" avec "nent" valeurs du cube (reel ou complex) et l'ecrit dans "pos"
|
---|
623 | {
|
---|
624 | if(ntname.size()<=0 || nent==0) return;
|
---|
625 | int nvar = 0;
|
---|
626 | if(array_type==1) nvar = 3;
|
---|
627 | else if(array_type==2) nvar = 4;
|
---|
628 | else return;
|
---|
629 | const char *vname[4] = {"t","z","re","im"};
|
---|
630 | float xnt[4];
|
---|
631 | NTuple nt(nvar,vname);
|
---|
632 |
|
---|
633 | if(array_type==1) {
|
---|
634 | unsigned long nmod = Nx_*Ny_*Nz_/nent; if(nmod==0) nmod=1;
|
---|
635 | unsigned long n=0;
|
---|
636 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
637 | if(n==nmod) {
|
---|
638 | int_8 ip = IndexR(i,j,l);
|
---|
639 | xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=data_[ip];
|
---|
640 | nt.Fill(xnt);
|
---|
641 | n=0;
|
---|
642 | }
|
---|
643 | n++;
|
---|
644 | }
|
---|
645 | } else {
|
---|
646 | unsigned long nmod = Nx_*Ny_*NCz_/nent; if(nmod==0) nmod=1;
|
---|
647 | unsigned long n=0;
|
---|
648 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<NCz_;l++) {
|
---|
649 | if(n==nmod) {
|
---|
650 | xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=T_(l,j,i).real(); xnt[3]=T_(l,j,i).imag();
|
---|
651 | nt.Fill(xnt);
|
---|
652 | n=0;
|
---|
653 | }
|
---|
654 | n++;
|
---|
655 | }
|
---|
656 | }
|
---|
657 |
|
---|
658 | pos.PutObject(nt,ntname);
|
---|
659 | }
|
---|
660 |
|
---|
661 | //-------------------------------------------------------
|
---|
662 | void GeneFluct3D::Print(void)
|
---|
663 | {
|
---|
664 | cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
|
---|
665 | cout<<"Space Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<Nz_<<" ("<<NTz_<<") size="
|
---|
666 | <<NRtot_<<endl;
|
---|
667 | cout<<" Resol: dx="<<Dx_<<" dy="<<Dy_<<" dz="<<Dz_<<" Mpc"
|
---|
668 | <<", dVol="<<dVol_<<", Vol="<<Vol_<<" Mpc^3"<<endl;
|
---|
669 | cout<<"Fourier Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<NCz_<<" ncoeff="<<NFcoef_<<endl;
|
---|
670 | cout<<" Resol: dkx="<<Dkx_<<" dky="<<Dky_<<" dkz="<<Dkz_<<" Mpc^-1"
|
---|
671 | <<", Dk3="<<Dk3_<<" Mpc^-3"<<endl;
|
---|
672 | cout<<" (2Pi/k: "<<2.*M_PI/Dkx_<<" "<<2.*M_PI/Dky_<<" "<<2.*M_PI/Dkz_<<" Mpc)"<<endl;
|
---|
673 | cout<<" Nyquist: kx="<<Knyqx_<<" ky="<<Knyqy_<<" kz="<<Knyqz_<<" Mpc^-1"
|
---|
674 | <<", Kmax="<<GetKmax()<<" Mpc^-1"<<endl;
|
---|
675 | cout<<" (2Pi/k: "<<2.*M_PI/Knyqx_<<" "<<2.*M_PI/Knyqy_<<" "<<2.*M_PI/Knyqz_<<" Mpc)"<<endl;
|
---|
676 | cout<<"Redshift "<<redsh_ref_<<" for z axe at k="<<kredsh_ref_<<endl;
|
---|
677 | }
|
---|
678 |
|
---|
679 | //-------------------------------------------------------
|
---|
680 | void GeneFluct3D::ComputeFourier0(PkSpectrumZ& pk_at_z)
|
---|
681 | // cf ComputeFourier() mais avec autre methode de realisation du spectre
|
---|
682 | // (attention on fait une fft pour realiser le spectre)
|
---|
683 | {
|
---|
684 | compute_pk_redsh_ref_ = pk_at_z.GetZ();
|
---|
685 | // --- realisation d'un tableau de tirage gaussiens
|
---|
686 | if(lp_>0) cout<<"--- ComputeFourier0 at z="<<compute_pk_redsh_ref_
|
---|
687 | <<": before gaussian filling ---"<<endl;
|
---|
688 | // On tient compte du pb de normalisation de FFTW3
|
---|
689 | double sntot = sqrt((double)NRtot_);
|
---|
690 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
691 | int_8 ip = IndexR(i,j,l);
|
---|
692 | data_[ip] = NorRand()/sntot;
|
---|
693 | }
|
---|
694 |
|
---|
695 | // --- realisation d'un tableau de tirage gaussiens
|
---|
696 | if(lp_>0) cout<<"...before fft real ---"<<endl;
|
---|
697 | GEN3D_FFTW_EXECUTE(pf_);
|
---|
698 |
|
---|
699 | // --- On remplit avec une realisation
|
---|
700 | if(lp_>0) cout<<"...before Fourier realization filling"<<endl;
|
---|
701 | T_(0,0,0) = complex<GEN3D_TYPE>(0.); // on coupe le continue et on l'initialise
|
---|
702 | long lmod = Nx_/20; if(lmod<1) lmod=1;
|
---|
703 | for(long i=0;i<Nx_;i++) {
|
---|
704 | double kx = AbsKx(i); kx *= kx;
|
---|
705 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<"\t nx-i="<<Nx_-i<<"\t kx="<<Kx(i)<<"\t pk="<<pk_at_z(kx)<<endl;
|
---|
706 | for(long j=0;j<Ny_;j++) {
|
---|
707 | double ky = AbsKy(j); ky *= ky;
|
---|
708 | for(long l=0;l<NCz_;l++) {
|
---|
709 | double kz = AbsKz(l); kz *= kz;
|
---|
710 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
711 | double k = sqrt(kx+ky+kz);
|
---|
712 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
713 | double pk = pk_at_z(k)/Vol_;
|
---|
714 | // ici pas de "/2" a cause de la remarque ci-dessus
|
---|
715 | T_(l,j,i) *= sqrt(pk);
|
---|
716 | }
|
---|
717 | }
|
---|
718 | }
|
---|
719 |
|
---|
720 | array_type = 2;
|
---|
721 |
|
---|
722 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
723 | double p = compute_power_carte();
|
---|
724 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
725 |
|
---|
726 | }
|
---|
727 |
|
---|
728 | //-------------------------------------------------------
|
---|
729 | void GeneFluct3D::ComputeFourier(PkSpectrumZ& pk_at_z)
|
---|
730 | // Calcule une realisation du spectre "pk_at_z"
|
---|
731 | // Attention: dans TArray le premier indice varie le + vite
|
---|
732 | // Explication normalisation: see Coles & Lucchin, Cosmology, p264-265
|
---|
733 | // FFTW3: on note N=Nx*Ny*Nz
|
---|
734 | // f --(FFT)--> F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f))
|
---|
735 | // sum(f(x_i)^2) = S
|
---|
736 | // sum(F(nu_i)^2) = S*N
|
---|
737 | // sum(fb(x_i)^2) = S*N^2
|
---|
738 | {
|
---|
739 | // --- RaZ du tableau
|
---|
740 | T_ = complex<GEN3D_TYPE>(0.);
|
---|
741 | compute_pk_redsh_ref_ = pk_at_z.GetZ();
|
---|
742 |
|
---|
743 | // --- On remplit avec une realisation
|
---|
744 | if(lp_>0) cout<<"--- ComputeFourier at z="<<compute_pk_redsh_ref_<<" ---"<<endl;
|
---|
745 | long lmod = Nx_/20; if(lmod<1) lmod=1;
|
---|
746 | for(long i=0;i<Nx_;i++) {
|
---|
747 | double kx = AbsKx(i); kx *= kx;
|
---|
748 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<"\t nx-i="<<Nx_-i<<"\t kx="<<Kx(i)<<"\t pk="<<pk_at_z(kx)<<endl;
|
---|
749 | for(long j=0;j<Ny_;j++) {
|
---|
750 | double ky = AbsKy(j); ky *= ky;
|
---|
751 | for(long l=0;l<NCz_;l++) {
|
---|
752 | double kz = AbsKz(l); kz *= kz;
|
---|
753 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
754 | double k = sqrt(kx+ky+kz);
|
---|
755 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
756 | double pk = pk_at_z(k)/Vol_;
|
---|
757 | // Explication de la division par 2: voir perandom.cc
|
---|
758 | // ou egalement Coles & Lucchin, Cosmology formula 13.7.2 p279
|
---|
759 | T_(l,j,i) = ComplexGaussianRand(sqrt(pk/2.));
|
---|
760 | }
|
---|
761 | }
|
---|
762 | }
|
---|
763 |
|
---|
764 | array_type = 2;
|
---|
765 | manage_coefficients(); // gros effet pour les spectres que l'on utilise !
|
---|
766 |
|
---|
767 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
768 | double p = compute_power_carte();
|
---|
769 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
770 |
|
---|
771 | }
|
---|
772 |
|
---|
773 | long GeneFluct3D::manage_coefficients(void)
|
---|
774 | // Take into account the real and complexe conjugate coefficients
|
---|
775 | // because we want a realization of a real data in real space
|
---|
776 | // On ecrit que: conj(P(k_x,k_y,k_z)) = P(-k_x,-k_y,-k_z)
|
---|
777 | // avec k_x = i, -k_x = N_x - i etc...
|
---|
778 | {
|
---|
779 | if(lp_>1) cout<<"...managing coefficients"<<endl;
|
---|
780 | check_array_alloc();
|
---|
781 |
|
---|
782 | // 1./ Le Continu et Nyquist sont reels
|
---|
783 | int_8 nreal = 0;
|
---|
784 | for(long kk=0;kk<2;kk++) {
|
---|
785 | long k=0; // continu
|
---|
786 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
787 | for(long jj=0;jj<2;jj++) {
|
---|
788 | long j=0;
|
---|
789 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
790 | for(long ii=0;ii<2;ii++) {
|
---|
791 | long i=0;
|
---|
792 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
793 | int_8 ip = IndexC(i,j,k);
|
---|
794 | //cout<<"i="<<i<<" j="<<j<<" k="<<k<<" = ("<<fdata_[ip][0]<<","<<fdata_[ip][1]<<")"<<endl;
|
---|
795 | fdata_[ip][1] = 0.; fdata_[ip][0] *= M_SQRT2;
|
---|
796 | nreal++;
|
---|
797 | }
|
---|
798 | }
|
---|
799 | }
|
---|
800 | if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
|
---|
801 |
|
---|
802 | // 2./ Les elements complexe conjugues (tous dans le plan k=0,Nyquist)
|
---|
803 |
|
---|
804 | // a./ les lignes et colonnes du continu et de nyquist
|
---|
805 | int_8 nconj1 = 0;
|
---|
806 | for(long kk=0;kk<2;kk++) {
|
---|
807 | long k=0; // continu
|
---|
808 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
809 | for(long jj=0;jj<2;jj++) { // selon j
|
---|
810 | long j=0;
|
---|
811 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
812 | for(long i=1;i<(Nx_+1)/2;i++) {
|
---|
813 | int_8 ip = IndexC(i,j,k);
|
---|
814 | int_8 ip1 = IndexC(Nx_-i,j,k);
|
---|
815 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
816 | nconj1++;
|
---|
817 | }
|
---|
818 | }
|
---|
819 | for(long ii=0;ii<2;ii++) {
|
---|
820 | long i=0;
|
---|
821 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
822 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
823 | int_8 ip = IndexC(i,j,k);
|
---|
824 | int_8 ip1 = IndexC(i,Ny_-j,k);
|
---|
825 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
826 | nconj1++;
|
---|
827 | }
|
---|
828 | }
|
---|
829 | }
|
---|
830 | if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
|
---|
831 |
|
---|
832 | // b./ les lignes et colonnes hors continu et de nyquist
|
---|
833 | int_8 nconj2 = 0;
|
---|
834 | for(long kk=0;kk<2;kk++) {
|
---|
835 | long k=0; // continu
|
---|
836 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
837 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
838 | if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
|
---|
839 | for(long i=1;i<Nx_;i++) {
|
---|
840 | if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
|
---|
841 | int_8 ip = IndexC(i,j,k);
|
---|
842 | int_8 ip1 = IndexC(Nx_-i,Ny_-j,k);
|
---|
843 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
844 | nconj2++;
|
---|
845 | }
|
---|
846 | }
|
---|
847 | }
|
---|
848 | if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
|
---|
849 |
|
---|
850 | if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(NFcoef_-nconj1-nconj2)-nreal<<endl;
|
---|
851 |
|
---|
852 | return nreal+nconj1+nconj2;
|
---|
853 | }
|
---|
854 |
|
---|
855 | double GeneFluct3D::compute_power_carte(void)
|
---|
856 | // Calcul la puissance de la realisation du spectre Pk
|
---|
857 | {
|
---|
858 | check_array_alloc();
|
---|
859 |
|
---|
860 | double s2 = 0.;
|
---|
861 | for(long l=0;l<NCz_;l++)
|
---|
862 | for(long j=0;j<Ny_;j++)
|
---|
863 | for(long i=0;i<Nx_;i++) s2 += MODULE2(T_(l,j,i));
|
---|
864 |
|
---|
865 | double s20 = 0.;
|
---|
866 | for(long j=0;j<Ny_;j++)
|
---|
867 | for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
|
---|
868 |
|
---|
869 | double s2n = 0.;
|
---|
870 | if(Nz_%2==0)
|
---|
871 | for(long j=0;j<Ny_;j++)
|
---|
872 | for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
|
---|
873 |
|
---|
874 | return 2.*s2 -s20 -s2n;
|
---|
875 | }
|
---|
876 |
|
---|
877 | //-------------------------------------------------------------------
|
---|
878 | void GeneFluct3D::FilterByPixel(void)
|
---|
879 | // Filtrage par la fonction fenetre du pixel (parallelepipede)
|
---|
880 | // TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
|
---|
881 | // e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
|
---|
882 | // = 2/(k_x*dx) * sin(k_x*dx/2) * (idem y) * (idem z)
|
---|
883 | // Gestion divergence en 0: sin(y)/y = 1 - y^2/6*(1-y^2/20)
|
---|
884 | // avec y = k_x*dx/2
|
---|
885 | {
|
---|
886 | if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
|
---|
887 | check_array_alloc();
|
---|
888 |
|
---|
889 | for(long i=0;i<Nx_;i++) {
|
---|
890 | double kx = AbsKx(i) *Dx_/2;
|
---|
891 | double pk_x = pixelfilter(kx);
|
---|
892 | for(long j=0;j<Ny_;j++) {
|
---|
893 | double ky = AbsKy(j) *Dy_/2;
|
---|
894 | double pk_y = pixelfilter(ky);
|
---|
895 | for(long l=0;l<NCz_;l++) {
|
---|
896 | double kz = AbsKz(l) *Dz_/2;
|
---|
897 | double pk_z = pixelfilter(kz);
|
---|
898 | T_(l,j,i) *= pk_x*pk_y*pk_z;
|
---|
899 | }
|
---|
900 | }
|
---|
901 | }
|
---|
902 |
|
---|
903 | }
|
---|
904 |
|
---|
905 | //-------------------------------------------------------------------
|
---|
906 | void GeneFluct3D::ToRedshiftSpace(void)
|
---|
907 | // Apply redshift distortion corrections
|
---|
908 | // ---- ATTENTION: Le spectre Pk doit etre (dRho/Rho)(k)
|
---|
909 | {
|
---|
910 | double zpk = compute_pk_redsh_ref_;
|
---|
911 | double beta = -(1.+zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
912 | if(lp_>0) cout<<"--- ToRedshiftSpace --- at z="<<zpk<<", Beta="<<beta
|
---|
913 | <<" (km/s)/Mpc (comp. for dz="<<good_dzinc_<<")"<<endl;
|
---|
914 | check_array_alloc();
|
---|
915 |
|
---|
916 | for(long i=0;i<Nx_;i++) {
|
---|
917 | double kx = Kx(i);
|
---|
918 | for(long j=0;j<Ny_;j++) {
|
---|
919 | double ky = Ky(j);
|
---|
920 | double kt2 = kx*kx + ky*ky;
|
---|
921 | for(long l=0;l<NCz_;l++) {
|
---|
922 | double kz = Kz(l);
|
---|
923 | if(l==0) continue;
|
---|
924 | T_(l,j,i) *= (1. + beta*kz*kz/(kt2 + kz*kz));
|
---|
925 | }
|
---|
926 | }
|
---|
927 | }
|
---|
928 | }
|
---|
929 |
|
---|
930 | //-------------------------------------------------------------------
|
---|
931 | void GeneFluct3D::ToVelLoS(void)
|
---|
932 | /*
|
---|
933 | ---- Vitesse particuliere
|
---|
934 | Un atome au redshift "z" emet a la longueur d'onde "Le" si il est au repos
|
---|
935 | Il a une vitesse (particuliere) "V" dans le referentiel comobile au redshift "z"
|
---|
936 | Dans le referentiel comobile au redshift "z", il emet a la longueur d'onde "Le*(1+V/C)"
|
---|
937 | L'observateur voit donc une longueur d'onde "Le*(1+V/C)*(1+z) = Le*[1+z+V/C*(1+z)]"
|
---|
938 | et croit que l'atome est au redshift "z' = z + V/C*(1+z)"
|
---|
939 | L'observateur observe donc une vitesse particuliere (1+z)*V/c
|
---|
940 | --> vitesse particuliere dans le repere comobile en z: V/C
|
---|
941 | vitesse particuliere dans le repere comobile de l'observateur: (1+z)*V/C
|
---|
942 | ---- ATTENTION: Le spectre Pk doit etre (dRho/Rho)(k)
|
---|
943 | ---- On genere la vitesse particuliere (dans le referentiel comobile au redshift "z")
|
---|
944 | Vlos(k) = i*Beta(z)*k(los)/k^2*(dRho/Rho)(k)
|
---|
945 | .avec Beta(z) = dln(D)/da = -(1+z)/D*dD/dz
|
---|
946 | .on convertit en vitesse en faisant Vlos(k)/H(z)
|
---|
947 | */
|
---|
948 | {
|
---|
949 | double zpk = compute_pk_redsh_ref_;
|
---|
950 | double dpsd = -cosmo_->H(zpk) * (1.+zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
951 | if(lp_>0) cout<<"--- ToVelLoS --- at z="<<zpk<<", -H*(1+z)*D'/D="<<dpsd
|
---|
952 | <<" (km/s)/Mpc (comp. for dz="<<good_dzinc_<<")"<<endl;
|
---|
953 | check_array_alloc();
|
---|
954 |
|
---|
955 | for(long i=0;i<Nx_;i++) {
|
---|
956 | double kx = Kx(i);
|
---|
957 | for(long j=0;j<Ny_;j++) {
|
---|
958 | double ky = Ky(j);
|
---|
959 | double kt2 = kx*kx + ky*ky;
|
---|
960 | for(long l=0;l<NCz_;l++) {
|
---|
961 | double kz = Kz(l);
|
---|
962 | double k2 = kt2 + kz*kz;
|
---|
963 | if(l==0) k2 = 0.; else k2 = dpsd*kz/k2;
|
---|
964 | T_(l,j,i) *= complex<double>(0.,k2);
|
---|
965 | }
|
---|
966 | }
|
---|
967 | }
|
---|
968 |
|
---|
969 | }
|
---|
970 |
|
---|
971 | //-------------------------------------------------------------------
|
---|
972 | void GeneFluct3D::ApplyGrowthFactor(int type_evol)
|
---|
973 | // Apply Growth to real space d(Rho)/Rho cube
|
---|
974 | // Using the correspondance between redshift and los comoving distance
|
---|
975 | // describe in vector "zred_" "loscom_"
|
---|
976 | // type_evol = 1 : evolution avec la distance a l'observateur
|
---|
977 | // 2 : evolution avec la distance du plan Z
|
---|
978 | // (tous les pixels d'un plan Z sont mis au meme redshift z que celui du milieu)
|
---|
979 | {
|
---|
980 | if(lp_>0) cout<<"--- ApplyGrowthFactor: evol="<<type_evol<<endl;
|
---|
981 | check_array_alloc();
|
---|
982 |
|
---|
983 | if(growth_ == NULL) {
|
---|
984 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
|
---|
985 | cout<<bla<<endl; throw ParmError(bla);
|
---|
986 | }
|
---|
987 | if(type_evol<1 || type_evol>2) {
|
---|
988 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
|
---|
989 | cout<<bla<<endl; throw ParmError(bla);
|
---|
990 | }
|
---|
991 |
|
---|
992 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
993 |
|
---|
994 | for(long i=0;i<Nx_;i++) {
|
---|
995 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
996 | for(long j=0;j<Ny_;j++) {
|
---|
997 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
998 | for(long l=0;l<Nz_;l++) {
|
---|
999 | double dz = DZcom(l);
|
---|
1000 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
1001 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
1002 | double z = interpinv(dz); // interpolation par morceau
|
---|
1003 | double dzgr = (*growth_)(z);
|
---|
1004 | int_8 ip = IndexR(i,j,l);
|
---|
1005 | data_[ip] *= dzgr;
|
---|
1006 | }
|
---|
1007 | }
|
---|
1008 | }
|
---|
1009 |
|
---|
1010 | }
|
---|
1011 |
|
---|
1012 | //-------------------------------------------------------------------
|
---|
1013 | void GeneFluct3D::ApplyRedshiftSpaceGrowthFactor(void)
|
---|
1014 | // Apply Growth(z) to redshift-distorted real space cube
|
---|
1015 | {
|
---|
1016 | const char *bla = "GeneFluct3D::ApplyRedshiftSpaceGrowthFactor_Error: redshift evolution impossible because factor depen on k";
|
---|
1017 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1018 | }
|
---|
1019 |
|
---|
1020 | //-------------------------------------------------------------------
|
---|
1021 | void GeneFluct3D::ApplyVelLosGrowthFactor(int type_evol)
|
---|
1022 | // Apply Growth(z) to transverse velocity real space cube
|
---|
1023 | {
|
---|
1024 | if(lp_>0) cout<<"--- ApplyVelLosGrowthFactor: evol="<<type_evol<<endl;
|
---|
1025 | check_array_alloc();
|
---|
1026 |
|
---|
1027 | if(growth_ == NULL) {
|
---|
1028 | const char *bla = "GeneFluct3D::ApplyVelLosGrowthFactor_Error: set GrowthFactor first";
|
---|
1029 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1030 | }
|
---|
1031 | if(type_evol<1 || type_evol>2) {
|
---|
1032 | const char *bla = "GeneFluct3D::ApplyVelLosGrowthFactor_Error: bad type_evol value";
|
---|
1033 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1034 | }
|
---|
1035 |
|
---|
1036 | double zpk = compute_pk_redsh_ref_;
|
---|
1037 | double dpsd_orig = - cosmo_->H(zpk) * (1.+zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
1038 |
|
---|
1039 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
1040 | InterpFunc interdpd(zredmin_dpsd_,zredmax_dpsd_,dpsdfrzred_);
|
---|
1041 |
|
---|
1042 | for(long i=0;i<Nx_;i++) {
|
---|
1043 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
1044 | for(long j=0;j<Ny_;j++) {
|
---|
1045 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
1046 | for(long l=0;l<Nz_;l++) {
|
---|
1047 | double dz = DZcom(l);
|
---|
1048 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
1049 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
1050 | double z = interpinv(dz); // interpolation par morceau
|
---|
1051 | double dpsd = interdpd(z);
|
---|
1052 | int_8 ip = IndexR(i,j,l);
|
---|
1053 | data_[ip] *= dpsd / dpsd_orig;
|
---|
1054 | }
|
---|
1055 | }
|
---|
1056 | }
|
---|
1057 |
|
---|
1058 | }
|
---|
1059 |
|
---|
1060 | //-------------------------------------------------------------------
|
---|
1061 | void GeneFluct3D::ComputeReal(void)
|
---|
1062 | // Calcule une realisation dans l'espace reel
|
---|
1063 | {
|
---|
1064 | if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
|
---|
1065 | check_array_alloc();
|
---|
1066 |
|
---|
1067 | // On fait la FFT
|
---|
1068 | GEN3D_FFTW_EXECUTE(pb_);
|
---|
1069 | array_type = 1;
|
---|
1070 | }
|
---|
1071 |
|
---|
1072 | //-------------------------------------------------------------------
|
---|
1073 | void GeneFluct3D::ReComputeFourier(void)
|
---|
1074 | {
|
---|
1075 | if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
|
---|
1076 | check_array_alloc();
|
---|
1077 |
|
---|
1078 | // On fait la FFT
|
---|
1079 | GEN3D_FFTW_EXECUTE(pf_);
|
---|
1080 | array_type = 2;
|
---|
1081 |
|
---|
1082 | // On corrige du pb de la normalisation de FFTW3
|
---|
1083 | complex<r_8> v((r_8)NRtot_,0.);
|
---|
1084 | for(long i=0;i<Nx_;i++)
|
---|
1085 | for(long j=0;j<Ny_;j++)
|
---|
1086 | for(long l=0;l<NCz_;l++) T_(l,j,i) /= v;
|
---|
1087 | }
|
---|
1088 |
|
---|
1089 | //-------------------------------------------------------------------
|
---|
1090 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
|
---|
1091 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
1092 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
1093 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
1094 | {
|
---|
1095 | if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
|
---|
1096 | check_array_alloc();
|
---|
1097 |
|
---|
1098 | if(herr.NBins()<0) return -1;
|
---|
1099 | herr.Zero();
|
---|
1100 |
|
---|
1101 | // Attention a l'ordre
|
---|
1102 | for(long i=0;i<Nx_;i++) {
|
---|
1103 | double kx = AbsKx(i); kx *= kx;
|
---|
1104 | for(long j=0;j<Ny_;j++) {
|
---|
1105 | double ky = AbsKy(j); ky *= ky;
|
---|
1106 | for(long l=0;l<NCz_;l++) {
|
---|
1107 | double kz = AbsKz(l);
|
---|
1108 | double k = sqrt(kx+ky+kz*kz);
|
---|
1109 | double pk = MODULE2(T_(l,j,i));
|
---|
1110 | herr.Add(k,pk);
|
---|
1111 | }
|
---|
1112 | }
|
---|
1113 | }
|
---|
1114 | herr.ToVariance();
|
---|
1115 |
|
---|
1116 | // renormalize to directly compare to original spectrum
|
---|
1117 | double norm = Vol_;
|
---|
1118 | herr *= norm;
|
---|
1119 |
|
---|
1120 | return 0;
|
---|
1121 | }
|
---|
1122 |
|
---|
1123 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
|
---|
1124 | {
|
---|
1125 | if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
|
---|
1126 | check_array_alloc();
|
---|
1127 |
|
---|
1128 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
1129 | herr.Zero();
|
---|
1130 |
|
---|
1131 | // Attention a l'ordre
|
---|
1132 | for(long i=0;i<Nx_;i++) {
|
---|
1133 | double kx = AbsKx(i); kx *= kx;
|
---|
1134 | for(long j=0;j<Ny_;j++) {
|
---|
1135 | double ky = AbsKy(j); ky *= ky;
|
---|
1136 | double kt = sqrt(kx+ky);
|
---|
1137 | for(long l=0;l<NCz_;l++) {
|
---|
1138 | double kz = AbsKz(l);
|
---|
1139 | double pk = MODULE2(T_(l,j,i));
|
---|
1140 | herr.Add(kt,kz,pk);
|
---|
1141 | }
|
---|
1142 | }
|
---|
1143 | }
|
---|
1144 | herr.ToVariance();
|
---|
1145 |
|
---|
1146 | // renormalize to directly compare to original spectrum
|
---|
1147 | double norm = Vol_;
|
---|
1148 | herr *= norm;
|
---|
1149 |
|
---|
1150 | return 0;
|
---|
1151 | }
|
---|
1152 |
|
---|
1153 | //-------------------------------------------------------------------
|
---|
1154 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor)
|
---|
1155 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
1156 | // AVEC la soustraction du niveau de bruit et la correction par filterpixel()
|
---|
1157 | // Si on ne fait pas ca, alors on obtient un spectre non-isotrope!
|
---|
1158 | //
|
---|
1159 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
1160 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
1161 | {
|
---|
1162 | if(lp_>0) cout<<"--- ComputeSpectrum: sigma="<<sigma<<endl;
|
---|
1163 | check_array_alloc();
|
---|
1164 |
|
---|
1165 | if(sigma<=0.) sigma = 0.;
|
---|
1166 | double sigma2 = sigma*sigma / (double)NRtot_;
|
---|
1167 |
|
---|
1168 | if(herr.NBins()<0) return -1;
|
---|
1169 | herr.Zero();
|
---|
1170 |
|
---|
1171 | TVector<r_8> vfz(NCz_);
|
---|
1172 | if(pixcor) // kz = l*Dkz_
|
---|
1173 | for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
|
---|
1174 |
|
---|
1175 | // Attention a l'ordre
|
---|
1176 | for(long i=0;i<Nx_;i++) {
|
---|
1177 | double kx = AbsKx(i);
|
---|
1178 | double fx = (pixcor) ? pixelfilter(kx*Dx_/2): 1.;
|
---|
1179 | kx *= kx; fx *= fx;
|
---|
1180 | for(long j=0;j<Ny_;j++) {
|
---|
1181 | double ky = AbsKy(j);
|
---|
1182 | double fy = (pixcor) ? pixelfilter(ky*Dy_/2): 1.;
|
---|
1183 | ky *= ky; fy *= fy;
|
---|
1184 | for(long l=0;l<NCz_;l++) {
|
---|
1185 | double kz = AbsKz(l);
|
---|
1186 | double k = sqrt(kx+ky+kz*kz);
|
---|
1187 | double pk = MODULE2(T_(l,j,i)) - sigma2;
|
---|
1188 | double fz = (pixcor) ? vfz(l): 1.;
|
---|
1189 | double f = fx*fy*fz;
|
---|
1190 | if(f>0.) herr.Add(k,pk/f);
|
---|
1191 | }
|
---|
1192 | }
|
---|
1193 | }
|
---|
1194 | herr.ToVariance();
|
---|
1195 | for(int i=0;i<herr.NBins();i++) herr(i) += sigma2;
|
---|
1196 |
|
---|
1197 | // renormalize to directly compare to original spectrum
|
---|
1198 | double norm = Vol_;
|
---|
1199 | herr *= norm;
|
---|
1200 |
|
---|
1201 | return 0;
|
---|
1202 | }
|
---|
1203 |
|
---|
1204 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor)
|
---|
1205 | // AVEC la soustraction du niveau de bruit et la correction par filterpixel()
|
---|
1206 | {
|
---|
1207 | if(lp_>0) cout<<"--- ComputeSpectrum2D: sigma="<<sigma<<endl;
|
---|
1208 | check_array_alloc();
|
---|
1209 |
|
---|
1210 | if(sigma<=0.) sigma = 0.;
|
---|
1211 | double sigma2 = sigma*sigma / (double)NRtot_;
|
---|
1212 |
|
---|
1213 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
1214 | herr.Zero();
|
---|
1215 |
|
---|
1216 | TVector<r_8> vfz(NCz_);
|
---|
1217 | if(pixcor) // kz = l*Dkz_
|
---|
1218 | for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
|
---|
1219 |
|
---|
1220 | // Attention a l'ordre
|
---|
1221 | for(long i=0;i<Nx_;i++) {
|
---|
1222 | double kx = AbsKx(i);
|
---|
1223 | double fx = (pixcor) ? pixelfilter(kx*Dx_/2) : 1.;
|
---|
1224 | kx *= kx; fx *= fx;
|
---|
1225 | for(long j=0;j<Ny_;j++) {
|
---|
1226 | double ky = AbsKy(j);
|
---|
1227 | double fy = (pixcor) ? pixelfilter(ky*Dy_/2) : 1.;
|
---|
1228 | ky *= ky; fy *= fy;
|
---|
1229 | double kt = sqrt(kx+ky);
|
---|
1230 | for(long l=0;l<NCz_;l++) {
|
---|
1231 | double kz = AbsKz(l);
|
---|
1232 | double pk = MODULE2(T_(l,j,i)) - sigma2;
|
---|
1233 | double fz = (pixcor) ? vfz(l): 1.;
|
---|
1234 | double f = fx*fy*fz;
|
---|
1235 | if(f>0.) herr.Add(kt,kz,pk/f);
|
---|
1236 | }
|
---|
1237 | }
|
---|
1238 | }
|
---|
1239 | herr.ToVariance();
|
---|
1240 | for(int i=0;i<herr.NBinX();i++)
|
---|
1241 | for(int j=0;j<herr.NBinY();j++) herr(i,j) += sigma2;
|
---|
1242 |
|
---|
1243 | // renormalize to directly compare to original spectrum
|
---|
1244 | double norm = Vol_;
|
---|
1245 | herr *= norm;
|
---|
1246 |
|
---|
1247 | return 0;
|
---|
1248 | }
|
---|
1249 |
|
---|
1250 | //-------------------------------------------------------
|
---|
1251 | int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
|
---|
1252 | // Recompute MASS variance in spherical top-hat (rayon=R)
|
---|
1253 | // Par definition: SigmaR^2 = <(M-<M>)^2>/<M>^2
|
---|
1254 | // ou M = masse dans sphere de rayon R
|
---|
1255 | // --- ATTENTION: la variance calculee a une tres grande dispersion
|
---|
1256 | // (surtout si le volume du cube est petit). Pour verifier
|
---|
1257 | // que le sigmaR calcule par cette methode est en accord avec
|
---|
1258 | // le sigmaR en input, il faut faire plusieurs simulations (~100)
|
---|
1259 | // et regarder la moyenne des sigmaR reconstruits
|
---|
1260 | {
|
---|
1261 | if(lp_>0) cout<<"--- VarianceFrReal R="<<R<<endl;
|
---|
1262 | check_array_alloc();
|
---|
1263 |
|
---|
1264 | long dnx = long(R/Dx_)+1; if(dnx<=0) dnx = 1;
|
---|
1265 | long dny = long(R/Dy_)+1; if(dny<=0) dny = 1;
|
---|
1266 | long dnz = long(R/Dz_)+1; if(dnz<=0) dnz = 1;
|
---|
1267 | if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
|
---|
1268 |
|
---|
1269 | double sum=0., sum2=0., sn=0., r2 = R*R;
|
---|
1270 | int_8 nsum=0;
|
---|
1271 |
|
---|
1272 | for(long i=dnx;i<Nx_-dnx;i+=2*dnx) {
|
---|
1273 | for(long j=dny;j<Ny_-dny;j+=2*dny) {
|
---|
1274 | for(long l=dnz;l<Nz_-dnz;l+=2*dnz) {
|
---|
1275 | double m=0.; int_8 n=0;
|
---|
1276 | for(long ii=i-dnx;ii<=i+dnx;ii++) {
|
---|
1277 | double x = (ii-i)*Dx_; x *= x;
|
---|
1278 | for(long jj=j-dny;jj<=j+dny;jj++) {
|
---|
1279 | double y = (jj-j)*Dy_; y *= y;
|
---|
1280 | for(long ll=l-dnz;ll<=l+dnz;ll++) {
|
---|
1281 | double z = (ll-l)*Dz_; z *= z;
|
---|
1282 | if(x+y+z>r2) continue;
|
---|
1283 | int_8 ip = IndexR(ii,jj,ll);
|
---|
1284 | m += 1.+data_[ip]; // 1+drho/rho
|
---|
1285 | n++;
|
---|
1286 | }
|
---|
1287 | }
|
---|
1288 | }
|
---|
1289 | if(n>0) {sum += m; sum2 += m*m; nsum++; sn += n;}
|
---|
1290 | //cout<<i<<","<<j<<","<<l<<" n="<<n<<" m="<<m<<" sum="<<sum<<" sum2="<<sum2<<endl;
|
---|
1291 | }
|
---|
1292 | }
|
---|
1293 | }
|
---|
1294 |
|
---|
1295 | if(nsum<=1) {var=0.; return nsum;}
|
---|
1296 | sum /= nsum;
|
---|
1297 | sum2 = sum2/nsum - sum*sum;
|
---|
1298 | sn /= nsum;
|
---|
1299 | if(lp_>0) cout<<"...<n>="<<sn<<", nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
|
---|
1300 | var = sum2/(sum*sum); // <dM^2>/<M>^2
|
---|
1301 | if(lp_>0) cout<<"...sigmaR^2 = <(M-<M>)^2>/<M>^2 = "<<var
|
---|
1302 | <<" -> sigmaR = "<<sqrt(var)<<endl;
|
---|
1303 |
|
---|
1304 | return nsum;
|
---|
1305 | }
|
---|
1306 |
|
---|
1307 | //-------------------------------------------------------
|
---|
1308 | int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
|
---|
1309 | // number of pixels outside of ]vmin,vmax[ extremites exclues
|
---|
1310 | // -> vmin and vmax are considered as bad
|
---|
1311 | {
|
---|
1312 | check_array_alloc();
|
---|
1313 |
|
---|
1314 | int_8 nbad = 0;
|
---|
1315 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1316 | int_8 ip = IndexR(i,j,l);
|
---|
1317 | double v = data_[ip];
|
---|
1318 | if(v<=vmin || v>=vmax) nbad++;
|
---|
1319 | }
|
---|
1320 |
|
---|
1321 | if(lp_>0) cout<<"--- NumberOfBad "<<nbad<<" px out of ]"<<vmin<<","<<vmax
|
---|
1322 | <<"[ i.e. frac="<<nbad/(double)NRtot_<<endl;
|
---|
1323 | return nbad;
|
---|
1324 | }
|
---|
1325 |
|
---|
1326 | int_8 GeneFluct3D::MinMax(double& xmin,double& xmax,double vmin,double vmax)
|
---|
1327 | // Calcul des valeurs xmin et xmax dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
1328 | {
|
---|
1329 | bool tstval = (vmax>vmin)? true: false;
|
---|
1330 | if(lp_>0) {
|
---|
1331 | cout<<"--- MinMax";
|
---|
1332 | if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
|
---|
1333 | cout<<endl;
|
---|
1334 | }
|
---|
1335 | check_array_alloc();
|
---|
1336 |
|
---|
1337 | int_8 n = 0;
|
---|
1338 | xmin = xmax = data_[0];
|
---|
1339 |
|
---|
1340 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1341 | int_8 ip = IndexR(i,j,l);
|
---|
1342 | double x = data_[ip];
|
---|
1343 | if(tstval && (x<=vmin || x>=vmax)) continue;
|
---|
1344 | if(x<xmin) xmin = x;
|
---|
1345 | if(x>xmax) xmax = x;
|
---|
1346 | n++;
|
---|
1347 | }
|
---|
1348 |
|
---|
1349 | if(lp_>0) cout<<" n="<<n<<" min="<<xmin<<" max="<<xmax<<endl;
|
---|
1350 |
|
---|
1351 | return n;
|
---|
1352 | }
|
---|
1353 |
|
---|
1354 | int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax
|
---|
1355 | ,bool useout,double vout)
|
---|
1356 | // Calcul de mean,sigma2 dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
1357 | // useout = false: ne pas utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
1358 | // true : utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
1359 | // en remplacant leurs valeurs par "vout"
|
---|
1360 | {
|
---|
1361 | bool tstval = (vmax>vmin)? true: false;
|
---|
1362 | if(lp_>0) {
|
---|
1363 | cout<<"--- MeanSigma2";
|
---|
1364 | if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
|
---|
1365 | if(useout) cout<<", useout="<<useout<<" vout="<<vout;
|
---|
1366 | cout<<endl;
|
---|
1367 | }
|
---|
1368 | check_array_alloc();
|
---|
1369 |
|
---|
1370 | int_8 n = 0;
|
---|
1371 | rm = rs2 = 0.;
|
---|
1372 |
|
---|
1373 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1374 | int_8 ip = IndexR(i,j,l);
|
---|
1375 | double v = data_[ip];
|
---|
1376 | if(tstval) {
|
---|
1377 | if(v<=vmin || v>=vmax) {if(useout) v=vout; else continue;}
|
---|
1378 | }
|
---|
1379 | rm += v;
|
---|
1380 | rs2 += v*v;
|
---|
1381 | n++;
|
---|
1382 | }
|
---|
1383 |
|
---|
1384 | if(n>1) {
|
---|
1385 | rm /= (double)n;
|
---|
1386 | rs2 = rs2/(double)n - rm*rm;
|
---|
1387 | }
|
---|
1388 |
|
---|
1389 | if(lp_>0) cout<<" n="<<n<<" m="<<rm<<" s2="<<rs2<<" s="<<sqrt(fabs(rs2))<<endl;
|
---|
1390 |
|
---|
1391 | return n;
|
---|
1392 | }
|
---|
1393 |
|
---|
1394 | int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
|
---|
1395 | // set to "val0" if out of range ]vmin,vmax[ extremites exclues
|
---|
1396 | // cad set to "val0" if in [vmin,vmax] -> vmin and vmax are set to val0
|
---|
1397 | {
|
---|
1398 | check_array_alloc();
|
---|
1399 |
|
---|
1400 | int_8 nbad = 0;
|
---|
1401 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1402 | int_8 ip = IndexR(i,j,l);
|
---|
1403 | double v = data_[ip];
|
---|
1404 | if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
|
---|
1405 | }
|
---|
1406 |
|
---|
1407 | if(lp_>0) cout<<"--- SetToVal "<<nbad<<" px set to="<<val0
|
---|
1408 | <<" because out of range=]"<<vmin<<","<<vmax<<"["<<endl;
|
---|
1409 | return nbad;
|
---|
1410 | }
|
---|
1411 |
|
---|
1412 | void GeneFluct3D::ScaleOffset(double scalecube,double offsetcube)
|
---|
1413 | // Replace "V" by "scalecube * ( V + offsetcube )"
|
---|
1414 | {
|
---|
1415 | if(lp_>0) cout<<"--- ScaleCube scale="<<scalecube<<" offset="<<offsetcube<<endl;
|
---|
1416 |
|
---|
1417 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1418 | int_8 ip = IndexR(i,j,l);
|
---|
1419 | data_[ip] = scalecube * ( data_[ip] + offsetcube );
|
---|
1420 | }
|
---|
1421 |
|
---|
1422 | return;
|
---|
1423 | }
|
---|
1424 |
|
---|
1425 | //-------------------------------------------------------
|
---|
1426 | void GeneFluct3D::TurnFluct2Mass(void)
|
---|
1427 | // d_rho/rho -> Mass (add one!)
|
---|
1428 | {
|
---|
1429 | if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
|
---|
1430 | check_array_alloc();
|
---|
1431 |
|
---|
1432 |
|
---|
1433 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1434 | int_8 ip = IndexR(i,j,l);
|
---|
1435 | data_[ip] += 1.;
|
---|
1436 | }
|
---|
1437 | }
|
---|
1438 |
|
---|
1439 | double GeneFluct3D::TurnFluct2MeanNumber(double val_by_mpc3)
|
---|
1440 | // ATTENTION: la gestion des pixels<0 proposee ici induit une perte de variance
|
---|
1441 | // sur la carte, le spectre Pk reconstruit sera plus faible!
|
---|
1442 | // L'effet sera d'autant plus grand que le nombre de pixels<0 sera grand.
|
---|
1443 | {
|
---|
1444 | if(lp_>0) cout<<"--- TurnFluct2MeanNumber : "<<val_by_mpc3<<" quantity (gal or mass)/Mpc^3"<<endl;
|
---|
1445 |
|
---|
1446 | // First convert dRho/Rho into 1+dRho/Rho
|
---|
1447 | int_8 nball = 0; double sumall = 0., sumall2 = 0.;
|
---|
1448 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1449 | int_8 ip = IndexR(i,j,l);
|
---|
1450 | data_[ip] += 1.;
|
---|
1451 | nball++; sumall += data_[ip]; sumall2 += data_[ip]*data_[ip];
|
---|
1452 | }
|
---|
1453 | if(nball>2) {
|
---|
1454 | sumall /= (double)nball;
|
---|
1455 | sumall2 = sumall2/(double)nball - sumall*sumall;
|
---|
1456 | if(lp_>0) cout<<"1+dRho/Rho: mean="<<sumall<<" variance="<<sumall2
|
---|
1457 | <<" -> "<<sqrt(fabs(sumall2))<<endl;
|
---|
1458 | }
|
---|
1459 |
|
---|
1460 | // Find contribution for positive pixels
|
---|
1461 | int_8 nbpos = 0; double sumpos = 0. , sumpos2 = 0.;
|
---|
1462 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1463 | int_8 ip = IndexR(i,j,l);
|
---|
1464 | double v = data_[ip];
|
---|
1465 | if(data_[ip]>0.) {nbpos++; sumpos += v; sumpos2 += v*v;}
|
---|
1466 | }
|
---|
1467 | if(nbpos<1) {
|
---|
1468 | cout<<"TurnFluct2MeanNumber_Error: nbpos<1"<<endl;
|
---|
1469 | throw RangeCheckError("TurnFluct2MeanNumber_Error: nbpos<1");
|
---|
1470 | }
|
---|
1471 | sumpos2 = sumpos2/nball - sumpos*sumpos/(nball*nball);
|
---|
1472 | if(lp_>0)
|
---|
1473 | cout<<"1+dRho/Rho with v<0 set to zero: mean="<<sumpos/nball
|
---|
1474 | <<" variance="<<sumpos2<<" -> "<<sqrt(fabs(sumpos2))<<endl;
|
---|
1475 | cout<<"Sum of positive values: sumpos="<<sumpos
|
---|
1476 | <<" (n(v>0) = "<<nbpos<<" frac(v>0)="<<nbpos/(double)NRtot_<<")"<<endl;
|
---|
1477 |
|
---|
1478 | // - Mettre exactement val_by_mpc3*Vol galaxies (ou Msol) dans notre survey
|
---|
1479 | // - Uniquement dans les pixels de masse >0.
|
---|
1480 | // - Mise a zero des pixels <0
|
---|
1481 | double dn = val_by_mpc3 * Vol_ / sumpos;
|
---|
1482 | if(lp_>0) cout<<"...density move from "
|
---|
1483 | <<val_by_mpc3*dVol_<<" to "<<dn<<" / pixel"<<endl;
|
---|
1484 |
|
---|
1485 | double sum = 0.;
|
---|
1486 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1487 | int_8 ip = IndexR(i,j,l);
|
---|
1488 | if(data_[ip]<=0.) data_[ip] = 0.;
|
---|
1489 | else {
|
---|
1490 | data_[ip] *= dn;
|
---|
1491 | sum += data_[ip];
|
---|
1492 | }
|
---|
1493 | }
|
---|
1494 |
|
---|
1495 | if(lp_>0) cout<<"...quantity put into survey "<<sum<<" / "<<val_by_mpc3*Vol_<<endl;
|
---|
1496 |
|
---|
1497 | return sum;
|
---|
1498 | }
|
---|
1499 |
|
---|
1500 | double GeneFluct3D::ApplyPoisson(void)
|
---|
1501 | // do NOT treate negative or nul mass -> let it as it is
|
---|
1502 | {
|
---|
1503 | if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
|
---|
1504 | check_array_alloc();
|
---|
1505 |
|
---|
1506 | double sum = 0.;
|
---|
1507 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1508 | int_8 ip = IndexR(i,j,l);
|
---|
1509 | double v = data_[ip];
|
---|
1510 | if(v>0.) {
|
---|
1511 | uint_8 dn = PoissonRand(v,10.);
|
---|
1512 | data_[ip] = (double)dn;
|
---|
1513 | sum += (double)dn;
|
---|
1514 | }
|
---|
1515 | }
|
---|
1516 | if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
|
---|
1517 |
|
---|
1518 | return sum;
|
---|
1519 | }
|
---|
1520 |
|
---|
1521 | double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
|
---|
1522 | // do NOT treate negative or nul mass -> let it as it is
|
---|
1523 | // INPUT:
|
---|
1524 | // massdist : distribution de masse (m*dn/dm)
|
---|
1525 | // axeslog = false : retourne la masse
|
---|
1526 | // = true : retourne le log10(mass)
|
---|
1527 | // RETURN la masse totale
|
---|
1528 | {
|
---|
1529 | if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
|
---|
1530 | check_array_alloc();
|
---|
1531 |
|
---|
1532 | double sum = 0.;
|
---|
1533 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1534 | int_8 ip = IndexR(i,j,l);
|
---|
1535 | double v = data_[ip];
|
---|
1536 | if(v>0.) {
|
---|
1537 | long ngal = long(v+0.1);
|
---|
1538 | data_[ip] = 0.;
|
---|
1539 | for(long i=0;i<ngal;i++) {
|
---|
1540 | double m = massdist.RandomInterp(); // massdist.Random();
|
---|
1541 | if(axeslog) m = pow(10.,m);
|
---|
1542 | data_[ip] += m;
|
---|
1543 | }
|
---|
1544 | sum += data_[ip];
|
---|
1545 | }
|
---|
1546 | }
|
---|
1547 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
1548 |
|
---|
1549 | return sum;
|
---|
1550 | }
|
---|
1551 |
|
---|
1552 | double GeneFluct3D::TurnNGal2MassQuick(SchechterMassDist& schmdist)
|
---|
1553 | // idem TurnNGal2Mass mais beaucoup plus rapide
|
---|
1554 | {
|
---|
1555 | if(lp_>0) cout<<"--- TurnNGal2MassQuick ---"<<endl;
|
---|
1556 | check_array_alloc();
|
---|
1557 |
|
---|
1558 | double sum = 0.;
|
---|
1559 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1560 | int_8 ip = IndexR(i,j,l);
|
---|
1561 | double v = data_[ip];
|
---|
1562 | if(v>0.) {
|
---|
1563 | long ngal = long(v+0.1);
|
---|
1564 | data_[ip] = schmdist.TirMass(ngal);
|
---|
1565 | sum += data_[ip];
|
---|
1566 | }
|
---|
1567 | }
|
---|
1568 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
1569 |
|
---|
1570 | return sum;
|
---|
1571 | }
|
---|
1572 |
|
---|
1573 | void GeneFluct3D::AddNoise2Real(double snoise,int type_evol)
|
---|
1574 | // add noise to every pixels (meme les <=0 !)
|
---|
1575 | // type_evol = 0 : pas d'evolution de la puissance du bruit
|
---|
1576 | // 1 : evolution de la puissance du bruit avec la distance a l'observateur
|
---|
1577 | // 2 : evolution de la puissance du bruit avec la distance du plan Z
|
---|
1578 | // (tous les plans Z sont mis au meme redshift z de leur milieu)
|
---|
1579 | {
|
---|
1580 | if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<" evol="<<type_evol<<endl;
|
---|
1581 | check_array_alloc();
|
---|
1582 |
|
---|
1583 | if(type_evol<0) type_evol = 0;
|
---|
1584 | if(type_evol>2) {
|
---|
1585 | const char *bla = "GeneFluct3D::AddNoise2Real_Error: bad type_evol value";
|
---|
1586 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1587 | }
|
---|
1588 |
|
---|
1589 | vector<double> correction;
|
---|
1590 | InterpFunc *intercor = NULL;
|
---|
1591 |
|
---|
1592 | if(type_evol>0) {
|
---|
1593 | // Sigma_Noise(en mass) :
|
---|
1594 | // Slim ~ 1/sqrt(DNu) * sqrt(nlobe) en W/m^2Hz
|
---|
1595 | // Flim ~ sqrt(DNu) * sqrt(nlobe) en W/m^2
|
---|
1596 | // Mlim ~ sqrt(DNu) * (Dlum)^2 * sqrt(nlobe) en Msol
|
---|
1597 | // nlobe ~ 1/Dtrcom^2
|
---|
1598 | // Mlim ~ sqrt(DNu) * (Dlum)^2 / Dtrcom
|
---|
1599 | if(cosmo_ == NULL || redsh_ref_<0. || loscom2zred_.size()<1) {
|
---|
1600 | const char *bla = "GeneFluct3D::AddNoise2Real_Error: set Observator and Cosmology first";
|
---|
1601 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1602 | }
|
---|
1603 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
1604 | long nsz = loscom2zred_.size(), nszmod=((nsz>10)? nsz/10: 1);
|
---|
1605 | for(long i=0;i<nsz;i++) {
|
---|
1606 | double d = interpinv.X(i);
|
---|
1607 | double zred = interpinv(d);
|
---|
1608 | double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
|
---|
1609 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
1610 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
1611 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
1612 | double corr = sqrt(dnu/dnu_ref_) * pow(dlum/dlum_ref_,2.) * dtrc_ref_/dtrc;
|
---|
1613 | if(lp_>0 && (i==0 || i==nsz-1 || i%nszmod==0))
|
---|
1614 | cout<<"i="<<i<<" d="<<d<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
|
---|
1615 | <<" dtrc="<<dtrc<<" dlum="<<dlum<<" -> cor="<<corr<<endl;
|
---|
1616 | correction.push_back(corr);
|
---|
1617 | }
|
---|
1618 | intercor = new InterpFunc(loscom2zred_min_,loscom2zred_max_,correction);
|
---|
1619 | }
|
---|
1620 |
|
---|
1621 | double corrlim[2] = {1.,1.};
|
---|
1622 | for(long i=0;i<Nx_;i++) {
|
---|
1623 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
1624 | for(long j=0;j<Ny_;j++) {
|
---|
1625 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
1626 | for(long l=0;l<Nz_;l++) {
|
---|
1627 | double corr = 1.;
|
---|
1628 | if(type_evol>0) {
|
---|
1629 | double dz = DZcom(l);
|
---|
1630 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
1631 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
1632 | corr = (*intercor)(dz);
|
---|
1633 | if(corr<corrlim[0]) corrlim[0]=corr; else if(corr>corrlim[1]) corrlim[1]=corr;
|
---|
1634 | }
|
---|
1635 | int_8 ip = IndexR(i,j,l);
|
---|
1636 | data_[ip] += snoise*corr*NorRand();
|
---|
1637 | }
|
---|
1638 | }
|
---|
1639 | }
|
---|
1640 | if(type_evol>0)
|
---|
1641 | cout<<"correction factor range: ["<<corrlim[0]<<","<<corrlim[1]<<"]"<<endl;
|
---|
1642 |
|
---|
1643 | if(intercor!=NULL) delete intercor;
|
---|
1644 | }
|
---|
1645 |
|
---|
1646 | } // Fin namespace SOPHYA
|
---|
1647 |
|
---|
1648 |
|
---|
1649 |
|
---|
1650 |
|
---|
1651 | /*********************************************************************
|
---|
1652 | void GeneFluct3D::AddAGN(double lfjy,double lsigma,double powlaw)
|
---|
1653 | // Add AGN flux into simulation:
|
---|
1654 | // --- Procedure:
|
---|
1655 | // 1. lancer "cmvdefsurv" avec les parametres du survey
|
---|
1656 | // (au redshift de reference du survey)
|
---|
1657 | // et recuperer l'angle solide "angsol sr" du pixel elementaire
|
---|
1658 | // au centre du cube.
|
---|
1659 | // 2. lancer "cmvtstagn" pour cet angle solide -> cmvtstagn.ppf
|
---|
1660 | // 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne
|
---|
1661 | // cad une moyenne <log10(S)> et un sigma "sig"
|
---|
1662 | // Attention: la distribution n'est pas gaussienne les "mean,sigma"
|
---|
1663 | // de l'histo ne sont pas vraiment ce que l'on veut
|
---|
1664 | // --- Limitations actuelle du code:
|
---|
1665 | // . les AGN sont supposes evoluer avec la meme loi de puissance pour tout theta,phi
|
---|
1666 | // . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee
|
---|
1667 | // . la distribution est approximee a une gaussienne
|
---|
1668 | // ... C'est une approximation pour un observateur loin du centre du cube
|
---|
1669 | // et pour un cube peu epais / distance observateur
|
---|
1670 | // --- Parametres de la routine:
|
---|
1671 | // llfy : c'est le <log10(S)> du flux depose par les AGN
|
---|
1672 | // dans l'angle solide du pixel elementaire de reference du cube
|
---|
1673 | // lsigma : c'est le sigma de la distribution des log10(S)
|
---|
1674 | // powlaw : c'est la pente de la distribution cad que le flux "lmsol"
|
---|
1675 | // et considere comme le flux a 1.4GHz et qu'on suppose une loi
|
---|
1676 | // F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz)
|
---|
1677 | // - Comme on est en echelle log10():
|
---|
1678 | // on tire log10(Msol) + X
|
---|
1679 | // ou X est une realisation sur une gaussienne de variance "sig^2"
|
---|
1680 | // La masse realisee est donc: Msol*10^X
|
---|
1681 | // - Pas de probleme de pixel negatif car on a une multiplication!
|
---|
1682 | {
|
---|
1683 | if(lp_>0) cout<<"--- AddAGN: <log10(S Jy)> = "<<lfjy<<" , sigma = "<<lsigma<<endl;
|
---|
1684 | check_array_alloc();
|
---|
1685 |
|
---|
1686 | if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) {
|
---|
1687 | char *bla = "GeneFluct3D::AddAGN_Error: set Observator and Cosmology first";
|
---|
1688 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1689 | }
|
---|
1690 |
|
---|
1691 | // Le flux des AGN en Jy et en mass solaire
|
---|
1692 | double fagnref = pow(10.,lfjy)*(dnu_ref_*1.e9); // Jy.Hz = W/m^2
|
---|
1693 | double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlum_ref_); // Msol
|
---|
1694 | if(lp_>0)
|
---|
1695 | cout<<"Au pixel de ref: fagnref="<<fagnref
|
---|
1696 | <<" Jy.Hz (a 1.4GHz), magnref="<<magnref<<" Msol"<<endl;
|
---|
1697 |
|
---|
1698 | if(powlaw!=0.) {
|
---|
1699 | // F(nu) = F(1.4GHz)*(nu GHz/1.4 Ghz)^p = F(1.4GHz)*(1/(1+z))^p , car nu = 1.4 GHz/(1+z)
|
---|
1700 | magnref *= pow(1/(1.+redsh_ref_),powlaw);
|
---|
1701 | if(lp_>0) cout<<" powlaw="<<powlaw<<" -> change magnref to "<<magnref<<" Msol"<<endl;
|
---|
1702 | }
|
---|
1703 |
|
---|
1704 | // Les infos en fonction de l'indice "l" selon Oz
|
---|
1705 | vector<double> correction;
|
---|
1706 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
1707 | long nzmod = ((Nz_>10)?Nz_/10:1);
|
---|
1708 | for(long l=0;l<Nz_;l++) {
|
---|
1709 | double z = fabs(DZcom(l));
|
---|
1710 | double zred = interpinv(z);
|
---|
1711 | double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
|
---|
1712 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
1713 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
1714 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
1715 | // on a: Mass ~ DNu * Dlum^2 / Dtrcom^2
|
---|
1716 | double corr = dnu/dnu_ref_*pow(dtrc_ref_/dtrc*dlum/dlum_ref_,2.);
|
---|
1717 | // F(nu) = F(1.4GHz)*(nu GHz/1.4 Ghz)^p = F(1.4GHz)*(1/(1+z))^p , car nu = 1.4 GHz/(1+z)
|
---|
1718 | if(powlaw!=0.) corr *= pow((1.+redsh_ref_)/(1.+zred),powlaw);
|
---|
1719 | correction.push_back(corr);
|
---|
1720 | if(lp_>0 && (l==0 || l==Nz_-1 || l%nzmod==0)) {
|
---|
1721 | cout<<"l="<<l<<" z="<<z<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
|
---|
1722 | <<" dtrc="<<dtrc<<" dlum="<<dlum
|
---|
1723 | <<" -> cor="<<corr<<endl;
|
---|
1724 | }
|
---|
1725 | }
|
---|
1726 |
|
---|
1727 | double sum=0., sum2=0., nsum=0.;
|
---|
1728 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) {
|
---|
1729 | double a = lsigma*NorRand();
|
---|
1730 | a = magnref*pow(10.,a);
|
---|
1731 | // On met le meme tirage le long de Oz (indice k)
|
---|
1732 | for(long l=0;l<Nz_;l++) {
|
---|
1733 | int_8 ip = IndexR(i,j,l);
|
---|
1734 | data_[ip] += a*correction[l];
|
---|
1735 | }
|
---|
1736 | sum += a; sum2 += a*a; nsum += 1.;
|
---|
1737 | }
|
---|
1738 |
|
---|
1739 | if(lp_>0 && nsum>1.) {
|
---|
1740 | sum /= nsum;
|
---|
1741 | sum2 = sum2/nsum - sum*sum;
|
---|
1742 | cout<<"...Mean mass="<<sum<<" Msol , s^2="<<sum2<<" s="<<sqrt(fabs(sum2))<<endl;
|
---|
1743 | }
|
---|
1744 |
|
---|
1745 | }
|
---|
1746 | *********************************************************************/
|
---|