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

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

intro ComputeSpectrum avec soustraction de bruit et deconvolution par la fct de transfert du pixel cmv 01/10/2007

File size: 46.5 KB
RevLine 
[3115]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
[3141]14#include "fabtcolread.h"
15#include "fabtwriter.h"
16#include "fioarr.h"
17
18#include "arrctcast.h"
19
[3115]20#include "constcosmo.h"
21#include "geneutils.h"
[3199]22#include "schechter.h"
[3115]23
24#include "genefluct3d.h"
25
[3281]26#define FFTW_THREAD
[3115]27
28#define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag()))
29
[3325]30namespace SOPHYA {
31
[3115]32//-------------------------------------------------------
[3141]33GeneFluct3D::GeneFluct3D(TArray< complex<r_8 > >& T)
[3271]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.)
[3199]42 , loscom2zred_min_(0.) , loscom2zred_max_(0.)
[3115]43{
[3157]44 xobs_[0] = xobs_[1] = xobs_[2] = 0.;
45 zred_.resize(0);
46 loscom_.resize(0);
[3199]47 loscom2zred_.resize(0);
[3115]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//-------------------------------------------------------
[3129]61void GeneFluct3D::SetSize(long nx,long ny,long nz,double dx,double dy,double dz)
[3115]62{
[3141]63 setsize(nx,ny,nz,dx,dy,dz);
64 setalloc();
65 setpointers(false);
[3154]66 init_fftw();
[3141]67}
68
[3154]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
[3267]77// Si kredshref<0 alors kredshref=nz/2 (milieu du cube)
[3157]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
[3154]80{
81 if(redshref<0.) redshref = 0.;
[3267]82 if(kredshref<0.) {
83 if(Nz_<=0) {
[3271]84 char *bla = "GeneFluct3D::SetObservator_Error: for kredsh_ref<0 SetSize should be called first";
[3267]85 cout<<bla<<endl; throw ParmError(bla);
86 }
87 kredshref = Nz_/2.;
88 }
[3271]89 redsh_ref_ = redshref;
90 kredsh_ref_ = kredshref;
[3199]91 if(lp_>0)
[3271]92 cout<<"--- GeneFluct3D::SetObservator zref="<<redsh_ref_<<" kref="<<kredsh_ref_<<endl;
[3154]93}
94
[3157]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
[3141]106void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz)
107{
[3155]108 if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<<nx<<","<<ny<<","<<nz
109 <<" D="<<dx<<","<<dy<<","<<dz<<endl;
[3141]110 if(nx<=0 || dx<=0.) {
[3267]111 char *bla = "GeneFluct3D::setsize_Error: bad value(s) for nn/dx";
[3199]112 cout<<bla<<endl; throw ParmError(bla);
[3115]113 }
114
[3141]115 // Les tailles des tableaux
[3115]116 Nx_ = nx;
117 Ny_ = ny; if(Ny_ <= 0) Ny_ = Nx_;
118 Nz_ = nz; if(Nz_ <= 0) Nz_ = Nx_;
[3141]119 N_.resize(0); N_.push_back(Nx_); N_.push_back(Ny_); N_.push_back(Nz_);
[3115]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_;
[3141]128 D_.resize(0); D_.push_back(Dx_); D_.push_back(Dy_); D_.push_back(Dz_);
[3115]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_);
[3141]136 Dk_.resize(0); Dk_.push_back(Dkx_); Dk_.push_back(Dky_); Dk_.push_back(Dkz_);
[3115]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_;
[3141]143 Knyq_.resize(0); Knyq_.push_back(Knyqx_); Knyq_.push_back(Knyqy_); Knyq_.push_back(Knyqz_);
144}
[3115]145
[3141]146void GeneFluct3D::setalloc(void)
147{
[3155]148 if(lp_>1) cout<<"--- GeneFluct3D::setalloc ---"<<endl;
[3141]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;
[3255]156 if(lp_>1) cout<<" allocating: "<<T_.Size()*sizeof(complex<r_8>)/1.e6<<" Mo"<<endl;
[3141]157 } catch (...) {
[3155]158 cout<<"GeneFluct3D::setalloc_Error: Problem allocating T_"<<endl;
[3141]159 }
160 T_.SetMemoryMapping(BaseArray::CMemoryMapping);
[3115]161}
162
[3141]163void GeneFluct3D::setpointers(bool from_real)
164{
[3155]165 if(lp_>1) cout<<"--- GeneFluct3D::setpointers ---"<<endl;
[3141]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");
[3199]179 cout<<bla<<endl; throw ParmError(bla);
[3141]180}
181
[3154]182void GeneFluct3D::init_fftw(void)
183{
184 // --- Initialisation de fftw3 (attention data est sur-ecrit a l'init)
[3155]185 if(lp_>1) cout<<"--- GeneFluct3D::init_fftw ---"<<endl;
[3154]186#ifdef FFTW_THREAD
187 if(nthread_>0) {
[3155]188 cout<<"...Computing with "<<nthread_<<" threads"<<endl;
[3154]189 fftw_init_threads();
190 fftw_plan_with_nthreads(nthread_);
191 }
192#endif
[3155]193 if(lp_>1) cout<<"...forward plan"<<endl;
[3154]194 pf_ = fftw_plan_dft_r2c_3d(Nx_,Ny_,Nz_,data_,fdata_,FFTW_ESTIMATE);
[3155]195 if(lp_>1) cout<<"...backward plan"<<endl;
[3154]196 pb_ = fftw_plan_dft_c2r_3d(Nx_,Ny_,Nz_,fdata_,data_,FFTW_ESTIMATE);
197}
[3141]198
[3157]199//-------------------------------------------------------
[3199]200long GeneFluct3D::LosComRedshift(double zinc,long npoints)
[3157]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
[3199]207// -- Input:
208// zinc : redshift increment for computation
209// npoints : number of points required for inverting loscom -> zred
[3157]210//
211{
[3199]212 if(lp_>0) cout<<"--- LosComRedshift: zinc="<<zinc<<" , npoints="<<npoints<<endl;
[3154]213
[3271]214 if(cosmo_ == NULL || redsh_ref_<0.) {
[3199]215 char *bla = "GeneFluct3D::LosComRedshift_Error: set Observator and Cosmology first";
216 cout<<bla<<endl; throw ParmError(bla);
[3157]217 }
218
[3271]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
[3199]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
[3157]241 xobs_[0] = Nx_/2.*Dx_;
242 xobs_[1] = Ny_/2.*Dy_;
[3271]243 xobs_[2] = kredsh_ref_*Dz_ - loscom_ref_;
[3157]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
[3271]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
[3157]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++) {
[3271]267 double dx2 = DXcom(i*(Nx_-1)); dx2 *= dx2;
[3157]268 for(long j=0;j<=1;j++) {
[3271]269 double dy2 = DYcom(j*(Ny_-1)); dy2 *= dy2;
[3157]270 for(long k=0;k<=1;k++) {
[3271]271 double dz2 = DZcom(k*(Nz_-1)); dz2 *= dz2;
[3157]272 dz2 = sqrt(dx2+dy2+dz2);
273 if(dz2>loscom_max_) loscom_max_ = dz2;
274 }
275 }
276 }
277 if(lp_>0) {
[3271]278 cout<<"...zref="<<redsh_ref_<<" kzref="<<kredsh_ref_<<" losref="<<loscom_ref_<<" Mpc\n"
[3157]279 <<" xobs="<<xobs_[0]<<" , "<<xobs_[1]<<" , "<<xobs_[2]<<" Mpc "
280 <<" in_cube="<<obs_in_cube
[3271]281 <<" loscom_min="<<loscom_min_<<" loscom_max="<<loscom_max_<<" Mpc (com)"<<endl;
[3157]282 }
283
[3199]284 // Fill the corresponding vectors for loscom and zred
[3267]285 // Be shure to have one dlc <loscom_min and one >loscom_max
[3199]286 if(zinc<=0.) zinc = 0.01;
[3157]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;
[3199]293 if(dlc>loscom_max_) break; // on sort apres avoir stoque un dlc>dlcmax
[3157]294 }
295
296 if(lp_>0) {
[3199]297 long n = zred_.size();
298 cout<<"...zred/loscom tables[zinc="<<zinc<<"]: n="<<n;
[3157]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
[3199]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++)
[3330]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;
[3199]327 }
328
329 return zred_.size();
[3157]330}
331
[3115]332//-------------------------------------------------------
[3141]333void GeneFluct3D::WriteFits(string cfname,int bitpix)
334{
[3155]335 cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
[3141]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");
[3271]347 fwrt.WriteKey("ZREF",redsh_ref_," reference redshift");
348 fwrt.WriteKey("KZREF",kredsh_ref_," reference redshift on z axe");
[3141]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{
[3155]362 cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
[3141]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");
[3154]372 double zref = fimg.ReadKey("ZREF");
373 double kzref = fimg.ReadKey("KZREF");
[3141]374 setsize(nx,ny,nz,dx,dy,dz);
375 setpointers(true);
[3154]376 init_fftw();
377 SetObservator(zref,kzref);
[3330]378 array_allocated_ = true;
[3141]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{
[3155]392 cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
[3141]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_;
[3271]400 R_.Info()["ZREF"] = (r_8)redsh_ref_;
401 R_.Info()["KZREF"] = (r_8)kredsh_ref_;
[3141]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{
[3155]417 cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
[3141]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) {
[3262]424 cout<<" ...reading spectrum into TArray<complex <r_8> >"<<endl;
[3141]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 }
[3154]431 setpointers(from_real); // a mettre ici pour relire les DVInfo
[3141]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"];
[3154]438 r_8 zref = R_.Info()["ZREF"];
439 r_8 kzref = R_.Info()["KZREF"];
[3141]440 setsize(nx,ny,nz,dx,dy,dz);
[3154]441 init_fftw();
442 SetObservator(zref,kzref);
[3330]443 array_allocated_ = true;
[3141]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
[3281]454void GeneFluct3D::WriteSlicePPF(string cfname)
455// On ecrit 3 tranches du cube selon chaque axe
456{
[3283]457 cout<<"--- GeneFluct3D::WriteSlicePPF: Writing Cube Slices "<<cfname<<endl;
[3281]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;
[3289]469 sprintf(str,"z%ld",l);
[3281]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;
[3289]478 sprintf(str,"y%ld",j);
[3281]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;
[3289]487 sprintf(str,"x%ld",i);
[3281]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
[3141]502//-------------------------------------------------------
[3115]503void GeneFluct3D::Print(void)
504{
[3141]505 cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
[3115]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;
[3271]517 cout<<"Redshift "<<redsh_ref_<<" for z axe at k="<<kredsh_ref_<<endl;
[3115]518}
519
520//-------------------------------------------------------
[3141]521void GeneFluct3D::ComputeFourier0(GenericFunc& pk_at_z)
[3115]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
[3155]527 if(lp_>0) cout<<"--- ComputeFourier0: before gaussian filling ---"<<endl;
[3115]528 // On tient compte du pb de normalisation de FFTW3
529 double sntot = sqrt((double)NRtot_);
[3129]530 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]531 int_8 ip = IndexR(i,j,l);
532 data_[ip] = NorRand()/sntot;
[3115]533 }
534
535 // --- realisation d'un tableau de tirage gaussiens
[3155]536 if(lp_>0) cout<<"...before fft real ---"<<endl;
[3115]537 fftw_execute(pf_);
538
539 // --- On remplit avec une realisation
[3157]540 if(lp_>0) cout<<"...before Fourier realization filling"<<endl;
[3115]541 T_(0,0,0) = complex<r_8>(0.); // on coupe le continue et on l'initialise
[3129]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;
[3115]545 double kx = ii*Dkx_; kx *= kx;
[3155]546 if(lp_>0 && i%lmod==0) cout<<"i="<<i<<" ii="<<ii<<endl;
[3129]547 for(long j=0;j<Ny_;j++) {
548 long jj = (j>Ny_/2) ? Ny_-j : j;
[3115]549 double ky = jj*Dky_; ky *= ky;
[3129]550 for(long l=0;l<NCz_;l++) {
[3115]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
[3141]555 double pk = pk_at_z(k)/Vol_;
[3115]556 // ici pas de "/2" a cause de la remarque ci-dessus
557 T_(l,j,i) *= sqrt(pk);
558 }
559 }
560 }
561
[3155]562 if(lp_>0) cout<<"...computing power"<<endl;
[3115]563 double p = compute_power_carte();
[3155]564 if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
[3115]565
566}
567
568//-------------------------------------------------------
[3141]569void GeneFluct3D::ComputeFourier(GenericFunc& pk_at_z)
570// Calcule une realisation du spectre "pk_at_z"
[3115]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
[3155]583 if(lp_>0) cout<<"--- ComputeFourier ---"<<endl;
[3129]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;
[3115]587 double kx = ii*Dkx_; kx *= kx;
[3155]588 if(lp_>0 && i%lmod==0) cout<<"i="<<i<<" ii="<<ii<<endl;
[3129]589 for(long j=0;j<Ny_;j++) {
590 long jj = (j>Ny_/2) ? Ny_-j : j;
[3115]591 double ky = jj*Dky_; ky *= ky;
[3129]592 for(long l=0;l<NCz_;l++) {
[3115]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
[3141]597 double pk = pk_at_z(k)/Vol_;
[3115]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
[3155]607 if(lp_>0) cout<<"...computing power"<<endl;
[3115]608 double p = compute_power_carte();
[3155]609 if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
[3115]610
611}
612
[3129]613long GeneFluct3D::manage_coefficients(void)
[3115]614// Take into account the real and complexe conjugate coefficients
615// because we want a realization of a real data in real space
616{
[3155]617 if(lp_>1) cout<<"...managing coefficients"<<endl;
[3141]618 check_array_alloc();
[3115]619
620 // 1./ Le Continu et Nyquist sont reels
[3129]621 long nreal = 0;
622 for(long kk=0;kk<2;kk++) {
623 long k=0; // continu
[3115]624 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
[3129]625 for(long jj=0;jj<2;jj++) {
626 long j=0;
[3115]627 if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
[3129]628 for(long ii=0;ii<2;ii++) {
629 long i=0;
[3115]630 if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
[3141]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;
[3115]634 nreal++;
635 }
636 }
637 }
[3155]638 if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
[3115]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
[3129]643 long nconj1 = 0;
644 for(long kk=0;kk<2;kk++) {
645 long k=0; // continu
[3115]646 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
[3129]647 for(long jj=0;jj<2;jj++) { // selon j
648 long j=0;
[3115]649 if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
[3129]650 for(long i=1;i<(Nx_+1)/2;i++) {
[3141]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];
[3115]654 nconj1++;
655 }
656 }
[3129]657 for(long ii=0;ii<2;ii++) {
658 long i=0;
[3115]659 if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
[3129]660 for(long j=1;j<(Ny_+1)/2;j++) {
[3141]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];
[3115]664 nconj1++;
665 }
666 }
667 }
[3155]668 if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
[3115]669
670 // b./ les lignes et colonnes hors continu et de nyquist
[3129]671 long nconj2 = 0;
672 for(long kk=0;kk<2;kk++) {
673 long k=0; // continu
[3115]674 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
[3129]675 for(long j=1;j<(Ny_+1)/2;j++) {
[3115]676 if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
[3129]677 for(long i=1;i<Nx_;i++) {
[3115]678 if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
[3141]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];
[3115]682 nconj2++;
683 }
684 }
685 }
[3155]686 if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
[3115]687
[3155]688 if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(Nx_*Ny_*NCz_-nconj1-nconj2)-8<<endl;
[3115]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{
[3141]696 check_array_alloc();
697
[3115]698 double s2 = 0.;
[3129]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));
[3115]702
703 double s20 = 0.;
[3129]704 for(long j=0;j<Ny_;j++)
705 for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
[3115]706
707 double s2n = 0.;
708 if(Nz_%2==0)
[3129]709 for(long j=0;j<Ny_;j++)
710 for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
[3115]711
712 return 2.*s2 -s20 -s2n;
713}
714
715//-------------------------------------------------------------------
716void GeneFluct3D::FilterByPixel(void)
717// Filtrage par la fonction fenetre du pixel (parallelepipede)
[3120]718// TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
[3115]719// e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
[3120]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
[3115]723{
[3155]724 if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
[3141]725 check_array_alloc();
726
[3129]727 for(long i=0;i<Nx_;i++) {
728 long ii = (i>Nx_/2) ? Nx_-i : i;
[3120]729 double kx = ii*Dkx_ *Dx_/2;
[3330]730 double pk_x = pixelfilter(kx);
[3129]731 for(long j=0;j<Ny_;j++) {
732 long jj = (j>Ny_/2) ? Ny_-j : j;
[3120]733 double ky = jj*Dky_ *Dy_/2;
[3330]734 double pk_y = pixelfilter(ky);
[3129]735 for(long l=0;l<NCz_;l++) {
[3120]736 double kz = l*Dkz_ *Dz_/2;
[3141]737 double pk_z = pixelfilter(kz);
738 T_(l,j,i) *= pk_x*pk_y*pk_z;
[3115]739 }
740 }
741 }
742
743}
744
745//-------------------------------------------------------------------
[3199]746void GeneFluct3D::ApplyGrowthFactor(void)
[3157]747// Apply Growth to real space
748// Using the correspondance between redshift and los comoving distance
749// describe in vector "zred_" "loscom_"
750{
[3199]751 if(lp_>0) cout<<"--- ApplyGrowthFactor ---"<<endl;
[3157]752 check_array_alloc();
753
754 if(growth_ == NULL) {
[3199]755 char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
756 cout<<bla<<endl; throw ParmError(bla);
[3157]757 }
758
[3199]759 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
[3157]760 unsigned short ok;
761
762 //CHECK: Histo hgr(0.9*zred_[0],1.1*zred_[n-1],1000);
763 for(long i=0;i<Nx_;i++) {
764 double dx2 = xobs_[0] - i*Dx_; dx2 *= dx2;
765 for(long j=0;j<Ny_;j++) {
766 double dy2 = xobs_[1] - j*Dy_; dy2 *= dy2;
767 for(long l=0;l<Nz_;l++) {
768 double dz2 = xobs_[2] - l*Dz_; dz2 *= dz2;
769 dz2 = sqrt(dx2+dy2+dz2);
770 double z = interpinv(dz2);
771 //CHECK: hgr.Add(z);
772 double dzgr = (*growth_)(z); // interpolation par morceau
773 //double dzgr = growth_->Linear(z,ok); // interpolation lineaire
774 //double dzgr = growth_->Parab(z,ok); // interpolation parabolique
775 int_8 ip = IndexR(i,j,l);
776 data_[ip] *= dzgr;
777 }
778 }
779 }
780
781 //CHECK: {POutPersist pos("applygrowth.ppf"); string tag="hgr"; pos.PutObject(hgr,tag);}
782
783}
784
785//-------------------------------------------------------------------
[3115]786void GeneFluct3D::ComputeReal(void)
787// Calcule une realisation dans l'espace reel
788{
[3155]789 if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
[3141]790 check_array_alloc();
[3115]791
792 // On fait la FFT
793 fftw_execute(pb_);
794}
795
796//-------------------------------------------------------------------
797void GeneFluct3D::ReComputeFourier(void)
798{
[3155]799 if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
[3141]800 check_array_alloc();
[3115]801
802 // On fait la FFT
803 fftw_execute(pf_);
804 // On corrige du pb de la normalisation de FFTW3
805 double v = (double)NRtot_;
[3129]806 for(long i=0;i<Nx_;i++)
807 for(long j=0;j<Ny_;j++)
808 for(long l=0;l<NCz_;l++) T_(l,j,i) /= complex<r_8>(v);
[3115]809
810}
811
812//-------------------------------------------------------------------
[3141]813int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
814// Compute spectrum from "T" and fill HistoErr "herr"
[3115]815// T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
816// cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
817{
[3155]818 if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
[3141]819 check_array_alloc();
[3115]820
[3141]821 if(herr.NBins()<0) return -1;
822 herr.Zero();
[3115]823
824 // Attention a l'ordre
[3129]825 for(long i=0;i<Nx_;i++) {
826 long ii = (i>Nx_/2) ? Nx_-i : i;
[3115]827 double kx = ii*Dkx_; kx *= kx;
[3129]828 for(long j=0;j<Ny_;j++) {
829 long jj = (j>Ny_/2) ? Ny_-j : j;
[3115]830 double ky = jj*Dky_; ky *= ky;
[3129]831 for(long l=0;l<NCz_;l++) {
[3330]832 double kz = l*Dkz_;
833 double k = sqrt(kx+ky+kz*kz);
[3115]834 double pk = MODULE2(T_(l,j,i));
[3141]835 herr.Add(k,pk);
[3115]836 }
837 }
838 }
[3150]839 herr.ToVariance();
[3115]840
841 // renormalize to directly compare to original spectrum
842 double norm = Vol_;
[3141]843 herr *= norm;
[3115]844
845 return 0;
846}
847
[3141]848int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
849{
[3155]850 if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
[3141]851 check_array_alloc();
852
853 if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
854 herr.Zero();
855
856 // Attention a l'ordre
857 for(long i=0;i<Nx_;i++) {
858 long ii = (i>Nx_/2) ? Nx_-i : i;
859 double kx = ii*Dkx_; kx *= kx;
860 for(long j=0;j<Ny_;j++) {
861 long jj = (j>Ny_/2) ? Ny_-j : j;
862 double ky = jj*Dky_; ky *= ky;
863 double kt = sqrt(kx+ky);
864 for(long l=0;l<NCz_;l++) {
865 double kz = l*Dkz_;
866 double pk = MODULE2(T_(l,j,i));
867 herr.Add(kt,kz,pk);
868 }
869 }
870 }
[3150]871 herr.ToVariance();
[3141]872
873 // renormalize to directly compare to original spectrum
874 double norm = Vol_;
875 herr *= norm;
876
877 return 0;
878}
879
[3330]880//-------------------------------------------------------------------
881int GeneFluct3D::ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor)
882// Compute spectrum from "T" and fill HistoErr "herr"
883// AVEC la soustraction du niveau de bruit et la correction par filterpixel()
884// Si on ne fait pas ca, alors on obtient un spectre non-isotrope!
885//
886// T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
887// cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
888{
889 if(lp_>0) cout<<"--- ComputeSpectrum: sigma="<<sigma<<endl;
890 check_array_alloc();
891
892 if(sigma<=0.) sigma = 0.;
893 double sigma2 = sigma*sigma / (double)NRtot_;
894
895 if(herr.NBins()<0) return -1;
896 herr.Zero();
897
898 TVector<r_8> vfz(NCz_);
899 if(pixcor) // kz = l*Dkz_
900 for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
901
902 // Attention a l'ordre
903 for(long i=0;i<Nx_;i++) {
904 long ii = (i>Nx_/2) ? Nx_-i : i;
905 double kx = ii*Dkx_;
906 double fx = (pixcor) ? pixelfilter(kx*Dx_/2): 1.;
907 kx *= kx; fx *= fx;
908 for(long j=0;j<Ny_;j++) {
909 long jj = (j>Ny_/2) ? Ny_-j : j;
910 double ky = jj*Dky_;
911 double fy = (pixcor) ? pixelfilter(ky*Dy_/2): 1.;
912 ky *= ky; fy *= fy;
913 for(long l=0;l<NCz_;l++) {
914 double kz = l*Dkz_;
915 double k = sqrt(kx+ky+kz*kz);
916 double pk = MODULE2(T_(l,j,i)) - sigma2;
917 double fz = (pixcor) ? vfz(l): 1.;
918 double f = fx*fy*fz;
919 if(f>0.) herr.Add(k,pk/f);
920 }
921 }
922 }
923 herr.ToVariance();
924
925 // renormalize to directly compare to original spectrum
926 double norm = Vol_;
927 herr *= norm;
928
929 return 0;
930}
931
932int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor)
933// AVEC la soustraction du niveau de bruit et la correction par filterpixel()
934{
935 if(lp_>0) cout<<"--- ComputeSpectrum2D: sigma="<<sigma<<endl;
936 check_array_alloc();
937
938 if(sigma<=0.) sigma = 0.;
939 double sigma2 = sigma*sigma / (double)NRtot_;
940
941 if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
942 herr.Zero();
943
944 TVector<r_8> vfz(NCz_);
945 if(pixcor) // kz = l*Dkz_
946 for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
947
948 // Attention a l'ordre
949 for(long i=0;i<Nx_;i++) {
950 long ii = (i>Nx_/2) ? Nx_-i : i;
951 double kx = ii*Dkx_;
952 double fx = (pixcor) ? pixelfilter(kx*Dx_/2) : 1.;
953 kx *= kx; fx *= fx;
954 for(long j=0;j<Ny_;j++) {
955 long jj = (j>Ny_/2) ? Ny_-j : j;
956 double ky = jj*Dky_;
957 double fy = (pixcor) ? pixelfilter(ky*Dy_/2) : 1.;
958 ky *= ky; fy *= fy;
959 double kt = sqrt(kx+ky);
960 for(long l=0;l<NCz_;l++) {
961 double kz = l*Dkz_;
962 double pk = MODULE2(T_(l,j,i)) - sigma2;
963 double fz = (pixcor) ? vfz(l): 1.;
964 double f = fx*fy*fz;
965 if(f>0.) herr.Add(kt,kz,pk/f);
966 }
967 }
968 }
969 herr.ToVariance();
970
971 // renormalize to directly compare to original spectrum
972 double norm = Vol_;
973 herr *= norm;
974
975 return 0;
976}
977
[3115]978//-------------------------------------------------------
[3134]979int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
[3115]980// Recompute MASS variance in spherical top-hat (rayon=R)
981{
[3262]982 if(lp_>0) cout<<"--- VarianceFrReal R="<<R<<endl;
[3141]983 check_array_alloc();
984
[3129]985 long dnx = long(R/Dx_+0.5); if(dnx<=0) dnx = 1;
986 long dny = long(R/Dy_+0.5); if(dny<=0) dny = 1;
987 long dnz = long(R/Dz_+0.5); if(dnz<=0) dnz = 1;
[3155]988 if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
[3115]989
[3134]990 double sum=0., sum2=0., r2 = R*R; int_8 nsum=0;
[3115]991
[3129]992 for(long i=dnx;i<Nx_-dnx;i+=dnx) {
993 for(long j=dny;j<Ny_-dny;j+=dny) {
994 for(long l=dnz;l<Nz_-dnz;l+=dnz) {
[3134]995 double s=0.; int_8 n=0;
[3129]996 for(long ii=i-dnx;ii<=i+dnx;ii++) {
[3115]997 double x = (ii-i)*Dx_; x *= x;
[3129]998 for(long jj=j-dny;jj<=j+dny;jj++) {
[3115]999 double y = (jj-j)*Dy_; y *= y;
[3129]1000 for(long ll=l-dnz;ll<=l+dnz;ll++) {
[3115]1001 double z = (ll-l)*Dz_; z *= z;
1002 if(x+y+z>r2) continue;
[3141]1003 int_8 ip = IndexR(ii,jj,ll);
1004 s += 1.+data_[ip];
[3115]1005 n++;
1006 }
1007 }
1008 }
1009 if(n>0) {sum += s; sum2 += s*s; nsum++;}
1010 //cout<<i<<","<<j<<","<<l<<" n="<<n<<" s="<<s<<" sum="<<sum<<" sum2="<<sum2<<endl;
1011 }
1012 }
1013 }
1014
1015 if(nsum<=1) {var=0.; return nsum;}
1016
1017 sum /= nsum;
1018 sum2 = sum2/nsum - sum*sum;
[3262]1019 if(lp_>0) cout<<"...nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
[3115]1020
1021 var = sum2/(sum*sum); // <dM>^2/<M>^2
[3262]1022 if(lp_>0) cout<<"...sigmaR^2="<<var<<" -> "<<sqrt(var)<<endl;
[3115]1023
1024 return nsum;
1025}
1026
1027//-------------------------------------------------------
[3134]1028int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
[3115]1029// number of pixels outside of ]vmin,vmax[ extremites exclues
1030// -> vmin and vmax are considered as bad
1031{
[3141]1032 check_array_alloc();
[3115]1033
[3134]1034 int_8 nbad = 0;
[3129]1035 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1036 int_8 ip = IndexR(i,j,l);
1037 double v = data_[ip];
[3115]1038 if(v<=vmin || v>=vmax) nbad++;
1039 }
1040
[3262]1041 if(lp_>0) cout<<"--- NumberOfBad "<<nbad<<" px out of ]"<<vmin<<","<<vmax<<"["<<endl;
[3115]1042 return nbad;
1043}
1044
[3320]1045int_8 GeneFluct3D::MinMax(double& xmin,double& xmax,double vmin,double vmax)
1046// Calcul des valeurs xmin et xmax dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
1047{
1048 bool tstval = (vmax>vmin)? true: false;
1049 if(lp_>0) {
1050 cout<<"--- MinMax";
1051 if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
1052 cout<<endl;
1053 }
1054 check_array_alloc();
1055
1056 int_8 n = 0;
1057 xmin = xmax = data_[0];
1058
1059 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1060 int_8 ip = IndexR(i,j,l);
1061 double x = data_[ip];
1062 if(tstval && (x<=vmin || x>=vmax)) continue;
1063 if(x<xmin) xmin = x;
1064 if(x>xmax) xmax = x;
1065 n++;
1066 }
1067
1068 if(lp_>0) cout<<" n="<<n<<" min="<<xmin<<" max="<<xmax<<endl;
1069
1070 return n;
1071}
1072
[3261]1073int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax
1074 ,bool useout,double vout)
1075// Calcul de mean,sigma2 dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
1076// useout = false: ne pas utiliser les pixels hors limites pour calculer mean,sigma2
1077// true : utiliser les pixels hors limites pour calculer mean,sigma2
1078// en remplacant leurs valeurs par "vout"
[3115]1079{
[3261]1080 bool tstval = (vmax>vmin)? true: false;
1081 if(lp_>0) {
[3262]1082 cout<<"--- MeanSigma2";
1083 if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
[3261]1084 if(useout) cout<<", useout="<<useout<<" vout="<<vout;
1085 cout<<endl;
1086 }
[3141]1087 check_array_alloc();
[3115]1088
[3134]1089 int_8 n = 0;
[3115]1090 rm = rs2 = 0.;
1091
[3129]1092 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1093 int_8 ip = IndexR(i,j,l);
1094 double v = data_[ip];
[3261]1095 if(tstval) {
1096 if(v<=vmin || v>=vmax) {if(useout) v=vout; else continue;}
1097 }
[3115]1098 rm += v;
1099 rs2 += v*v;
1100 n++;
1101 }
1102
1103 if(n>1) {
1104 rm /= (double)n;
1105 rs2 = rs2/(double)n - rm*rm;
1106 }
1107
[3261]1108 if(lp_>0) cout<<" n="<<n<<" m="<<rm<<" s2="<<rs2<<" s="<<sqrt(fabs(rs2))<<endl;
1109
[3115]1110 return n;
1111}
1112
[3134]1113int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
[3115]1114// set to "val0" if out of range ]vmin,vmax[ extremites exclues
[3261]1115// cad set to "val0" if in [vmin,vmax] -> vmin and vmax are set to val0
[3115]1116{
[3141]1117 check_array_alloc();
[3115]1118
[3134]1119 int_8 nbad = 0;
[3129]1120 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1121 int_8 ip = IndexR(i,j,l);
1122 double v = data_[ip];
1123 if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
[3115]1124 }
1125
[3262]1126 if(lp_>0) cout<<"--- SetToVal "<<nbad<<" px set to="<<val0
1127 <<" because out of range=]"<<vmin<<","<<vmax<<"["<<endl;
[3115]1128 return nbad;
1129}
1130
[3283]1131void GeneFluct3D::ScaleOffset(double scalecube,double offsetcube)
1132// Replace "V" by "scalecube * ( V + offsetcube )"
1133{
[3284]1134 if(lp_>0) cout<<"--- ScaleCube scale="<<scalecube<<" offset="<<offsetcube<<endl;
[3283]1135
1136 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1137 int_8 ip = IndexR(i,j,l);
1138 data_[ip] = scalecube * ( data_[ip] + offsetcube );
1139 }
1140
1141 return;
1142}
1143
[3115]1144//-------------------------------------------------------
1145void GeneFluct3D::TurnFluct2Mass(void)
1146// d_rho/rho -> Mass (add one!)
1147{
[3155]1148 if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
[3141]1149 check_array_alloc();
1150
[3115]1151
[3129]1152 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1153 int_8 ip = IndexR(i,j,l);
1154 data_[ip] += 1.;
[3115]1155 }
1156}
1157
[3329]1158double GeneFluct3D::TurnMass2HIMass(double m_by_mpc3)
1159{
1160 if(lp_>0) cout<<"--- TurnMass2HIMass : "<<m_by_mpc3<<" Msol/Mpc^3"<<endl;
1161
1162 double mall=0., mgood=0.;
1163 int_8 nall=0, ngood=0;
1164 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1165 int_8 ip = IndexR(i,j,l);
1166 mall += data_[ip]; nall++;
1167 if(data_[ip]>0.) {mgood += data_[ip]; ngood++;}
1168 }
1169 if(ngood>0) mgood /= (double)ngood;
1170 if(nall>0) mall /= (double)nall;
1171 if(lp_>0) cout<<"...ngood="<<ngood<<" mgood="<<mgood
1172 <<", nall="<<nall<<" mall="<<mall<<endl;
1173 if(ngood<=0 || mall<=0.) {
1174 cout<<"TurnMass2HIMass_Error: ngood="<<ngood<<" <=0 || mall="<<mall<<" <=0"<<endl;
1175 throw RangeCheckError("TurnMass2HIMass_Error: ngood<=0 || mall<=0");
1176 }
1177
1178 // On doit mettre m*Vol masse de HI dans notre survey
1179 // On en met uniquement dans les pixels de masse >0.
1180 // On MET a zero les pixels <0
1181 // On renormalise sur les pixels>0 pour qu'on ait m*Vol masse de HI
1182 // comme on ne prend que les pixels >0, on doit normaliser
1183 // a la moyenne de <1+d_rho/rho> sur ces pixels
1184 // (rappel sur tout les pixels <1+d_rho/rho>=1)
1185 // masse de HI a mettre ds 1 px:
1186 double dm = m_by_mpc3*Vol_/ (mgood/mall) /(double)ngood;
1187 if(lp_>0) cout<<"...HI mass density move from "
1188 <<m_by_mpc3*Vol_/double(NRtot_)<<" to "<<dm<<" / pixel"<<endl;
1189
1190 double sum = 0.;
1191 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1192 int_8 ip = IndexR(i,j,l);
1193 if(data_[ip]<=0.) data_[ip] = 0.;
1194 else {
1195 data_[ip] *= dm; // par coherence on multiplie aussi les <=0
1196 sum += data_[ip]; // mais on ne les compte pas
1197 }
1198 }
1199
1200 if(lp_>0) cout<<sum<<"...mass HI put into survey / "<<m_by_mpc3*Vol_<<endl;
1201
1202 return sum;
1203}
1204
[3115]1205double GeneFluct3D::TurnMass2MeanNumber(double n_by_mpc3)
1206// do NOT treate negative or nul values
1207{
[3329]1208 if(lp_>0) cout<<"--- TurnMass2MeanNumber : "<<n_by_mpc3<<" gal/Mpc^3"<<endl;
[3115]1209
[3262]1210 double mall=0., mgood=0.;
1211 int_8 nall=0, ngood=0;
1212 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1213 int_8 ip = IndexR(i,j,l);
1214 mall += data_[ip]; nall++;
1215 if(data_[ip]>0.) {mgood += data_[ip]; ngood++;}
1216 }
1217 if(ngood>0) mgood /= (double)ngood;
1218 if(nall>0) mall /= (double)nall;
1219 if(lp_>0) cout<<"...ngood="<<ngood<<" mgood="<<mgood
1220 <<", nall="<<nall<<" mall="<<mall<<endl;
1221 if(ngood<=0 || mall<=0.) {
1222 cout<<"TurnMass2MeanNumber_Error: ngood="<<ngood<<" <=0 || mall="<<mall<<" <=0"<<endl;
1223 throw RangeCheckError("TurnMass2MeanNumber_Error: ngood<=0 || mall<=0");
1224 }
[3115]1225
1226 // On doit mettre n*Vol galaxies dans notre survey
1227 // On en met uniquement dans les pixels de masse >0.
1228 // On NE met PAS a zero les pixels <0
1229 // On renormalise sur les pixels>0 pour qu'on ait n*Vol galaxies
1230 // comme on ne prend que les pixels >0, on doit normaliser
1231 // a la moyenne de <1+d_rho/rho> sur ces pixels
1232 // (rappel sur tout les pixels <1+d_rho/rho>=1)
[3262]1233 // nb de gal a mettre ds 1 px:
1234 double dn = n_by_mpc3*Vol_/ (mgood/mall) /(double)ngood;
[3155]1235 if(lp_>0) cout<<"...galaxy density move from "
1236 <<n_by_mpc3*Vol_/double(NRtot_)<<" to "<<dn<<" / pixel"<<endl;
[3262]1237
[3115]1238 double sum = 0.;
[3129]1239 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1240 int_8 ip = IndexR(i,j,l);
1241 data_[ip] *= dn; // par coherence on multiplie aussi les <=0
1242 if(data_[ip]>0.) sum += data_[ip]; // mais on ne les compte pas
[3115]1243 }
[3262]1244
[3155]1245 if(lp_>0) cout<<sum<<"...galaxies put into survey / "<<n_by_mpc3*Vol_<<endl;
[3115]1246
1247 return sum;
1248}
1249
1250double GeneFluct3D::ApplyPoisson(void)
1251// do NOT treate negative or nul mass -> let it as it is
1252{
[3155]1253 if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
[3141]1254 check_array_alloc();
1255
[3115]1256 double sum = 0.;
[3129]1257 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1258 int_8 ip = IndexR(i,j,l);
1259 double v = data_[ip];
[3115]1260 if(v>0.) {
1261 unsigned long dn = PoissRandLimit(v,10.);
[3141]1262 data_[ip] = (double)dn;
[3115]1263 sum += (double)dn;
1264 }
1265 }
[3155]1266 if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
[3115]1267
1268 return sum;
1269}
1270
1271double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
1272// do NOT treate negative or nul mass -> let it as it is
1273// INPUT:
1274// massdist : distribution de masse (m*dn/dm)
1275// axeslog = false : retourne la masse
1276// = true : retourne le log10(mass)
1277// RETURN la masse totale
1278{
[3155]1279 if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
[3141]1280 check_array_alloc();
1281
[3115]1282 double sum = 0.;
[3129]1283 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1284 int_8 ip = IndexR(i,j,l);
1285 double v = data_[ip];
[3115]1286 if(v>0.) {
[3129]1287 long ngal = long(v+0.1);
[3141]1288 data_[ip] = 0.;
[3129]1289 for(long i=0;i<ngal;i++) {
[3115]1290 double m = massdist.RandomInterp(); // massdist.Random();
1291 if(axeslog) m = pow(10.,m);
[3141]1292 data_[ip] += m;
[3115]1293 }
[3141]1294 sum += data_[ip];
[3115]1295 }
1296 }
[3155]1297 if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
[3115]1298
1299 return sum;
1300}
1301
[3320]1302double GeneFluct3D::TurnNGal2MassQuick(SchechterMassDist& schmdist)
1303// idem TurnNGal2Mass mais beaucoup plus rapide
1304{
1305 if(lp_>0) cout<<"--- TurnNGal2MassQuick ---"<<endl;
1306 check_array_alloc();
1307
1308 double sum = 0.;
1309 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1310 int_8 ip = IndexR(i,j,l);
1311 double v = data_[ip];
1312 if(v>0.) {
1313 long ngal = long(v+0.1);
1314 data_[ip] = schmdist.TirMass(ngal);
1315 sum += data_[ip];
1316 }
1317 }
1318 if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
1319
1320 return sum;
1321}
1322
[3199]1323void GeneFluct3D::AddAGN(double lfjy,double lsigma,double powlaw)
[3196]1324// Add AGN flux into simulation:
1325// --- Procedure:
1326// 1. lancer "cmvdefsurv" avec les parametres du survey
[3199]1327// (au redshift de reference du survey)
[3196]1328// et recuperer l'angle solide "angsol sr" du pixel elementaire
1329// au centre du cube.
1330// 2. lancer "cmvtstagn" pour cet angle solide -> cmvtstagn.ppf
1331// 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne
1332// cad une moyenne <log10(S)> et un sigma "sig"
[3199]1333// Attention: la distribution n'est pas gaussienne les "mean,sigma"
1334// de l'histo ne sont pas vraiment ce que l'on veut
[3196]1335// --- Limitations actuelle du code:
[3271]1336// . les AGN sont supposes evoluer avec la meme loi de puissance pour tout theta,phi
[3199]1337// . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee
1338// . la distribution est approximee a une gaussienne
1339// ... C'est une approximation pour un observateur loin du centre du cube
1340// et pour un cube peu epais / distance observateur
[3196]1341// --- Parametres de la routine:
[3271]1342// llfy : c'est le <log10(S)> du flux depose par les AGN
1343// dans l'angle solide du pixel elementaire de reference du cube
1344// lsigma : c'est le sigma de la distribution des log10(S)
1345// powlaw : c'est la pente de la distribution cad que le flux "lmsol"
[3199]1346// et considere comme le flux a 1.4GHz et qu'on suppose une loi
1347// F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz)
[3196]1348// - Comme on est en echelle log10():
1349// on tire log10(Msol) + X
1350// ou X est une realisation sur une gaussienne de variance "sig^2"
1351// La masse realisee est donc: Msol*10^X
1352// - Pas de probleme de pixel negatif car on a une multiplication!
1353{
[3199]1354 if(lp_>0) cout<<"--- AddAGN: <log10(S Jy)> = "<<lfjy<<" , sigma = "<<lsigma<<endl;
[3196]1355 check_array_alloc();
1356
[3271]1357 if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) {
[3199]1358 char *bla = "GeneFluct3D::AddAGN_Error: set Observator and Cosmology first";
1359 cout<<bla<<endl; throw ParmError(bla);
1360 }
[3196]1361
[3271]1362 // Le flux des AGN en Jy et en mass solaire
1363 double fagnref = pow(10.,lfjy)*(dnu_ref_*1.e9); // Jy.Hz = W/m^2
1364 double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlum_ref_); // Msol
1365 if(lp_>0)
1366 cout<<"Au pixel de ref: fagnref="<<fagnref
1367 <<" Jy.Hz (a 1.4GHz), magnref="<<magnref<<" Msol"<<endl;
[3196]1368
[3199]1369 if(powlaw!=0.) {
[3271]1370 // 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)
1371 magnref *= pow(1/(1.+redsh_ref_),powlaw);
[3199]1372 if(lp_>0) cout<<" powlaw="<<powlaw<<" -> change magnref to "<<magnref<<" Msol"<<endl;
1373 }
1374
1375 // Les infos en fonction de l'indice "l" selon Oz
1376 vector<double> correction;
1377 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
[3271]1378 long nzmod = ((Nz_>10)?Nz_/10:1);
[3199]1379 for(long l=0;l<Nz_;l++) {
[3271]1380 double z = fabs(DZcom(l));
[3199]1381 double zred = interpinv(z);
[3271]1382 double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
[3199]1383 double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
1384 double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
1385 double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
[3271]1386 // on a: Mass ~ DNu * Dlum^2 / Dtrcom^2
1387 double corr = dnu/dnu_ref_*pow(dtrc_ref_/dtrc*dlum/dlum_ref_,2.);
1388 // 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)
1389 if(powlaw!=0.) corr *= pow((1.+redsh_ref_)/(1.+zred),powlaw);
[3199]1390 correction.push_back(corr);
[3271]1391 if(lp_>0 && (l==0 || l==Nz_-1 || l%nzmod==0)) {
1392 cout<<"l="<<l<<" z="<<z<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
1393 <<" dtrc="<<dtrc<<" dlum="<<dlum
1394 <<" -> cor="<<corr<<endl;
[3199]1395 }
1396 }
1397
1398 double sum=0., sum2=0., nsum=0.;
1399 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) {
1400 double a = lsigma*NorRand();
1401 a = magnref*pow(10.,a);
1402 // On met le meme tirage le long de Oz (indice k)
1403 for(long l=0;l<Nz_;l++) {
1404 int_8 ip = IndexR(i,j,l);
1405 data_[ip] += a*correction[l];
1406 }
1407 sum += a; sum2 += a*a; nsum += 1.;
1408 }
1409
1410 if(lp_>0 && nsum>1.) {
[3196]1411 sum /= nsum;
1412 sum2 = sum2/nsum - sum*sum;
1413 cout<<"...Mean mass="<<sum<<" Msol , s^2="<<sum2<<" s="<<sqrt(fabs(sum2))<<endl;
1414 }
1415
1416}
1417
[3267]1418void GeneFluct3D::AddNoise2Real(double snoise,bool with_evol)
[3115]1419// add noise to every pixels (meme les <=0 !)
1420{
[3267]1421 if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<" evol="<<with_evol<<endl;
[3141]1422 check_array_alloc();
1423
[3271]1424 vector<double> correction;
1425 InterpFunc *intercor = NULL;
1426
1427 if(with_evol) {
1428 // Sigma_Noise(en mass) :
1429 // Slim ~ 1/sqrt(DNu) * sqrt(nlobe) en W/m^2Hz
1430 // Flim ~ sqrt(DNu) * sqrt(nlobe) en W/m^2
1431 // Mlim ~ sqrt(DNu) * (Dlum)^2 * sqrt(nlobe) en Msol
1432 // nlobe ~ 1/Dtrcom^2
1433 // Mlim ~ sqrt(DNu) * (Dlum)^2 / Dtrcom
1434 if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) {
1435 char *bla = "GeneFluct3D::AddNoise2Real_Error: set Observator and Cosmology first";
1436 cout<<bla<<endl; throw ParmError(bla);
1437 }
1438 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
1439 long nsz = loscom2zred_.size(), nszmod=((nsz>10)? nsz/10: 1);
1440 for(long i=0;i<nsz;i++) {
1441 double d = interpinv.X(i);
1442 double zred = interpinv(d);
1443 double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
1444 double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
1445 double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
1446 double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
1447 double corr = sqrt(dnu/dnu_ref_) * pow(dlum/dlum_ref_,2.) * dtrc_ref_/dtrc;
1448 if(lp_>0 && (i==0 || i==nsz-1 || i%nszmod==0))
1449 cout<<"i="<<i<<" d="<<d<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
1450 <<" dtrc="<<dtrc<<" dlum="<<dlum<<" -> cor="<<corr<<endl;
1451 correction.push_back(corr);
1452 }
1453 intercor = new InterpFunc(loscom2zred_min_,loscom2zred_max_,correction);
1454 }
1455
1456 double dx2=0., dy2=0., dz2=0.;
[3199]1457 for(long i=0;i<Nx_;i++) {
[3271]1458 dx2 = DXcom(i); dx2 *= dx2;
[3199]1459 for(long j=0;j<Ny_;j++) {
[3271]1460 dy2 = DYcom(j); dy2 *= dy2;
[3267]1461 for(long l=0;l<Nz_;l++) {
[3271]1462 double corr = 1.;
1463 if(with_evol) {
1464 dz2 = DZcom(l); dz2 *= dz2; dz2 = sqrt(dx2+dy2+dz2);
1465 corr = (*intercor)(dz2);
1466 }
[3267]1467 int_8 ip = IndexR(i,j,l);
[3271]1468 data_[ip] += snoise*corr*NorRand();
[3199]1469 }
1470 }
1471 }
1472
[3271]1473 if(intercor!=NULL) delete intercor;
[3199]1474}
[3325]1475
1476} // Fin namespace SOPHYA
Note: See TracBrowser for help on using the repository browser.