source: Sophya/trunk/Cosmo/SimLSS/genefluct3d.h@ 3805

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

Refonte totale de la structure des classes InitialSpectrum, TransfertFunction, GrowthFactor, PkSpectrum et classes tabulate pour lecture des fichiers CAMB, cmv 24/07/2010

File size: 7.6 KB
RevLine 
[3115]1#ifndef GENEFLUCT3D_SEEN
2#define GENEFLUCT3D_SEEN
3
4#include "machdefs.h"
[3289]5#include <math.h>
[3115]6#include "genericfunc.h"
7#include "tarray.h"
[3141]8#include "histerr.h"
9#include "hist2err.h"
[3115]10#include "perandom.h"
11
[3141]12#include "FFTW/fftw3.h"
13#include "FitsIO/fitsio.h"
14
[3115]15#include <vector>
[3619]16#include <algorithm>
[3115]17
[3157]18#include "cosmocalc.h"
[3115]19#include "pkspectrum.h"
[3768]20#include "schechter.h"
[3115]21
[3518]22#define WITH_FFTW_THREAD
[3115]23
[3770]24//NE PAS DECOMMENTER, UTILISEZ LA MAKEFILE ("g++ -c -DGEN3D_FLOAT ...")
[3518]25#if defined(GEN3D_FLOAT)
26#define GEN3D_TYPE r_4
27#define GEN3D_FFTW_PLAN fftwf_plan
28#define GEN3D_FFTW_COMPLEX fftwf_complex
29#else
30#define GEN3D_TYPE r_8
31#define GEN3D_FFTW_PLAN fftw_plan
32#define GEN3D_FFTW_COMPLEX fftw_complex
33#endif
34
[3115]35namespace SOPHYA {
36
37//-----------------------------------------------------------------------------------
38class GeneFluct3D {
39public:
[3349]40 GeneFluct3D(long nx,long ny,long nz,double dx,double dy,double dz,unsigned short nthread=0,int lp=0); // Mpc
41 GeneFluct3D(unsigned short nthread=0);
[3115]42 virtual ~GeneFluct3D(void);
43
[3271]44 // Distance los comobile a l'observateur
[3157]45 void SetObservator(double redshref=0.,double kredshref=0.);
[3271]46 inline double DXcom(long i) {return i*Dx_ - xobs_[0];}
47 inline double DYcom(long j) {return j*Dy_ - xobs_[1];}
48 inline double DZcom(long k) {return k*Dz_ - xobs_[2];}
49 inline double Dcom(long i,long j,long k) {
50 double dx=DXcom(i), dy=DYcom(j), dz=DZcom(k);
51 return sqrt(dx*dx+dy*dy+dz*dz);
52 }
[3157]53 void SetCosmology(CosmoCalc& cosmo);
54 void SetGrowthFactor(GrowthFactor& growth);
[3199]55 long LosComRedshift(double zinc=0.001,long npoints=-1);
[3115]56
[3518]57 TArray< complex<GEN3D_TYPE> >& GetComplexArray(void) {return T_;}
58 GEN3D_FFTW_COMPLEX * GetComplexPointer(void) {return fdata_;}
59 TArray<GEN3D_TYPE>& GetRealArray(void) {return R_;}
60 GEN3D_TYPE* GetRealPointer(void) {return data_;}
[3141]61
62 // Pour adressage data_[ip]
63 inline int_8 IndexR(long i,long j,long k) {return (int_8)(k+NTz_*(j+Ny_*i));}
64 // Pour adressage fdata_[ip][0-1]
65 inline int_8 IndexC(long i,long j,long k) {return (int_8)(k+NCz_*(j+Ny_*i));}
[3330]66 // On peut aussi adresser:
67 // TArray< complex<r_8> >& pk = gf3d.GetComplexArray();
68 // pk(k,j,i) avec k=[0,NCz_[ j=[0,Ny_[ i=[0,Nx_[
69 // pk[IndexC(i,j,k)]
70 // TArray<r_8>& rgen = gf3d.GetRealArray();
71 // rgen(k,j,i) avec k=[0,NTz_[ j=[0,Ny_[ i=[0,Nx_[
72 // mais seul k=[0,Nz_[ est utile
73 // rgen[IndexR(i,j,k)]
74 // ATTENTION: TArray adresse en memoire a l'envers du C !
75 // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
[3141]76
77 vector<long> GetNpix(void) {return N_;}
78 int_8 NPix(void) {return NRtot_;}
[3781]79 int_8 NCoeff(void) {return NFcoef_;}
[3768]80 long Nx(void) {return Nx_;}
81 long Ny(void) {return Ny_;}
82 long Nz(void) {return Nz_;}
83 double Dx(void) {return Dx_;}
84 double Dy(void) {return Dy_;}
85 double Dz(void) {return Dz_;}
[3141]86
[3290]87 // Return |K_i| module relative to pixel indices
[3770]88 inline r_8 AbsKx(long i) {long ii=(i>Nx_/2)? Nx_-i :i; return ii*Dkx_;}
89 inline r_8 AbsKy(long j) {long jj=(j>Ny_/2)? Ny_-j :j; return jj*Dky_;}
90 inline r_8 AbsKz(long l) {return l*Dkz_;}
91 inline r_8 Kx(long i) {long ii=(i>Nx_/2)? i-Nx_ :i; return ii*Dkx_;}
92 inline r_8 Ky(long j) {long jj=(j>Ny_/2)? j-Ny_ :j; return jj*Dky_;}
[3290]93 inline r_8 Kz(long l) {return l*Dkz_;}
94
[3141]95 vector<r_8> GetDinc(void) {return D_;}
96 double GetDVol(void) {return dVol_;}
[3115]97 double GetVol(void) {return Vol_;}
[3141]98
99 vector<r_8> GetKinc(void) {return Dk_;}
100 vector<r_8> GetKnyq(void) {return Knyq_;}
[3115]101 double GetKmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_+Knyqz_*Knyqz_);}
[3141]102 double GetKTmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_);}
103 double GetKincMin(void)
104 {vector<r_8>::const_iterator it = min_element(Dk_.begin(), Dk_.end()); return *it;}
[3289]105 double GetKincMax(void)
106 {vector<r_8>::const_iterator it = max_element(Dk_.begin(), Dk_.end()); return *it;}
[3141]107 double GetKTincMin(void) {return min(Dk_[0],Dk_[1]);}
[3289]108 double GetKTincMax(void) {return max(Dk_[0],Dk_[1]);}
[3115]109
[3768]110 double Zref(void) {return redsh_ref_;}
111 double ZrefPk(void) {return compute_pk_redsh_ref_;}
112 double Dref(void) {return growth_ref_;}
113 double DrefPk(void)
114 {return (growth_!=NULL && compute_pk_redsh_ref_>=0.) ? (*growth_)(compute_pk_redsh_ref_): -999.;}
[3770]115 double Href(void) {return h_ref_;}
[3768]116
[3805]117 void ComputeFourier0(PkSpectrum& pk_at_z);
118 void ComputeFourier(PkSpectrum& pk_at_z);
[3115]119 void FilterByPixel(void);
[3800]120 void ToRedshiftSpace(void);
[3770]121 void ToVelLoS(void);
[3141]122
[3115]123 void ComputeReal(void);
[3331]124 void ApplyGrowthFactor(int type_evol=1);
[3800]125 void ApplyRedshiftSpaceGrowthFactor(void);
126 void ApplyVelLosGrowthFactor(int type_evol=1);
[3141]127
[3115]128 void ReComputeFourier(void);
129
[3141]130 int ComputeSpectrum(HistoErr& herr);
131 int ComputeSpectrum2D(Histo2DErr& herr);
[3330]132 int ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor);
133 int ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor);
[3141]134
[3134]135 int_8 VarianceFrReal(double R,double& var);
[3261]136 int_8 MeanSigma2(double& rm,double& rs2,double vmin=1.,double vmax=-1.
137 ,bool useout=false,double vout=0.);
[3320]138 int_8 MinMax(double& xmin,double& xmax,double vmin=1.,double vmax=-1.);
[3134]139 int_8 NumberOfBad(double vmin=-1.e+150,double vmax=1.e+150);
140 int_8 SetToVal(double vmin, double vmax,double val0=0.);
[3283]141 void ScaleOffset(double scalecube=1.,double offsetcube=0.);
[3115]142
143 void TurnFluct2Mass(void);
[3358]144 double TurnFluct2MeanNumber(double val_by_mpc3);
[3115]145 double ApplyPoisson(void);
146 double TurnNGal2Mass(FunRan& massdist,bool axeslog=false);
[3320]147 double TurnNGal2MassQuick(SchechterMassDist& schmdist);
[3115]148 double TurnMass2Flux(void);
[3349]149 //void AddAGN(double lfjy,double lsigma,double powlaw=0.);
[3331]150 void AddNoise2Real(double snoise,int type_evol=0);
[3115]151
[3141]152 void WriteFits(string cfname,int bitpix=FLOAT_IMG);
153 void ReadFits(string cfname);
154
155 void WritePPF(string cfname,bool write_real=true);
156 void ReadPPF(string cfname);
[3281]157 void WriteSlicePPF(string cfname);
[3524]158 void NTupleCheck(POutPersist &pos,string ntname,unsigned long nent);
[3141]159
[3150]160 void SetPrtLevel(int lp=0) {lp_ = lp;}
[3115]161 void Print(void);
162
[3199]163//-------------------------------------------------------------------
164
[3115]165protected:
[3349]166 void init_default(void);
[3141]167 void setsize(long nx,long ny,long nz,double dx,double dy,double dz);
168 void setalloc(void);
169 void setpointers(bool from_real);
[3154]170 void init_fftw(void);
[3349]171 void delete_fftw(void);
[3129]172 long manage_coefficients(void);
[3115]173 double compute_power_carte(void);
[3141]174 void check_array_alloc(void);
[3662]175 inline double pixelfilter(double x) // ATTENTION: seulement pour x>0
[3120]176 {return (x<0.025) ? 1.-x*x/6.*(1.-x*x/20.): sin(x)/x;}
[3115]177
[3154]178 // valeurs dans l'espace reel
[3141]179 long Nx_,Ny_,Nz_; vector<long> N_;
[3129]180 long NCz_,NTz_;
[3781]181 int_8 NRtot_,NFcoef_;
[3141]182
183 double Dx_,Dy_,Dz_; vector<double> D_;
184
[3154]185 // valeurs dans l'espace des K
[3141]186 double Dkx_,Dky_,Dkz_; vector<double> Dk_;
187 double Knyqx_,Knyqy_,Knyqz_; vector<double> Knyq_;
188 double Dk3_;
[3115]189 double dVol_, Vol_;
190
[3154]191 // la gestion de la FFT
[3518]192 bool is_set_fft_plan;
193 GEN3D_FFTW_PLAN pf_,pb_;
[3115]194 unsigned short nthread_;
[3150]195 int lp_;
[3115]196
[3154]197 // le stockage du Cube de donnees et les pointeurs
[3141]198 bool array_allocated_; // true if array has been allocated
[3524]199 unsigned short array_type; // 0=empty, 1=real, 2=complex
[3518]200 TArray< complex<GEN3D_TYPE> > T_;
201 GEN3D_FFTW_COMPLEX *fdata_;
202 TArray<GEN3D_TYPE> R_;
203 GEN3D_TYPE *data_;
[3154]204
205 // l'observateur
[3157]206 CosmoCalc *cosmo_;
207 GrowthFactor *growth_;
[3768]208 double good_dzinc_;
[3271]209 double redsh_ref_,kredsh_ref_,dred_ref_;
[3768]210 double compute_pk_redsh_ref_;
[3770]211 double h_ref_, loscom_ref_,dtrc_ref_, dlum_ref_, dang_ref_;
[3768]212 double growth_ref_, dsdz_growth_ref_;
[3271]213 double nu_ref_, dnu_ref_ ;
[3157]214 double xobs_[3];
[3768]215
[3271]216 double loscom_min_, loscom_max_;
[3157]217 vector<double> zred_, loscom_;
[3199]218 double loscom2zred_min_, loscom2zred_max_;
219 vector<double> loscom2zred_;
220
[3768]221 double zredmin_dpsd_, zredmax_dpsd_;
222 vector<double> dpsdfrzred_;
223
[3115]224};
225
[3325]226} // Fin du namespace SOPHYA
[3115]227
228#endif
Note: See TracBrowser for help on using the repository browser.