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

Last change on this file since 3621 was 3619, checked in by cmv, 16 years ago

add various #include<> for g++ 4.3 (jaunty 9.04), cmv 05/05/2009

File size: 6.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"
20
[3518]21#define WITH_FFTW_THREAD
[3115]22
[3521]23//NE PAS DECOMMENTER, UTILISEZ LA MAKEFILE #define GEN3D_FLOAT
24
[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_;}
[3330]79 long GetNx(void) {return Nx_;}
80 long GetNy(void) {return Ny_;}
81 long GetNz(void) {return Nz_;}
[3141]82
[3290]83 // Return |K_i| module relative to pixel indices
84 inline r_8 Kx(long i) {long ii=(i>Nx_/2)? Nx_-i :i; return ii*Dkx_;}
85 inline r_8 Ky(long j) {long jj=(j>Ny_/2)? Ny_-j :j; return jj*Dky_;}
86 inline r_8 Kz(long l) {return l*Dkz_;}
87
[3141]88 vector<r_8> GetDinc(void) {return D_;}
89 double GetDVol(void) {return dVol_;}
[3115]90 double GetVol(void) {return Vol_;}
[3141]91
92 vector<r_8> GetKinc(void) {return Dk_;}
93 vector<r_8> GetKnyq(void) {return Knyq_;}
[3115]94 double GetKmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_+Knyqz_*Knyqz_);}
[3141]95 double GetKTmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_);}
96 double GetKincMin(void)
97 {vector<r_8>::const_iterator it = min_element(Dk_.begin(), Dk_.end()); return *it;}
[3289]98 double GetKincMax(void)
99 {vector<r_8>::const_iterator it = max_element(Dk_.begin(), Dk_.end()); return *it;}
[3141]100 double GetKTincMin(void) {return min(Dk_[0],Dk_[1]);}
[3289]101 double GetKTincMax(void) {return max(Dk_[0],Dk_[1]);}
[3115]102
[3141]103 void ComputeFourier0(GenericFunc& pk_at_z);
104 void ComputeFourier(GenericFunc& pk_at_z);
[3115]105 void FilterByPixel(void);
[3141]106
[3115]107 void ComputeReal(void);
[3331]108 void ApplyGrowthFactor(int type_evol=1);
[3141]109
[3115]110 void ReComputeFourier(void);
111
[3141]112 int ComputeSpectrum(HistoErr& herr);
113 int ComputeSpectrum2D(Histo2DErr& herr);
[3330]114 int ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor);
115 int ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor);
[3141]116
[3134]117 int_8 VarianceFrReal(double R,double& var);
[3261]118 int_8 MeanSigma2(double& rm,double& rs2,double vmin=1.,double vmax=-1.
119 ,bool useout=false,double vout=0.);
[3320]120 int_8 MinMax(double& xmin,double& xmax,double vmin=1.,double vmax=-1.);
[3134]121 int_8 NumberOfBad(double vmin=-1.e+150,double vmax=1.e+150);
122 int_8 SetToVal(double vmin, double vmax,double val0=0.);
[3283]123 void ScaleOffset(double scalecube=1.,double offsetcube=0.);
[3115]124
125 void TurnFluct2Mass(void);
[3358]126 double TurnFluct2MeanNumber(double val_by_mpc3);
[3115]127 double ApplyPoisson(void);
128 double TurnNGal2Mass(FunRan& massdist,bool axeslog=false);
[3320]129 double TurnNGal2MassQuick(SchechterMassDist& schmdist);
[3115]130 double TurnMass2Flux(void);
[3349]131 //void AddAGN(double lfjy,double lsigma,double powlaw=0.);
[3331]132 void AddNoise2Real(double snoise,int type_evol=0);
[3115]133
[3141]134 void WriteFits(string cfname,int bitpix=FLOAT_IMG);
135 void ReadFits(string cfname);
136
137 void WritePPF(string cfname,bool write_real=true);
138 void ReadPPF(string cfname);
[3281]139 void WriteSlicePPF(string cfname);
[3524]140 void NTupleCheck(POutPersist &pos,string ntname,unsigned long nent);
[3141]141
[3150]142 void SetPrtLevel(int lp=0) {lp_ = lp;}
[3115]143 void Print(void);
144
[3199]145//-------------------------------------------------------------------
146
[3115]147protected:
[3349]148 void init_default(void);
[3141]149 void setsize(long nx,long ny,long nz,double dx,double dy,double dz);
150 void setalloc(void);
151 void setpointers(bool from_real);
[3154]152 void init_fftw(void);
[3349]153 void delete_fftw(void);
[3129]154 long manage_coefficients(void);
[3115]155 double compute_power_carte(void);
[3141]156 void check_array_alloc(void);
[3120]157 inline double pixelfilter(double x)
158 {return (x<0.025) ? 1.-x*x/6.*(1.-x*x/20.): sin(x)/x;}
[3115]159
[3154]160 // valeurs dans l'espace reel
[3141]161 long Nx_,Ny_,Nz_; vector<long> N_;
[3129]162 long NCz_,NTz_;
[3134]163 int_8 NRtot_;
[3141]164
165 double Dx_,Dy_,Dz_; vector<double> D_;
166
[3154]167 // valeurs dans l'espace des K
[3141]168 double Dkx_,Dky_,Dkz_; vector<double> Dk_;
169 double Knyqx_,Knyqy_,Knyqz_; vector<double> Knyq_;
170 double Dk3_;
[3115]171 double dVol_, Vol_;
172
[3154]173 // la gestion de la FFT
[3518]174 bool is_set_fft_plan;
175 GEN3D_FFTW_PLAN pf_,pb_;
[3115]176 unsigned short nthread_;
[3150]177 int lp_;
[3115]178
[3154]179 // le stockage du Cube de donnees et les pointeurs
[3141]180 bool array_allocated_; // true if array has been allocated
[3524]181 unsigned short array_type; // 0=empty, 1=real, 2=complex
[3518]182 TArray< complex<GEN3D_TYPE> > T_;
183 GEN3D_FFTW_COMPLEX *fdata_;
184 TArray<GEN3D_TYPE> R_;
185 GEN3D_TYPE *data_;
[3154]186
187 // l'observateur
[3157]188 CosmoCalc *cosmo_;
189 GrowthFactor *growth_;
[3271]190 double redsh_ref_,kredsh_ref_,dred_ref_;
191 double loscom_ref_,dtrc_ref_, dlum_ref_, dang_ref_;
192 double nu_ref_, dnu_ref_ ;
[3157]193 double xobs_[3];
[3271]194 double loscom_min_, loscom_max_;
[3157]195 vector<double> zred_, loscom_;
[3199]196 double loscom2zred_min_, loscom2zred_max_;
197 vector<double> loscom2zred_;
198
[3115]199};
200
[3325]201} // Fin du namespace SOPHYA
[3115]202
203#endif
Note: See TracBrowser for help on using the repository browser.