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

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

suite de la mise au point pour lecture fichiers CAMB, cmv 24/07/2010

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