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

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

etude des cubes generes par le GSM, cmv 28/06/2010

File size: 54.1 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;
423 double v = -cosmo_->H(z) * growth_->DsDz(z,good_dzinc_) / (*growth_)(z);
424 dpsdfrzred_.push_back(v);
425 if(lp_ && (i%nmod==0 || i==nptd-1)) cout<<" z="<<z<<" D'/D="<<v<<" km/s/Mpc"<<endl;
426 }
427
[3199]428 return zred_.size();
[3157]429}
430
[3115]431//-------------------------------------------------------
[3141]432void GeneFluct3D::WriteFits(string cfname,int bitpix)
433{
[3155]434 cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
[3141]435 try {
436 FitsImg3DWriter fwrt(cfname.c_str(),bitpix,5);
[3790]437 fwrt.WriteKey("NX",Nx_," axe transverse 1"); // NAXIS3
438 fwrt.WriteKey("NY",Ny_," axe transverse 2"); // NAXIS2
439 fwrt.WriteKey("NZ",Nz_," axe longitudinal (redshift)"); // NAXIS1
[3141]440 fwrt.WriteKey("DX",Dx_," Mpc");
441 fwrt.WriteKey("DY",Dy_," Mpc");
442 fwrt.WriteKey("DZ",Dz_," Mpc");
443 fwrt.WriteKey("DKX",Dkx_," Mpc^-1");
444 fwrt.WriteKey("DKY",Dky_," Mpc^-1");
445 fwrt.WriteKey("DKZ",Dkz_," Mpc^-1");
[3768]446 fwrt.WriteKey("ZREFPK",compute_pk_redsh_ref_," Pk computed redshift");
[3271]447 fwrt.WriteKey("ZREF",redsh_ref_," reference redshift");
448 fwrt.WriteKey("KZREF",kredsh_ref_," reference redshift on z axe");
[3770]449 fwrt.WriteKey("HREF",h_ref_," ref H(z)");
[3768]450 fwrt.WriteKey("Growth",Dref()," growth at reference redshift");
451 fwrt.WriteKey("GrowthPk",DrefPk()," growth at Pk computed redshift");
[3141]452 fwrt.Write(R_);
453 } catch (PThrowable & exc) {
454 cout<<"Exception : "<<(string)typeid(exc).name()
455 <<" - Msg= "<<exc.Msg()<<endl;
456 return;
457 } catch (...) {
458 cout<<" some other exception was caught !"<<endl;
459 return;
460 }
461}
462
463void GeneFluct3D::ReadFits(string cfname)
464{
[3155]465 cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
[3141]466 try {
467 FitsImg3DRead fimg(cfname.c_str(),0,5);
468 fimg.Read(R_);
[3790]469 long nx = fimg.ReadKeyL("NX"); // NAXIS3
470 long ny = fimg.ReadKeyL("NY"); // NAXIS2
471 long nz = fimg.ReadKeyL("NZ"); // NAXIS1
[3141]472 double dx = fimg.ReadKey("DX");
473 double dy = fimg.ReadKey("DY");
474 double dz = fimg.ReadKey("DZ");
[3768]475 double pkzref = fimg.ReadKey("ZREFPK");
[3154]476 double zref = fimg.ReadKey("ZREF");
477 double kzref = fimg.ReadKey("KZREF");
[3141]478 setsize(nx,ny,nz,dx,dy,dz);
479 setpointers(true);
[3154]480 init_fftw();
481 SetObservator(zref,kzref);
[3330]482 array_allocated_ = true;
[3768]483 compute_pk_redsh_ref_ = pkzref;
[3141]484 } catch (PThrowable & exc) {
485 cout<<"Exception : "<<(string)typeid(exc).name()
486 <<" - Msg= "<<exc.Msg()<<endl;
487 return;
488 } catch (...) {
489 cout<<" some other exception was caught !"<<endl;
490 return;
491 }
492}
493
494void GeneFluct3D::WritePPF(string cfname,bool write_real)
495// On ecrit soit le TArray<r_8> ou le TArray<complex <r_8> >
496{
[3155]497 cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
[3141]498 try {
499 R_.Info()["NX"] = (int_8)Nx_;
500 R_.Info()["NY"] = (int_8)Ny_;
501 R_.Info()["NZ"] = (int_8)Nz_;
502 R_.Info()["DX"] = (r_8)Dx_;
503 R_.Info()["DY"] = (r_8)Dy_;
504 R_.Info()["DZ"] = (r_8)Dz_;
[3768]505 R_.Info()["ZREFPK"] = (r_8)compute_pk_redsh_ref_;
[3271]506 R_.Info()["ZREF"] = (r_8)redsh_ref_;
507 R_.Info()["KZREF"] = (r_8)kredsh_ref_;
[3770]508 R_.Info()["HREF"] = (r_8)h_ref_;
[3768]509 R_.Info()["Growth"] = (r_8)Dref();
510 R_.Info()["GrowthPk"] = (r_8)DrefPk();
[3141]511 POutPersist pos(cfname.c_str());
512 if(write_real) pos << PPFNameTag("rgen") << R_;
513 else pos << PPFNameTag("pkgen") << T_;
514 } catch (PThrowable & exc) {
515 cout<<"Exception : "<<(string)typeid(exc).name()
516 <<" - Msg= "<<exc.Msg()<<endl;
517 return;
518 } catch (...) {
519 cout<<" some other exception was caught !"<<endl;
520 return;
521 }
522}
523
524void GeneFluct3D::ReadPPF(string cfname)
525{
[3155]526 cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
[3141]527 try {
528 bool from_real = true;
529 PInPersist pis(cfname.c_str());
530 string name_tag_k = "pkgen";
531 bool found_tag_k = pis.GotoNameTag("pkgen");
532 if(found_tag_k) {
[3262]533 cout<<" ...reading spectrum into TArray<complex <r_8> >"<<endl;
[3141]534 pis >> PPFNameTag("pkgen") >> T_;
535 from_real = false;
536 } else {
537 cout<<" ...reading space into TArray<r_8>"<<endl;
538 pis >> PPFNameTag("rgen") >> R_;
539 }
[3154]540 setpointers(from_real); // a mettre ici pour relire les DVInfo
[3141]541 int_8 nx = R_.Info()["NX"];
542 int_8 ny = R_.Info()["NY"];
543 int_8 nz = R_.Info()["NZ"];
544 r_8 dx = R_.Info()["DX"];
545 r_8 dy = R_.Info()["DY"];
546 r_8 dz = R_.Info()["DZ"];
[3768]547 r_8 pkzref = R_.Info()["ZREFPK"];
[3154]548 r_8 zref = R_.Info()["ZREF"];
549 r_8 kzref = R_.Info()["KZREF"];
[3141]550 setsize(nx,ny,nz,dx,dy,dz);
[3154]551 init_fftw();
552 SetObservator(zref,kzref);
[3330]553 array_allocated_ = true;
[3768]554 compute_pk_redsh_ref_ = pkzref;
[3141]555 } catch (PThrowable & exc) {
556 cout<<"Exception : "<<(string)typeid(exc).name()
557 <<" - Msg= "<<exc.Msg()<<endl;
558 return;
559 } catch (...) {
560 cout<<" some other exception was caught !"<<endl;
561 return;
562 }
563}
564
[3281]565void GeneFluct3D::WriteSlicePPF(string cfname)
566// On ecrit 3 tranches du cube selon chaque axe
567{
[3283]568 cout<<"--- GeneFluct3D::WriteSlicePPF: Writing Cube Slices "<<cfname<<endl;
[3281]569 try {
570
571 POutPersist pos(cfname.c_str());
572 TMatrix<r_4> S;
573 char str[16];
574 long i,j,l;
575
576 // Tranches en Z
577 for(int s=0;s<3;s++) {
578 S.ReSize(Nx_,Ny_);
579 if(s==0) l=0; else if(s==1) l=(Nz_+1)/2; else l=Nz_-1;
[3289]580 sprintf(str,"z%ld",l);
[3281]581 for(i=0;i<Nx_;i++) for(j=0;j<Ny_;j++) S(i,j)=data_[IndexR(i,j,l)];
582 pos<<PPFNameTag(str)<<S; S.RenewObjId();
583 }
584
585 // Tranches en Y
586 for(int s=0;s<3;s++) {
587 S.ReSize(Nz_,Nx_);
588 if(s==0) j=0; else if(s==1) j=(Ny_+1)/2; else j=Ny_-1;
[3289]589 sprintf(str,"y%ld",j);
[3281]590 for(i=0;i<Nx_;i++) for(l=0;l<Nz_;l++) S(l,i)=data_[IndexR(i,j,l)];
591 pos<<PPFNameTag(str)<<S; S.RenewObjId();
592 }
593
594 // Tranches en X
595 for(int s=0;s<3;s++) {
596 S.ReSize(Nz_,Ny_);
597 if(s==0) i=0; else if(s==1) i=(Nx_+1)/2; else i=Nx_-1;
[3289]598 sprintf(str,"x%ld",i);
[3281]599 for(j=0;j<Ny_;j++) for(l=0;l<Nz_;l++) S(l,j)=data_[IndexR(i,j,l)];
600 pos<<PPFNameTag(str)<<S; S.RenewObjId();
601 }
602
603 } catch (PThrowable & exc) {
604 cout<<"Exception : "<<(string)typeid(exc).name()
605 <<" - Msg= "<<exc.Msg()<<endl;
606 return;
607 } catch (...) {
608 cout<<" some other exception was caught !"<<endl;
609 return;
610 }
611}
612
[3141]613//-------------------------------------------------------
[3524]614void GeneFluct3D::NTupleCheck(POutPersist &pos,string ntname,unsigned long nent)
615// Remplit le NTuple "ntname" avec "nent" valeurs du cube (reel ou complex) et l'ecrit dans "pos"
616{
617 if(ntname.size()<=0 || nent==0) return;
618 int nvar = 0;
619 if(array_type==1) nvar = 3;
620 else if(array_type==2) nvar = 4;
621 else return;
[3572]622 const char *vname[4] = {"t","z","re","im"};
[3524]623 float xnt[4];
624 NTuple nt(nvar,vname);
625
626 if(array_type==1) {
627 unsigned long nmod = Nx_*Ny_*Nz_/nent; if(nmod==0) nmod=1;
628 unsigned long n=0;
629 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
630 if(n==nmod) {
631 int_8 ip = IndexR(i,j,l);
632 xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=data_[ip];
633 nt.Fill(xnt);
634 n=0;
635 }
636 n++;
637 }
638 } else {
639 unsigned long nmod = Nx_*Ny_*NCz_/nent; if(nmod==0) nmod=1;
640 unsigned long n=0;
641 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<NCz_;l++) {
642 if(n==nmod) {
643 xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=T_(l,j,i).real(); xnt[3]=T_(l,j,i).imag();
644 nt.Fill(xnt);
645 n=0;
646 }
647 n++;
648 }
649 }
650
651 pos.PutObject(nt,ntname);
652}
653
654//-------------------------------------------------------
[3115]655void GeneFluct3D::Print(void)
656{
[3141]657 cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
[3115]658 cout<<"Space Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<Nz_<<" ("<<NTz_<<") size="
659 <<NRtot_<<endl;
660 cout<<" Resol: dx="<<Dx_<<" dy="<<Dy_<<" dz="<<Dz_<<" Mpc"
661 <<", dVol="<<dVol_<<", Vol="<<Vol_<<" Mpc^3"<<endl;
[3781]662 cout<<"Fourier Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<NCz_<<" ncoeff="<<NFcoef_<<endl;
[3115]663 cout<<" Resol: dkx="<<Dkx_<<" dky="<<Dky_<<" dkz="<<Dkz_<<" Mpc^-1"
664 <<", Dk3="<<Dk3_<<" Mpc^-3"<<endl;
665 cout<<" (2Pi/k: "<<2.*M_PI/Dkx_<<" "<<2.*M_PI/Dky_<<" "<<2.*M_PI/Dkz_<<" Mpc)"<<endl;
666 cout<<" Nyquist: kx="<<Knyqx_<<" ky="<<Knyqy_<<" kz="<<Knyqz_<<" Mpc^-1"
667 <<", Kmax="<<GetKmax()<<" Mpc^-1"<<endl;
668 cout<<" (2Pi/k: "<<2.*M_PI/Knyqx_<<" "<<2.*M_PI/Knyqy_<<" "<<2.*M_PI/Knyqz_<<" Mpc)"<<endl;
[3271]669 cout<<"Redshift "<<redsh_ref_<<" for z axe at k="<<kredsh_ref_<<endl;
[3115]670}
671
672//-------------------------------------------------------
[3768]673void GeneFluct3D::ComputeFourier0(PkSpectrumZ& pk_at_z)
[3115]674// cf ComputeFourier() mais avec autre methode de realisation du spectre
675// (attention on fait une fft pour realiser le spectre)
676{
[3768]677 compute_pk_redsh_ref_ = pk_at_z.GetZ();
[3115]678 // --- realisation d'un tableau de tirage gaussiens
[3768]679 if(lp_>0) cout<<"--- ComputeFourier0 at z="<<compute_pk_redsh_ref_
680 <<": before gaussian filling ---"<<endl;
[3115]681 // On tient compte du pb de normalisation de FFTW3
682 double sntot = sqrt((double)NRtot_);
[3129]683 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]684 int_8 ip = IndexR(i,j,l);
685 data_[ip] = NorRand()/sntot;
[3115]686 }
687
688 // --- realisation d'un tableau de tirage gaussiens
[3155]689 if(lp_>0) cout<<"...before fft real ---"<<endl;
[3518]690 GEN3D_FFTW_EXECUTE(pf_);
[3115]691
692 // --- On remplit avec une realisation
[3157]693 if(lp_>0) cout<<"...before Fourier realization filling"<<endl;
[3518]694 T_(0,0,0) = complex<GEN3D_TYPE>(0.); // on coupe le continue et on l'initialise
[3768]695 long lmod = Nx_/20; if(lmod<1) lmod=1;
[3129]696 for(long i=0;i<Nx_;i++) {
[3770]697 double kx = AbsKx(i); kx *= kx;
698 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]699 for(long j=0;j<Ny_;j++) {
[3770]700 double ky = AbsKy(j); ky *= ky;
[3129]701 for(long l=0;l<NCz_;l++) {
[3770]702 double kz = AbsKz(l); kz *= kz;
[3115]703 if(i==0 && j==0 && l==0) continue; // Suppression du continu
704 double k = sqrt(kx+ky+kz);
705 // cf normalisation: Peacock, Cosmology, formule 16.38 p504
[3141]706 double pk = pk_at_z(k)/Vol_;
[3115]707 // ici pas de "/2" a cause de la remarque ci-dessus
708 T_(l,j,i) *= sqrt(pk);
709 }
710 }
711 }
712
[3524]713 array_type = 2;
714
[3155]715 if(lp_>0) cout<<"...computing power"<<endl;
[3115]716 double p = compute_power_carte();
[3155]717 if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
[3115]718
719}
720
721//-------------------------------------------------------
[3768]722void GeneFluct3D::ComputeFourier(PkSpectrumZ& pk_at_z)
[3141]723// Calcule une realisation du spectre "pk_at_z"
[3115]724// Attention: dans TArray le premier indice varie le + vite
725// Explication normalisation: see Coles & Lucchin, Cosmology, p264-265
726// FFTW3: on note N=Nx*Ny*Nz
727// f --(FFT)--> F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f))
728// sum(f(x_i)^2) = S
729// sum(F(nu_i)^2) = S*N
730// sum(fb(x_i)^2) = S*N^2
731{
732 // --- RaZ du tableau
[3518]733 T_ = complex<GEN3D_TYPE>(0.);
[3768]734 compute_pk_redsh_ref_ = pk_at_z.GetZ();
[3115]735
736 // --- On remplit avec une realisation
[3768]737 if(lp_>0) cout<<"--- ComputeFourier at z="<<compute_pk_redsh_ref_<<" ---"<<endl;
738 long lmod = Nx_/20; if(lmod<1) lmod=1;
[3129]739 for(long i=0;i<Nx_;i++) {
[3770]740 double kx = AbsKx(i); kx *= kx;
741 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]742 for(long j=0;j<Ny_;j++) {
[3770]743 double ky = AbsKy(j); ky *= ky;
[3129]744 for(long l=0;l<NCz_;l++) {
[3770]745 double kz = AbsKz(l); kz *= kz;
[3115]746 if(i==0 && j==0 && l==0) continue; // Suppression du continu
747 double k = sqrt(kx+ky+kz);
748 // cf normalisation: Peacock, Cosmology, formule 16.38 p504
[3141]749 double pk = pk_at_z(k)/Vol_;
[3115]750 // Explication de la division par 2: voir perandom.cc
751 // ou egalement Coles & Lucchin, Cosmology formula 13.7.2 p279
[3617]752 T_(l,j,i) = ComplexGaussianRand(sqrt(pk/2.));
[3115]753 }
754 }
755 }
756
[3524]757 array_type = 2;
[3115]758 manage_coefficients(); // gros effet pour les spectres que l'on utilise !
759
[3155]760 if(lp_>0) cout<<"...computing power"<<endl;
[3115]761 double p = compute_power_carte();
[3155]762 if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
[3115]763
764}
765
[3129]766long GeneFluct3D::manage_coefficients(void)
[3115]767// Take into account the real and complexe conjugate coefficients
768// because we want a realization of a real data in real space
[3773]769// On ecrit que: conj(P(k_x,k_y,k_z)) = P(-k_x,-k_y,-k_z)
[3768]770// avec k_x = i, -k_x = N_x - i etc...
[3115]771{
[3155]772 if(lp_>1) cout<<"...managing coefficients"<<endl;
[3141]773 check_array_alloc();
[3115]774
775 // 1./ Le Continu et Nyquist sont reels
[3773]776 int_8 nreal = 0;
[3129]777 for(long kk=0;kk<2;kk++) {
778 long k=0; // continu
[3115]779 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
[3129]780 for(long jj=0;jj<2;jj++) {
781 long j=0;
[3115]782 if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
[3129]783 for(long ii=0;ii<2;ii++) {
784 long i=0;
[3115]785 if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
[3141]786 int_8 ip = IndexC(i,j,k);
787 //cout<<"i="<<i<<" j="<<j<<" k="<<k<<" = ("<<fdata_[ip][0]<<","<<fdata_[ip][1]<<")"<<endl;
788 fdata_[ip][1] = 0.; fdata_[ip][0] *= M_SQRT2;
[3115]789 nreal++;
790 }
791 }
792 }
[3155]793 if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
[3115]794
795 // 2./ Les elements complexe conjugues (tous dans le plan k=0,Nyquist)
796
797 // a./ les lignes et colonnes du continu et de nyquist
[3773]798 int_8 nconj1 = 0;
[3129]799 for(long kk=0;kk<2;kk++) {
800 long k=0; // continu
[3115]801 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
[3129]802 for(long jj=0;jj<2;jj++) { // selon j
803 long j=0;
[3115]804 if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
[3129]805 for(long i=1;i<(Nx_+1)/2;i++) {
[3141]806 int_8 ip = IndexC(i,j,k);
807 int_8 ip1 = IndexC(Nx_-i,j,k);
808 fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
[3115]809 nconj1++;
810 }
811 }
[3129]812 for(long ii=0;ii<2;ii++) {
813 long i=0;
[3115]814 if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
[3129]815 for(long j=1;j<(Ny_+1)/2;j++) {
[3141]816 int_8 ip = IndexC(i,j,k);
817 int_8 ip1 = IndexC(i,Ny_-j,k);
818 fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
[3115]819 nconj1++;
820 }
821 }
822 }
[3155]823 if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
[3115]824
825 // b./ les lignes et colonnes hors continu et de nyquist
[3773]826 int_8 nconj2 = 0;
[3129]827 for(long kk=0;kk<2;kk++) {
828 long k=0; // continu
[3115]829 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
[3129]830 for(long j=1;j<(Ny_+1)/2;j++) {
[3115]831 if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
[3129]832 for(long i=1;i<Nx_;i++) {
[3115]833 if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
[3141]834 int_8 ip = IndexC(i,j,k);
835 int_8 ip1 = IndexC(Nx_-i,Ny_-j,k);
836 fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
[3115]837 nconj2++;
838 }
839 }
840 }
[3155]841 if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
[3115]842
[3781]843 if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(NFcoef_-nconj1-nconj2)-nreal<<endl;
[3115]844
845 return nreal+nconj1+nconj2;
846}
847
848double GeneFluct3D::compute_power_carte(void)
849// Calcul la puissance de la realisation du spectre Pk
850{
[3141]851 check_array_alloc();
852
[3115]853 double s2 = 0.;
[3129]854 for(long l=0;l<NCz_;l++)
855 for(long j=0;j<Ny_;j++)
856 for(long i=0;i<Nx_;i++) s2 += MODULE2(T_(l,j,i));
[3115]857
858 double s20 = 0.;
[3129]859 for(long j=0;j<Ny_;j++)
860 for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
[3115]861
862 double s2n = 0.;
863 if(Nz_%2==0)
[3129]864 for(long j=0;j<Ny_;j++)
865 for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
[3115]866
867 return 2.*s2 -s20 -s2n;
868}
869
870//-------------------------------------------------------------------
871void GeneFluct3D::FilterByPixel(void)
872// Filtrage par la fonction fenetre du pixel (parallelepipede)
[3120]873// TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
[3115]874// e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
[3120]875// = 2/(k_x*dx) * sin(k_x*dx/2) * (idem y) * (idem z)
876// Gestion divergence en 0: sin(y)/y = 1 - y^2/6*(1-y^2/20)
877// avec y = k_x*dx/2
[3115]878{
[3155]879 if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
[3141]880 check_array_alloc();
881
[3129]882 for(long i=0;i<Nx_;i++) {
[3770]883 double kx = AbsKx(i) *Dx_/2;
[3330]884 double pk_x = pixelfilter(kx);
[3129]885 for(long j=0;j<Ny_;j++) {
[3770]886 double ky = AbsKy(j) *Dy_/2;
[3330]887 double pk_y = pixelfilter(ky);
[3129]888 for(long l=0;l<NCz_;l++) {
[3770]889 double kz = AbsKz(l) *Dz_/2;
[3141]890 double pk_z = pixelfilter(kz);
891 T_(l,j,i) *= pk_x*pk_y*pk_z;
[3115]892 }
893 }
894 }
895
896}
897
898//-------------------------------------------------------------------
[3770]899void GeneFluct3D::ToVelLoS(void)
[3773]900// Le spectre Pk doit etre (dRho/Rho)(k)
[3768]901{
902 double zpk = compute_pk_redsh_ref_;
903 double dpsd = -cosmo_->H(zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
[3770]904 if(lp_>0) cout<<"--- ToVelLoS --- at z="<<zpk<<", D'/D="<<dpsd
[3768]905 <<" (km/s)/Mpc (comp. for dz="<<good_dzinc_<<")"<<endl;
906 check_array_alloc();
907
908 for(long i=0;i<Nx_;i++) {
909 double kx = Kx(i);
910 for(long j=0;j<Ny_;j++) {
911 double ky = Ky(j);
912 double kt2 = kx*kx + ky*ky;
913 for(long l=0;l<NCz_;l++) {
914 double kz = Kz(l);
915 double k2 = kt2 + kz*kz;
[3773]916 if(l==0) k2 = 0.; else k2 = dpsd*kz/k2;
917 T_(l,j,i) *= complex<double>(0.,k2);
[3768]918 }
919 }
920 }
[3773]921
[3768]922}
923
924//-------------------------------------------------------------------
[3331]925void GeneFluct3D::ApplyGrowthFactor(int type_evol)
[3157]926// Apply Growth to real space
927// Using the correspondance between redshift and los comoving distance
928// describe in vector "zred_" "loscom_"
[3516]929// type_evol = 1 : evolution avec la distance a l'observateur
930// 2 : evolution avec la distance du plan Z
[3331]931// (tous les pixels d'un plan Z sont mis au meme redshift z que celui du milieu)
[3157]932{
[3331]933 if(lp_>0) cout<<"--- ApplyGrowthFactor: evol="<<type_evol<<endl;
[3157]934 check_array_alloc();
935
936 if(growth_ == NULL) {
[3572]937 const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
[3199]938 cout<<bla<<endl; throw ParmError(bla);
[3157]939 }
[3331]940 if(type_evol<1 || type_evol>2) {
[3572]941 const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
[3331]942 cout<<bla<<endl; throw ParmError(bla);
943 }
[3157]944
[3199]945 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
[3157]946
947 for(long i=0;i<Nx_;i++) {
[3331]948 double dx2 = DXcom(i); dx2 *= dx2;
[3157]949 for(long j=0;j<Ny_;j++) {
[3331]950 double dy2 = DYcom(j); dy2 *= dy2;
[3157]951 for(long l=0;l<Nz_;l++) {
[3331]952 double dz = DZcom(l);
953 if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
954 else dz = fabs(dz); // tous les plans Z au meme redshift
[3768]955 double z = interpinv(dz); // interpolation par morceau
956 double dzgr = (*growth_)(z);
[3157]957 int_8 ip = IndexR(i,j,l);
958 data_[ip] *= dzgr;
959 }
960 }
961 }
962
[3768]963}
[3157]964
[3768]965//-------------------------------------------------------------------
966void GeneFluct3D::ApplyDerGrowthFactor(int type_evol)
967// Apply Conformal derivative of Growth to real space for transverse velocity cube
968{
969 if(lp_>0) cout<<"--- ApplyDerGrowthFactor: evol="<<type_evol<<endl;
970 check_array_alloc();
971
972 if(growth_ == NULL) {
973 const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
974 cout<<bla<<endl; throw ParmError(bla);
975 }
976 if(type_evol<1 || type_evol>2) {
977 const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
978 cout<<bla<<endl; throw ParmError(bla);
979 }
980
981 double zpk = compute_pk_redsh_ref_;
982 double dpsd_orig = - cosmo_->H(zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
983
984 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
985 InterpFunc interdpd(zredmin_dpsd_,zredmax_dpsd_,dpsdfrzred_);
986
987 for(long i=0;i<Nx_;i++) {
988 double dx2 = DXcom(i); dx2 *= dx2;
989 for(long j=0;j<Ny_;j++) {
990 double dy2 = DYcom(j); dy2 *= dy2;
991 for(long l=0;l<Nz_;l++) {
992 double dz = DZcom(l);
993 if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
994 else dz = fabs(dz); // tous les plans Z au meme redshift
995 double z = interpinv(dz); // interpolation par morceau
996 double dpsd = interdpd(z);
997 int_8 ip = IndexR(i,j,l);
998 data_[ip] *= dpsd / dpsd_orig;
999 }
1000 }
1001 }
1002
[3157]1003}
1004
1005//-------------------------------------------------------------------
[3115]1006void GeneFluct3D::ComputeReal(void)
1007// Calcule une realisation dans l'espace reel
1008{
[3155]1009 if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
[3141]1010 check_array_alloc();
[3115]1011
1012 // On fait la FFT
[3518]1013 GEN3D_FFTW_EXECUTE(pb_);
[3524]1014 array_type = 1;
[3115]1015}
1016
1017//-------------------------------------------------------------------
1018void GeneFluct3D::ReComputeFourier(void)
1019{
[3155]1020 if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
[3141]1021 check_array_alloc();
[3115]1022
1023 // On fait la FFT
[3518]1024 GEN3D_FFTW_EXECUTE(pf_);
[3524]1025 array_type = 2;
1026
[3115]1027 // On corrige du pb de la normalisation de FFTW3
[3594]1028 complex<r_8> v((r_8)NRtot_,0.);
[3129]1029 for(long i=0;i<Nx_;i++)
1030 for(long j=0;j<Ny_;j++)
[3594]1031 for(long l=0;l<NCz_;l++) T_(l,j,i) /= v;
[3115]1032}
1033
1034//-------------------------------------------------------------------
[3141]1035int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
1036// Compute spectrum from "T" and fill HistoErr "herr"
[3115]1037// T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
1038// cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
1039{
[3155]1040 if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
[3141]1041 check_array_alloc();
[3115]1042
[3141]1043 if(herr.NBins()<0) return -1;
1044 herr.Zero();
[3115]1045
1046 // Attention a l'ordre
[3129]1047 for(long i=0;i<Nx_;i++) {
[3770]1048 double kx = AbsKx(i); kx *= kx;
[3129]1049 for(long j=0;j<Ny_;j++) {
[3770]1050 double ky = AbsKy(j); ky *= ky;
[3129]1051 for(long l=0;l<NCz_;l++) {
[3770]1052 double kz = AbsKz(l);
[3330]1053 double k = sqrt(kx+ky+kz*kz);
[3115]1054 double pk = MODULE2(T_(l,j,i));
[3141]1055 herr.Add(k,pk);
[3115]1056 }
1057 }
1058 }
[3150]1059 herr.ToVariance();
[3115]1060
1061 // renormalize to directly compare to original spectrum
1062 double norm = Vol_;
[3141]1063 herr *= norm;
[3115]1064
1065 return 0;
1066}
1067
[3141]1068int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
1069{
[3155]1070 if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
[3141]1071 check_array_alloc();
1072
1073 if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
1074 herr.Zero();
1075
1076 // Attention a l'ordre
1077 for(long i=0;i<Nx_;i++) {
[3770]1078 double kx = AbsKx(i); kx *= kx;
[3141]1079 for(long j=0;j<Ny_;j++) {
[3770]1080 double ky = AbsKy(j); ky *= ky;
[3141]1081 double kt = sqrt(kx+ky);
1082 for(long l=0;l<NCz_;l++) {
[3770]1083 double kz = AbsKz(l);
[3141]1084 double pk = MODULE2(T_(l,j,i));
1085 herr.Add(kt,kz,pk);
1086 }
1087 }
1088 }
[3150]1089 herr.ToVariance();
[3141]1090
1091 // renormalize to directly compare to original spectrum
1092 double norm = Vol_;
1093 herr *= norm;
1094
1095 return 0;
1096}
1097
[3330]1098//-------------------------------------------------------------------
1099int GeneFluct3D::ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor)
1100// Compute spectrum from "T" and fill HistoErr "herr"
1101// AVEC la soustraction du niveau de bruit et la correction par filterpixel()
1102// Si on ne fait pas ca, alors on obtient un spectre non-isotrope!
1103//
1104// T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
1105// cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
1106{
1107 if(lp_>0) cout<<"--- ComputeSpectrum: sigma="<<sigma<<endl;
1108 check_array_alloc();
1109
1110 if(sigma<=0.) sigma = 0.;
1111 double sigma2 = sigma*sigma / (double)NRtot_;
1112
1113 if(herr.NBins()<0) return -1;
1114 herr.Zero();
1115
1116 TVector<r_8> vfz(NCz_);
1117 if(pixcor) // kz = l*Dkz_
1118 for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
1119
1120 // Attention a l'ordre
1121 for(long i=0;i<Nx_;i++) {
[3770]1122 double kx = AbsKx(i);
[3330]1123 double fx = (pixcor) ? pixelfilter(kx*Dx_/2): 1.;
1124 kx *= kx; fx *= fx;
1125 for(long j=0;j<Ny_;j++) {
[3770]1126 double ky = AbsKy(j);
[3330]1127 double fy = (pixcor) ? pixelfilter(ky*Dy_/2): 1.;
1128 ky *= ky; fy *= fy;
1129 for(long l=0;l<NCz_;l++) {
[3770]1130 double kz = AbsKz(l);
[3330]1131 double k = sqrt(kx+ky+kz*kz);
1132 double pk = MODULE2(T_(l,j,i)) - sigma2;
1133 double fz = (pixcor) ? vfz(l): 1.;
1134 double f = fx*fy*fz;
1135 if(f>0.) herr.Add(k,pk/f);
1136 }
1137 }
1138 }
1139 herr.ToVariance();
[3351]1140 for(int i=0;i<herr.NBins();i++) herr(i) += sigma2;
[3330]1141
1142 // renormalize to directly compare to original spectrum
1143 double norm = Vol_;
1144 herr *= norm;
1145
1146 return 0;
1147}
1148
1149int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor)
1150// AVEC la soustraction du niveau de bruit et la correction par filterpixel()
1151{
1152 if(lp_>0) cout<<"--- ComputeSpectrum2D: sigma="<<sigma<<endl;
1153 check_array_alloc();
1154
1155 if(sigma<=0.) sigma = 0.;
1156 double sigma2 = sigma*sigma / (double)NRtot_;
1157
1158 if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
1159 herr.Zero();
1160
1161 TVector<r_8> vfz(NCz_);
1162 if(pixcor) // kz = l*Dkz_
1163 for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
1164
1165 // Attention a l'ordre
1166 for(long i=0;i<Nx_;i++) {
[3770]1167 double kx = AbsKx(i);
[3330]1168 double fx = (pixcor) ? pixelfilter(kx*Dx_/2) : 1.;
1169 kx *= kx; fx *= fx;
1170 for(long j=0;j<Ny_;j++) {
[3770]1171 double ky = AbsKy(j);
[3330]1172 double fy = (pixcor) ? pixelfilter(ky*Dy_/2) : 1.;
1173 ky *= ky; fy *= fy;
1174 double kt = sqrt(kx+ky);
1175 for(long l=0;l<NCz_;l++) {
[3770]1176 double kz = AbsKz(l);
[3330]1177 double pk = MODULE2(T_(l,j,i)) - sigma2;
1178 double fz = (pixcor) ? vfz(l): 1.;
1179 double f = fx*fy*fz;
1180 if(f>0.) herr.Add(kt,kz,pk/f);
1181 }
1182 }
1183 }
1184 herr.ToVariance();
[3351]1185 for(int i=0;i<herr.NBinX();i++)
1186 for(int j=0;j<herr.NBinY();j++) herr(i,j) += sigma2;
[3330]1187
1188 // renormalize to directly compare to original spectrum
1189 double norm = Vol_;
1190 herr *= norm;
1191
1192 return 0;
1193}
1194
[3115]1195//-------------------------------------------------------
[3134]1196int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
[3115]1197// Recompute MASS variance in spherical top-hat (rayon=R)
[3353]1198// Par definition: SigmaR^2 = <(M-<M>)^2>/<M>^2
1199// ou M = masse dans sphere de rayon R
[3354]1200// --- ATTENTION: la variance calculee a une tres grande dispersion
1201// (surtout si le volume du cube est petit). Pour verifier
1202// que le sigmaR calcule par cette methode est en accord avec
1203// le sigmaR en input, il faut faire plusieurs simulations (~100)
1204// et regarder la moyenne des sigmaR reconstruits
[3115]1205{
[3262]1206 if(lp_>0) cout<<"--- VarianceFrReal R="<<R<<endl;
[3141]1207 check_array_alloc();
1208
[3353]1209 long dnx = long(R/Dx_)+1; if(dnx<=0) dnx = 1;
1210 long dny = long(R/Dy_)+1; if(dny<=0) dny = 1;
1211 long dnz = long(R/Dz_)+1; if(dnz<=0) dnz = 1;
[3155]1212 if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
[3115]1213
[3353]1214 double sum=0., sum2=0., sn=0., r2 = R*R;
1215 int_8 nsum=0;
[3115]1216
[3353]1217 for(long i=dnx;i<Nx_-dnx;i+=2*dnx) {
1218 for(long j=dny;j<Ny_-dny;j+=2*dny) {
1219 for(long l=dnz;l<Nz_-dnz;l+=2*dnz) {
1220 double m=0.; int_8 n=0;
[3129]1221 for(long ii=i-dnx;ii<=i+dnx;ii++) {
[3115]1222 double x = (ii-i)*Dx_; x *= x;
[3129]1223 for(long jj=j-dny;jj<=j+dny;jj++) {
[3115]1224 double y = (jj-j)*Dy_; y *= y;
[3129]1225 for(long ll=l-dnz;ll<=l+dnz;ll++) {
[3115]1226 double z = (ll-l)*Dz_; z *= z;
1227 if(x+y+z>r2) continue;
[3141]1228 int_8 ip = IndexR(ii,jj,ll);
[3353]1229 m += 1.+data_[ip]; // 1+drho/rho
[3115]1230 n++;
1231 }
1232 }
1233 }
[3353]1234 if(n>0) {sum += m; sum2 += m*m; nsum++; sn += n;}
1235 //cout<<i<<","<<j<<","<<l<<" n="<<n<<" m="<<m<<" sum="<<sum<<" sum2="<<sum2<<endl;
[3115]1236 }
1237 }
1238 }
1239
1240 if(nsum<=1) {var=0.; return nsum;}
1241 sum /= nsum;
1242 sum2 = sum2/nsum - sum*sum;
[3353]1243 sn /= nsum;
1244 if(lp_>0) cout<<"...<n>="<<sn<<", nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
1245 var = sum2/(sum*sum); // <dM^2>/<M>^2
1246 if(lp_>0) cout<<"...sigmaR^2 = <(M-<M>)^2>/<M>^2 = "<<var
1247 <<" -> sigmaR = "<<sqrt(var)<<endl;
[3115]1248
1249 return nsum;
1250}
1251
1252//-------------------------------------------------------
[3134]1253int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
[3115]1254// number of pixels outside of ]vmin,vmax[ extremites exclues
1255// -> vmin and vmax are considered as bad
1256{
[3141]1257 check_array_alloc();
[3115]1258
[3134]1259 int_8 nbad = 0;
[3129]1260 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1261 int_8 ip = IndexR(i,j,l);
1262 double v = data_[ip];
[3115]1263 if(v<=vmin || v>=vmax) nbad++;
1264 }
1265
[3358]1266 if(lp_>0) cout<<"--- NumberOfBad "<<nbad<<" px out of ]"<<vmin<<","<<vmax
1267 <<"[ i.e. frac="<<nbad/(double)NRtot_<<endl;
[3115]1268 return nbad;
1269}
1270
[3320]1271int_8 GeneFluct3D::MinMax(double& xmin,double& xmax,double vmin,double vmax)
1272// Calcul des valeurs xmin et xmax dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
1273{
1274 bool tstval = (vmax>vmin)? true: false;
1275 if(lp_>0) {
1276 cout<<"--- MinMax";
1277 if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
1278 cout<<endl;
1279 }
1280 check_array_alloc();
1281
1282 int_8 n = 0;
1283 xmin = xmax = data_[0];
1284
1285 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1286 int_8 ip = IndexR(i,j,l);
1287 double x = data_[ip];
1288 if(tstval && (x<=vmin || x>=vmax)) continue;
1289 if(x<xmin) xmin = x;
1290 if(x>xmax) xmax = x;
1291 n++;
1292 }
1293
1294 if(lp_>0) cout<<" n="<<n<<" min="<<xmin<<" max="<<xmax<<endl;
1295
1296 return n;
1297}
1298
[3261]1299int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax
1300 ,bool useout,double vout)
1301// Calcul de mean,sigma2 dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
1302// useout = false: ne pas utiliser les pixels hors limites pour calculer mean,sigma2
1303// true : utiliser les pixels hors limites pour calculer mean,sigma2
1304// en remplacant leurs valeurs par "vout"
[3115]1305{
[3261]1306 bool tstval = (vmax>vmin)? true: false;
1307 if(lp_>0) {
[3262]1308 cout<<"--- MeanSigma2";
1309 if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
[3261]1310 if(useout) cout<<", useout="<<useout<<" vout="<<vout;
1311 cout<<endl;
1312 }
[3141]1313 check_array_alloc();
[3115]1314
[3134]1315 int_8 n = 0;
[3115]1316 rm = rs2 = 0.;
1317
[3129]1318 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1319 int_8 ip = IndexR(i,j,l);
1320 double v = data_[ip];
[3261]1321 if(tstval) {
1322 if(v<=vmin || v>=vmax) {if(useout) v=vout; else continue;}
1323 }
[3115]1324 rm += v;
1325 rs2 += v*v;
1326 n++;
1327 }
1328
1329 if(n>1) {
1330 rm /= (double)n;
1331 rs2 = rs2/(double)n - rm*rm;
1332 }
1333
[3261]1334 if(lp_>0) cout<<" n="<<n<<" m="<<rm<<" s2="<<rs2<<" s="<<sqrt(fabs(rs2))<<endl;
1335
[3115]1336 return n;
1337}
1338
[3134]1339int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
[3115]1340// set to "val0" if out of range ]vmin,vmax[ extremites exclues
[3261]1341// cad set to "val0" if in [vmin,vmax] -> vmin and vmax are set to val0
[3115]1342{
[3141]1343 check_array_alloc();
[3115]1344
[3134]1345 int_8 nbad = 0;
[3129]1346 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1347 int_8 ip = IndexR(i,j,l);
1348 double v = data_[ip];
1349 if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
[3115]1350 }
1351
[3262]1352 if(lp_>0) cout<<"--- SetToVal "<<nbad<<" px set to="<<val0
1353 <<" because out of range=]"<<vmin<<","<<vmax<<"["<<endl;
[3115]1354 return nbad;
1355}
1356
[3283]1357void GeneFluct3D::ScaleOffset(double scalecube,double offsetcube)
1358// Replace "V" by "scalecube * ( V + offsetcube )"
1359{
[3284]1360 if(lp_>0) cout<<"--- ScaleCube scale="<<scalecube<<" offset="<<offsetcube<<endl;
[3283]1361
1362 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1363 int_8 ip = IndexR(i,j,l);
1364 data_[ip] = scalecube * ( data_[ip] + offsetcube );
1365 }
1366
1367 return;
1368}
1369
[3115]1370//-------------------------------------------------------
1371void GeneFluct3D::TurnFluct2Mass(void)
1372// d_rho/rho -> Mass (add one!)
1373{
[3155]1374 if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
[3141]1375 check_array_alloc();
1376
[3115]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 data_[ip] += 1.;
[3115]1381 }
1382}
1383
[3358]1384double GeneFluct3D::TurnFluct2MeanNumber(double val_by_mpc3)
[3365]1385// ATTENTION: la gestion des pixels<0 proposee ici induit une perte de variance
1386// sur la carte, le spectre Pk reconstruit sera plus faible!
1387// L'effet sera d'autant plus grand que le nombre de pixels<0 sera grand.
[3329]1388{
[3358]1389 if(lp_>0) cout<<"--- TurnFluct2MeanNumber : "<<val_by_mpc3<<" quantity (gal or mass)/Mpc^3"<<endl;
[3329]1390
[3358]1391 // First convert dRho/Rho into 1+dRho/Rho
1392 int_8 nball = 0; double sumall = 0., sumall2 = 0.;
[3329]1393 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1394 int_8 ip = IndexR(i,j,l);
[3358]1395 data_[ip] += 1.;
1396 nball++; sumall += data_[ip]; sumall2 += data_[ip]*data_[ip];
[3329]1397 }
[3358]1398 if(nball>2) {
1399 sumall /= (double)nball;
1400 sumall2 = sumall2/(double)nball - sumall*sumall;
1401 if(lp_>0) cout<<"1+dRho/Rho: mean="<<sumall<<" variance="<<sumall2
1402 <<" -> "<<sqrt(fabs(sumall2))<<endl;
[3329]1403 }
1404
[3358]1405 // Find contribution for positive pixels
1406 int_8 nbpos = 0; double sumpos = 0. , sumpos2 = 0.;
1407 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1408 int_8 ip = IndexR(i,j,l);
1409 double v = data_[ip];
1410 if(data_[ip]>0.) {nbpos++; sumpos += v; sumpos2 += v*v;}
1411 }
1412 if(nbpos<1) {
1413 cout<<"TurnFluct2MeanNumber_Error: nbpos<1"<<endl;
1414 throw RangeCheckError("TurnFluct2MeanNumber_Error: nbpos<1");
1415 }
1416 sumpos2 = sumpos2/nball - sumpos*sumpos/(nball*nball);
1417 if(lp_>0)
1418 cout<<"1+dRho/Rho with v<0 set to zero: mean="<<sumpos/nball
1419 <<" variance="<<sumpos2<<" -> "<<sqrt(fabs(sumpos2))<<endl;
1420 cout<<"Sum of positive values: sumpos="<<sumpos
1421 <<" (n(v>0) = "<<nbpos<<" frac(v>0)="<<nbpos/(double)NRtot_<<")"<<endl;
[3329]1422
[3358]1423 // - Mettre exactement val_by_mpc3*Vol galaxies (ou Msol) dans notre survey
1424 // - Uniquement dans les pixels de masse >0.
1425 // - Mise a zero des pixels <0
1426 double dn = val_by_mpc3 * Vol_ / sumpos;
1427 if(lp_>0) cout<<"...density move from "
1428 <<val_by_mpc3*dVol_<<" to "<<dn<<" / pixel"<<endl;
1429
[3329]1430 double sum = 0.;
1431 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1432 int_8 ip = IndexR(i,j,l);
1433 if(data_[ip]<=0.) data_[ip] = 0.;
1434 else {
[3349]1435 data_[ip] *= dn;
1436 sum += data_[ip];
[3329]1437 }
1438 }
1439
[3358]1440 if(lp_>0) cout<<"...quantity put into survey "<<sum<<" / "<<val_by_mpc3*Vol_<<endl;
[3329]1441
1442 return sum;
1443}
1444
[3115]1445double GeneFluct3D::ApplyPoisson(void)
1446// do NOT treate negative or nul mass -> let it as it is
1447{
[3155]1448 if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
[3141]1449 check_array_alloc();
1450
[3115]1451 double sum = 0.;
[3129]1452 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1453 int_8 ip = IndexR(i,j,l);
1454 double v = data_[ip];
[3115]1455 if(v>0.) {
[3617]1456 uint_8 dn = PoissonRand(v,10.);
[3141]1457 data_[ip] = (double)dn;
[3115]1458 sum += (double)dn;
1459 }
1460 }
[3155]1461 if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
[3115]1462
1463 return sum;
1464}
1465
1466double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
1467// do NOT treate negative or nul mass -> let it as it is
1468// INPUT:
1469// massdist : distribution de masse (m*dn/dm)
1470// axeslog = false : retourne la masse
1471// = true : retourne le log10(mass)
1472// RETURN la masse totale
1473{
[3155]1474 if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
[3141]1475 check_array_alloc();
1476
[3115]1477 double sum = 0.;
[3129]1478 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
[3141]1479 int_8 ip = IndexR(i,j,l);
1480 double v = data_[ip];
[3115]1481 if(v>0.) {
[3129]1482 long ngal = long(v+0.1);
[3141]1483 data_[ip] = 0.;
[3129]1484 for(long i=0;i<ngal;i++) {
[3115]1485 double m = massdist.RandomInterp(); // massdist.Random();
1486 if(axeslog) m = pow(10.,m);
[3141]1487 data_[ip] += m;
[3115]1488 }
[3141]1489 sum += data_[ip];
[3115]1490 }
1491 }
[3155]1492 if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
[3115]1493
1494 return sum;
1495}
1496
[3320]1497double GeneFluct3D::TurnNGal2MassQuick(SchechterMassDist& schmdist)
1498// idem TurnNGal2Mass mais beaucoup plus rapide
1499{
1500 if(lp_>0) cout<<"--- TurnNGal2MassQuick ---"<<endl;
1501 check_array_alloc();
1502
1503 double sum = 0.;
1504 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
1505 int_8 ip = IndexR(i,j,l);
1506 double v = data_[ip];
1507 if(v>0.) {
1508 long ngal = long(v+0.1);
1509 data_[ip] = schmdist.TirMass(ngal);
1510 sum += data_[ip];
1511 }
1512 }
1513 if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
1514
1515 return sum;
1516}
1517
[3349]1518void GeneFluct3D::AddNoise2Real(double snoise,int type_evol)
1519// add noise to every pixels (meme les <=0 !)
1520// type_evol = 0 : pas d'evolution de la puissance du bruit
1521// 1 : evolution de la puissance du bruit avec la distance a l'observateur
1522// 2 : evolution de la puissance du bruit avec la distance du plan Z
1523// (tous les plans Z sont mis au meme redshift z de leur milieu)
1524{
1525 if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<" evol="<<type_evol<<endl;
1526 check_array_alloc();
1527
1528 if(type_evol<0) type_evol = 0;
1529 if(type_evol>2) {
[3572]1530 const char *bla = "GeneFluct3D::AddNoise2Real_Error: bad type_evol value";
[3349]1531 cout<<bla<<endl; throw ParmError(bla);
1532 }
1533
1534 vector<double> correction;
1535 InterpFunc *intercor = NULL;
1536
1537 if(type_evol>0) {
1538 // Sigma_Noise(en mass) :
1539 // Slim ~ 1/sqrt(DNu) * sqrt(nlobe) en W/m^2Hz
1540 // Flim ~ sqrt(DNu) * sqrt(nlobe) en W/m^2
1541 // Mlim ~ sqrt(DNu) * (Dlum)^2 * sqrt(nlobe) en Msol
1542 // nlobe ~ 1/Dtrcom^2
1543 // Mlim ~ sqrt(DNu) * (Dlum)^2 / Dtrcom
[3662]1544 if(cosmo_ == NULL || redsh_ref_<0. || loscom2zred_.size()<1) {
[3572]1545 const char *bla = "GeneFluct3D::AddNoise2Real_Error: set Observator and Cosmology first";
[3349]1546 cout<<bla<<endl; throw ParmError(bla);
1547 }
1548 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
1549 long nsz = loscom2zred_.size(), nszmod=((nsz>10)? nsz/10: 1);
1550 for(long i=0;i<nsz;i++) {
1551 double d = interpinv.X(i);
1552 double zred = interpinv(d);
1553 double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
1554 double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
1555 double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
1556 double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
1557 double corr = sqrt(dnu/dnu_ref_) * pow(dlum/dlum_ref_,2.) * dtrc_ref_/dtrc;
1558 if(lp_>0 && (i==0 || i==nsz-1 || i%nszmod==0))
1559 cout<<"i="<<i<<" d="<<d<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
1560 <<" dtrc="<<dtrc<<" dlum="<<dlum<<" -> cor="<<corr<<endl;
1561 correction.push_back(corr);
1562 }
1563 intercor = new InterpFunc(loscom2zred_min_,loscom2zred_max_,correction);
1564 }
1565
1566 double corrlim[2] = {1.,1.};
1567 for(long i=0;i<Nx_;i++) {
1568 double dx2 = DXcom(i); dx2 *= dx2;
1569 for(long j=0;j<Ny_;j++) {
1570 double dy2 = DYcom(j); dy2 *= dy2;
1571 for(long l=0;l<Nz_;l++) {
1572 double corr = 1.;
1573 if(type_evol>0) {
1574 double dz = DZcom(l);
1575 if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
1576 else dz = fabs(dz); // tous les plans Z au meme redshift
1577 corr = (*intercor)(dz);
1578 if(corr<corrlim[0]) corrlim[0]=corr; else if(corr>corrlim[1]) corrlim[1]=corr;
1579 }
1580 int_8 ip = IndexR(i,j,l);
1581 data_[ip] += snoise*corr*NorRand();
1582 }
1583 }
1584 }
1585 if(type_evol>0)
1586 cout<<"correction factor range: ["<<corrlim[0]<<","<<corrlim[1]<<"]"<<endl;
1587
1588 if(intercor!=NULL) delete intercor;
1589}
1590
1591} // Fin namespace SOPHYA
1592
1593
1594
1595
1596/*********************************************************************
[3199]1597void GeneFluct3D::AddAGN(double lfjy,double lsigma,double powlaw)
[3196]1598// Add AGN flux into simulation:
1599// --- Procedure:
1600// 1. lancer "cmvdefsurv" avec les parametres du survey
[3199]1601// (au redshift de reference du survey)
[3196]1602// et recuperer l'angle solide "angsol sr" du pixel elementaire
1603// au centre du cube.
1604// 2. lancer "cmvtstagn" pour cet angle solide -> cmvtstagn.ppf
1605// 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne
1606// cad une moyenne <log10(S)> et un sigma "sig"
[3199]1607// Attention: la distribution n'est pas gaussienne les "mean,sigma"
1608// de l'histo ne sont pas vraiment ce que l'on veut
[3196]1609// --- Limitations actuelle du code:
[3271]1610// . les AGN sont supposes evoluer avec la meme loi de puissance pour tout theta,phi
[3199]1611// . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee
1612// . la distribution est approximee a une gaussienne
1613// ... C'est une approximation pour un observateur loin du centre du cube
1614// et pour un cube peu epais / distance observateur
[3196]1615// --- Parametres de la routine:
[3271]1616// llfy : c'est le <log10(S)> du flux depose par les AGN
1617// dans l'angle solide du pixel elementaire de reference du cube
1618// lsigma : c'est le sigma de la distribution des log10(S)
1619// powlaw : c'est la pente de la distribution cad que le flux "lmsol"
[3199]1620// et considere comme le flux a 1.4GHz et qu'on suppose une loi
1621// F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz)
[3196]1622// - Comme on est en echelle log10():
1623// on tire log10(Msol) + X
1624// ou X est une realisation sur une gaussienne de variance "sig^2"
1625// La masse realisee est donc: Msol*10^X
1626// - Pas de probleme de pixel negatif car on a une multiplication!
1627{
[3199]1628 if(lp_>0) cout<<"--- AddAGN: <log10(S Jy)> = "<<lfjy<<" , sigma = "<<lsigma<<endl;
[3196]1629 check_array_alloc();
1630
[3271]1631 if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) {
[3199]1632 char *bla = "GeneFluct3D::AddAGN_Error: set Observator and Cosmology first";
1633 cout<<bla<<endl; throw ParmError(bla);
1634 }
[3196]1635
[3271]1636 // Le flux des AGN en Jy et en mass solaire
1637 double fagnref = pow(10.,lfjy)*(dnu_ref_*1.e9); // Jy.Hz = W/m^2
1638 double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlum_ref_); // Msol
1639 if(lp_>0)
1640 cout<<"Au pixel de ref: fagnref="<<fagnref
1641 <<" Jy.Hz (a 1.4GHz), magnref="<<magnref<<" Msol"<<endl;
[3196]1642
[3199]1643 if(powlaw!=0.) {
[3271]1644 // 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)
1645 magnref *= pow(1/(1.+redsh_ref_),powlaw);
[3199]1646 if(lp_>0) cout<<" powlaw="<<powlaw<<" -> change magnref to "<<magnref<<" Msol"<<endl;
1647 }
1648
1649 // Les infos en fonction de l'indice "l" selon Oz
1650 vector<double> correction;
1651 InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
[3271]1652 long nzmod = ((Nz_>10)?Nz_/10:1);
[3199]1653 for(long l=0;l<Nz_;l++) {
[3271]1654 double z = fabs(DZcom(l));
[3199]1655 double zred = interpinv(z);
[3271]1656 double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
[3199]1657 double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
1658 double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
1659 double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
[3271]1660 // on a: Mass ~ DNu * Dlum^2 / Dtrcom^2
1661 double corr = dnu/dnu_ref_*pow(dtrc_ref_/dtrc*dlum/dlum_ref_,2.);
1662 // 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)
1663 if(powlaw!=0.) corr *= pow((1.+redsh_ref_)/(1.+zred),powlaw);
[3199]1664 correction.push_back(corr);
[3271]1665 if(lp_>0 && (l==0 || l==Nz_-1 || l%nzmod==0)) {
1666 cout<<"l="<<l<<" z="<<z<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
1667 <<" dtrc="<<dtrc<<" dlum="<<dlum
1668 <<" -> cor="<<corr<<endl;
[3199]1669 }
1670 }
1671
1672 double sum=0., sum2=0., nsum=0.;
1673 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) {
1674 double a = lsigma*NorRand();
1675 a = magnref*pow(10.,a);
1676 // On met le meme tirage le long de Oz (indice k)
1677 for(long l=0;l<Nz_;l++) {
1678 int_8 ip = IndexR(i,j,l);
1679 data_[ip] += a*correction[l];
1680 }
1681 sum += a; sum2 += a*a; nsum += 1.;
1682 }
1683
1684 if(lp_>0 && nsum>1.) {
[3196]1685 sum /= nsum;
1686 sum2 = sum2/nsum - sum*sum;
1687 cout<<"...Mean mass="<<sum<<" Msol , s^2="<<sum2<<" s="<<sqrt(fabs(sum2))<<endl;
1688 }
1689
1690}
[3349]1691*********************************************************************/
Note: See TracBrowser for help on using the repository browser.