source: Sophya/trunk/Cosmo/SimLSS/genefluct3d.cc@ 3773

Last change on this file since 3773 was 3773, checked in by cmv, 15 years ago

speed-up FITS file read, cmv 09/05/2010

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