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

Last change on this file since 4047 was 4047, checked in by cmv, 14 years ago

details..., cmv 12/01/2012

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