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

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

intro choix dabs calcul de beta pour les redshift-distorsions, cmv 10/12/2010

File size: 7.7 KB
Line 
1#ifndef GENEFLUCT3D_SEEN
2#define GENEFLUCT3D_SEEN
3
4#include "machdefs.h"
5#include <math.h>
6#include "genericfunc.h"
7#include "tarray.h"
8#include "histerr.h"
9#include "hist2err.h"
10#include "perandom.h"
11
12#include "FFTW/fftw3.h"
13#include "FitsIO/fitsio.h"
14
15#include <vector>
16#include <algorithm>
17
18#include "cosmocalc.h"
19#include "pkspectrum.h"
20#include "schechter.h"
21
22#define WITH_FFTW_THREAD
23
24//NE PAS DECOMMENTER, UTILISEZ LA MAKEFILE ("g++ -c -DGEN3D_FLOAT ...")
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
35namespace SOPHYA {
36
37//-----------------------------------------------------------------------------------
38class GeneFluct3D {
39public:
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);
42 virtual ~GeneFluct3D(void);
43
44 // Distance los comobile a l'observateur
45 void SetObservator(double redshref=0.,double kredshref=0.);
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 }
53 void SetCosmology(CosmoCalc& cosmo);
54 void SetGrowthFactor(GrowthFactor& growth);
55 long LosComRedshift(double zinc=0.001,long npoints=-1);
56
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_;}
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));}
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]
76
77 vector<long> GetNpix(void) {return N_;}
78 int_8 NPix(void) {return NRtot_;}
79 int_8 NCoeff(void) {return NFcoef_;}
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_;}
86
87 // Return |K_i| module relative to pixel indices
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_;}
93 inline r_8 Kz(long l) {return l*Dkz_;}
94
95 vector<r_8> GetDinc(void) {return D_;}
96 double GetDVol(void) {return dVol_;}
97 double GetVol(void) {return Vol_;}
98
99 vector<r_8> GetKinc(void) {return Dk_;}
100 vector<r_8> GetKnyq(void) {return Knyq_;}
101 double GetKmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_+Knyqz_*Knyqz_);}
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;}
105 double GetKincMax(void)
106 {vector<r_8>::const_iterator it = max_element(Dk_.begin(), Dk_.end()); return *it;}
107 double GetKTincMin(void) {return min(Dk_[0],Dk_[1]);}
108 double GetKTincMax(void) {return max(Dk_[0],Dk_[1]);}
109
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.;}
115 double Href(void) {return h_ref_;}
116
117 void ComputeFourier0(PkSpectrum& pk_at_z);
118 void ComputeFourier(PkSpectrum& pk_at_z);
119 void FilterByPixel(void);
120 void ToRedshiftSpace(int type_beta=0,double fixed_beta=0.);
121 void ToVelLoS(int type_beta=0,double fixed_beta=0.);
122
123 void ComputeReal(void);
124 void ApplyGrowthFactor(int type_evol=1);
125 void ApplyRedshiftSpaceGrowthFactor(void);
126 void ApplyVelLosGrowthFactor(int type_evol=1,int type_beta=0,double fixed_beta=0.);
127
128 void ReComputeFourier(void);
129
130 int ComputeSpectrum(HistoErr& herr);
131 int ComputeSpectrum2D(Histo2DErr& herr);
132 int ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor);
133 int ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor);
134
135 int_8 VarianceFrReal(double R,double& var);
136 int_8 MeanSigma2(double& rm,double& rs2,double vmin=1.,double vmax=-1.
137 ,bool useout=false,double vout=0.);
138 int_8 MinMax(double& xmin,double& xmax,double vmin=1.,double vmax=-1.);
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.);
141 void ScaleOffset(double scalecube=1.,double offsetcube=0.);
142
143 void TurnFluct2Mass(void);
144 double TurnFluct2MeanNumber(double val_by_mpc3);
145 double ApplyPoisson(void);
146 double TurnNGal2Mass(FunRan& massdist,bool axeslog=false);
147 double TurnNGal2MassQuick(SchechterMassDist& schmdist);
148 double TurnMass2Flux(void);
149 //void AddAGN(double lfjy,double lsigma,double powlaw=0.);
150 void AddNoise2Real(double snoise,int type_evol=0);
151
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);
157 void WriteSlicePPF(string cfname);
158 void NTupleCheck(POutPersist &pos,string ntname,unsigned long nent);
159
160 void SetPrtLevel(int lp=0) {lp_ = lp;}
161 void Print(void);
162
163//-------------------------------------------------------------------
164
165protected:
166 void init_default(void);
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);
170 void init_fftw(void);
171 void delete_fftw(void);
172 long manage_coefficients(void);
173 double compute_power_carte(void);
174 void check_array_alloc(void);
175 inline double pixelfilter(double x) // ATTENTION: seulement pour x>0
176 {return (x<0.025) ? 1.-x*x/6.*(1.-x*x/20.): sin(x)/x;}
177
178 // valeurs dans l'espace reel
179 long Nx_,Ny_,Nz_; vector<long> N_;
180 long NCz_,NTz_;
181 int_8 NRtot_,NFcoef_;
182
183 double Dx_,Dy_,Dz_; vector<double> D_;
184
185 // valeurs dans l'espace des K
186 double Dkx_,Dky_,Dkz_; vector<double> Dk_;
187 double Knyqx_,Knyqy_,Knyqz_; vector<double> Knyq_;
188 double Dk3_;
189 double dVol_, Vol_;
190
191 // la gestion de la FFT
192 bool is_set_fft_plan;
193 GEN3D_FFTW_PLAN pf_,pb_;
194 unsigned short nthread_;
195 int lp_;
196
197 // le stockage du Cube de donnees et les pointeurs
198 bool array_allocated_; // true if array has been allocated
199 unsigned short array_type; // 0=empty, 1=real, 2=complex
200 TArray< complex<GEN3D_TYPE> > T_;
201 GEN3D_FFTW_COMPLEX *fdata_;
202 TArray<GEN3D_TYPE> R_;
203 GEN3D_TYPE *data_;
204
205 // l'observateur
206 CosmoCalc *cosmo_;
207 GrowthFactor *growth_;
208 double good_dzinc_;
209 double redsh_ref_,kredsh_ref_,dred_ref_;
210 double compute_pk_redsh_ref_;
211 double h_ref_, loscom_ref_,dtrc_ref_, dlum_ref_, dang_ref_;
212 double growth_ref_, dsdz_growth_ref_;
213 double nu_ref_, dnu_ref_ ;
214 double xobs_[3];
215
216 double loscom_min_, loscom_max_;
217 vector<double> zred_, loscom_;
218 double loscom2zred_min_, loscom2zred_max_;
219 vector<double> loscom2zred_;
220
221 double zredmin_dpsd_, zredmax_dpsd_;
222 vector<double> dpsdfrzred_;
223
224};
225
226} // Fin du namespace SOPHYA
227
228#endif
Note: See TracBrowser for help on using the repository browser.