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

Last change on this file since 3336 was 3331, checked in by cmv, 18 years ago

mise en place evolution seulement sur distance au plan Z , cmv 02/10/2007

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