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

Last change on this file since 3350 was 3349, checked in by cmv, 18 years ago
  • gros changements dans la structure de la classe GeneFluct3D (constructeur et logique d'aloocation memoire, init_fftw etc...)
  • suppression des valeurs de masse<0 mises a -999. directement mises a zero
  • suppression de TurnMass2HIMass qui fait maintenant la meme chose que TurnMass2MeanNumber
  • legere restructuration de cmvobserv3d.cc pour compat.

cmv 11/10/2007

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