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

Last change on this file since 3330 was 3330, checked in by cmv, 18 years ago

intro ComputeSpectrum avec soustraction de bruit et deconvolution par la fct de transfert du pixel cmv 01/10/2007

File size: 6.1 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
17#include "cosmocalc.h"
18#include "pkspectrum.h"
19
20
21namespace SOPHYA {
22
23//-----------------------------------------------------------------------------------
24class GeneFluct3D {
25public:
26 GeneFluct3D(TArray< complex<r_8 > >& T);
27 virtual ~GeneFluct3D(void);
28
29 void SetNThread(unsigned short nthread=0) {nthread_ = nthread;}
30 void SetSize(long nx,long ny,long nz,double dx,double dy,double dz); // Mpc
31 // Distance los comobile a l'observateur
32 void SetObservator(double redshref=0.,double kredshref=0.);
33 inline double DXcom(long i) {return i*Dx_ - xobs_[0];}
34 inline double DYcom(long j) {return j*Dy_ - xobs_[1];}
35 inline double DZcom(long k) {return k*Dz_ - xobs_[2];}
36 inline double Dcom(long i,long j,long k) {
37 double dx=DXcom(i), dy=DYcom(j), dz=DZcom(k);
38 return sqrt(dx*dx+dy*dy+dz*dz);
39 }
40 void SetCosmology(CosmoCalc& cosmo);
41 void SetGrowthFactor(GrowthFactor& growth);
42 long LosComRedshift(double zinc=0.001,long npoints=-1);
43
44 TArray< complex<r_8> >& GetComplexArray(void) {return T_;}
45 fftw_complex* GetComplexPointer(void) {return fdata_;}
46 TArray<r_8>& GetRealArray(void) {return R_;}
47 r_8* GetRealPointer(void) {return data_;}
48
49 // Pour adressage data_[ip]
50 inline int_8 IndexR(long i,long j,long k) {return (int_8)(k+NTz_*(j+Ny_*i));}
51 // Pour adressage fdata_[ip][0-1]
52 inline int_8 IndexC(long i,long j,long k) {return (int_8)(k+NCz_*(j+Ny_*i));}
53 // On peut aussi adresser:
54 // TArray< complex<r_8> >& pk = gf3d.GetComplexArray();
55 // pk(k,j,i) avec k=[0,NCz_[ j=[0,Ny_[ i=[0,Nx_[
56 // pk[IndexC(i,j,k)]
57 // TArray<r_8>& rgen = gf3d.GetRealArray();
58 // rgen(k,j,i) avec k=[0,NTz_[ j=[0,Ny_[ i=[0,Nx_[
59 // mais seul k=[0,Nz_[ est utile
60 // rgen[IndexR(i,j,k)]
61 // ATTENTION: TArray adresse en memoire a l'envers du C !
62 // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
63
64 vector<long> GetNpix(void) {return N_;}
65 int_8 NPix(void) {return NRtot_;}
66 long GetNx(void) {return Nx_;}
67 long GetNy(void) {return Ny_;}
68 long GetNz(void) {return Nz_;}
69
70 // Return |K_i| module relative to pixel indices
71 inline r_8 Kx(long i) {long ii=(i>Nx_/2)? Nx_-i :i; return ii*Dkx_;}
72 inline r_8 Ky(long j) {long jj=(j>Ny_/2)? Ny_-j :j; return jj*Dky_;}
73 inline r_8 Kz(long l) {return l*Dkz_;}
74
75 vector<r_8> GetDinc(void) {return D_;}
76 double GetDVol(void) {return dVol_;}
77 double GetVol(void) {return Vol_;}
78
79 vector<r_8> GetKinc(void) {return Dk_;}
80 vector<r_8> GetKnyq(void) {return Knyq_;}
81 double GetKmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_+Knyqz_*Knyqz_);}
82 double GetKTmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_);}
83 double GetKincMin(void)
84 {vector<r_8>::const_iterator it = min_element(Dk_.begin(), Dk_.end()); return *it;}
85 double GetKincMax(void)
86 {vector<r_8>::const_iterator it = max_element(Dk_.begin(), Dk_.end()); return *it;}
87 double GetKTincMin(void) {return min(Dk_[0],Dk_[1]);}
88 double GetKTincMax(void) {return max(Dk_[0],Dk_[1]);}
89
90 void ComputeFourier0(GenericFunc& pk_at_z);
91 void ComputeFourier(GenericFunc& pk_at_z);
92 void FilterByPixel(void);
93
94 void ComputeReal(void);
95 void ApplyGrowthFactor(void);
96
97 void ReComputeFourier(void);
98
99 int ComputeSpectrum(HistoErr& herr);
100 int ComputeSpectrum2D(Histo2DErr& herr);
101 int ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor);
102 int ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor);
103
104 int_8 VarianceFrReal(double R,double& var);
105 int_8 MeanSigma2(double& rm,double& rs2,double vmin=1.,double vmax=-1.
106 ,bool useout=false,double vout=0.);
107 int_8 MinMax(double& xmin,double& xmax,double vmin=1.,double vmax=-1.);
108 int_8 NumberOfBad(double vmin=-1.e+150,double vmax=1.e+150);
109 int_8 SetToVal(double vmin, double vmax,double val0=0.);
110 void ScaleOffset(double scalecube=1.,double offsetcube=0.);
111
112 void TurnFluct2Mass(void);
113 double TurnMass2HIMass(double m_by_mpc3);
114 double TurnMass2MeanNumber(double n_by_mpc3);
115 double ApplyPoisson(void);
116 double TurnNGal2Mass(FunRan& massdist,bool axeslog=false);
117 double TurnNGal2MassQuick(SchechterMassDist& schmdist);
118 double TurnMass2Flux(void);
119 void AddAGN(double lfjy,double lsigma,double powlaw=0.);
120 void AddNoise2Real(double snoise,bool with_evol=false);
121
122 void WriteFits(string cfname,int bitpix=FLOAT_IMG);
123 void ReadFits(string cfname);
124
125 void WritePPF(string cfname,bool write_real=true);
126 void ReadPPF(string cfname);
127 void WriteSlicePPF(string cfname);
128
129 void SetPrtLevel(int lp=0) {lp_ = lp;}
130 void Print(void);
131
132//-------------------------------------------------------------------
133
134protected:
135 void setsize(long nx,long ny,long nz,double dx,double dy,double dz);
136 void setalloc(void);
137 void setpointers(bool from_real);
138 void init_fftw(void);
139 long manage_coefficients(void);
140 double compute_power_carte(void);
141 void check_array_alloc(void);
142 inline double pixelfilter(double x)
143 {return (x<0.025) ? 1.-x*x/6.*(1.-x*x/20.): sin(x)/x;}
144
145 // valeurs dans l'espace reel
146 long Nx_,Ny_,Nz_; vector<long> N_;
147 long NCz_,NTz_;
148 int_8 NRtot_;
149
150 double Dx_,Dy_,Dz_; vector<double> D_;
151
152 // valeurs dans l'espace des K
153 double Dkx_,Dky_,Dkz_; vector<double> Dk_;
154 double Knyqx_,Knyqy_,Knyqz_; vector<double> Knyq_;
155 double Dk3_;
156 double dVol_, Vol_;
157
158 // la gestion de la FFT
159 fftw_plan pf_,pb_;
160 unsigned short nthread_;
161 int lp_;
162
163 // le stockage du Cube de donnees et les pointeurs
164 bool array_allocated_; // true if array has been allocated
165 TArray< complex<r_8> >& T_;
166 fftw_complex *fdata_;
167 TArray<r_8> R_;
168 double *data_;
169
170 // l'observateur
171 CosmoCalc *cosmo_;
172 GrowthFactor *growth_;
173 double redsh_ref_,kredsh_ref_,dred_ref_;
174 double loscom_ref_,dtrc_ref_, dlum_ref_, dang_ref_;
175 double nu_ref_, dnu_ref_ ;
176 double xobs_[3];
177 double loscom_min_, loscom_max_;
178 vector<double> zred_, loscom_;
179 double loscom2zred_min_, loscom2zred_max_;
180 vector<double> loscom2zred_;
181
182};
183
184} // Fin du namespace SOPHYA
185
186#endif
Note: See TracBrowser for help on using the repository browser.