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

Last change on this file since 3743 was 3743, checked in by cmv, 16 years ago

on vire le SetMemoryMapping(CMemoryMapping), cmv 08/02/2010

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