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

Last change on this file since 3154 was 3154, checked in by cmv, 19 years ago

modifs d'options de cmvobserv3d.cc
creation de init_fftw() et deplacement a l'initialisation
creation de SetObservator pour positionner l'observateur + ecriture FITS et PPF

(mais le codage des distance/redshift reste a coder)

cmv 19/01/2007

File size: 3.9 KB
Line 
1#ifndef GENEFLUCT3D_SEEN
2#define GENEFLUCT3D_SEEN
3
4#include "machdefs.h"
5#include "genericfunc.h"
6#include "tarray.h"
7#include "histerr.h"
8#include "hist2err.h"
9#include "perandom.h"
10
11#include "FFTW/fftw3.h"
12#include "FitsIO/fitsio.h"
13
14#include <vector>
15
16#include "pkspectrum.h"
17
18
19namespace SOPHYA {
20
21//-----------------------------------------------------------------------------------
22class GeneFluct3D {
23public:
24 GeneFluct3D(TArray< complex<r_8 > >& T);
25 virtual ~GeneFluct3D(void);
26
27 void SetNThread(unsigned short nthread=0) {nthread_ = nthread;}
28 void SetSize(long nx,long ny,long nz,double dx,double dy,double dz); // Mpc
29 void SetObservator(double redshref=0.,double kredshref=-1.);
30
31 TArray< complex<r_8> >& GetComplexArray(void) {return T_;}
32 fftw_complex* GetComplexPointer(void) {return fdata_;}
33 TArray<r_8>& GetRealArray(void) {return R_;}
34 r_8* GetRealPointer(void) {return data_;}
35
36 // Pour adressage data_[ip]
37 inline int_8 IndexR(long i,long j,long k) {return (int_8)(k+NTz_*(j+Ny_*i));}
38 // Pour adressage fdata_[ip][0-1]
39 inline int_8 IndexC(long i,long j,long k) {return (int_8)(k+NCz_*(j+Ny_*i));}
40 // - On peut aussi adresser:
41 // TArray< complex<r_8> >& pk = gf3d.GetComplexArray(); pk(i,j,k) = ...;
42 // TArray<r_8>& rgen = gf3d.GetRealArray(); rgen(i,j,k) = ...;
43
44 vector<long> GetNpix(void) {return N_;}
45 int_8 NPix(void) {return NRtot_;}
46
47 vector<r_8> GetDinc(void) {return D_;}
48 double GetDVol(void) {return dVol_;}
49 double GetVol(void) {return Vol_;}
50
51 vector<r_8> GetKinc(void) {return Dk_;}
52 vector<r_8> GetKnyq(void) {return Knyq_;}
53 double GetKmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_+Knyqz_*Knyqz_);}
54 double GetKTmax(void) {return sqrt(Knyqx_*Knyqx_+Knyqy_*Knyqy_);}
55 double GetKincMin(void)
56 {vector<r_8>::const_iterator it = min_element(Dk_.begin(), Dk_.end()); return *it;}
57 double GetKTincMin(void) {return min(Dk_[0],Dk_[1]);}
58
59 void ComputeFourier0(GenericFunc& pk_at_z);
60 void ComputeFourier(GenericFunc& pk_at_z);
61 void FilterByPixel(void);
62
63 void ComputeReal(void);
64
65 void ReComputeFourier(void);
66
67 int ComputeSpectrum(HistoErr& herr);
68 int ComputeSpectrum2D(Histo2DErr& herr);
69
70 int_8 VarianceFrReal(double R,double& var);
71 int_8 MeanSigma2(double& rm,double& rs2,double vmin=-1.e+150,double vmax=1.e+150);
72 int_8 NumberOfBad(double vmin=-1.e+150,double vmax=1.e+150);
73 int_8 SetToVal(double vmin, double vmax,double val0=0.);
74
75 void TurnFluct2Mass(void);
76 double TurnMass2MeanNumber(double n_by_mpc3);
77 double ApplyPoisson(void);
78 double TurnNGal2Mass(FunRan& massdist,bool axeslog=false);
79 double TurnMass2Flux(void);
80 void AddNoise2Real(double snoise);
81
82 void WriteFits(string cfname,int bitpix=FLOAT_IMG);
83 void ReadFits(string cfname);
84
85 void WritePPF(string cfname,bool write_real=true);
86 void ReadPPF(string cfname);
87
88 void SetPrtLevel(int lp=0) {lp_ = lp;}
89 void Print(void);
90
91protected:
92 void setsize(long nx,long ny,long nz,double dx,double dy,double dz);
93 void setalloc(void);
94 void setpointers(bool from_real);
95 void init_fftw(void);
96 long manage_coefficients(void);
97 double compute_power_carte(void);
98 void check_array_alloc(void);
99 inline double pixelfilter(double x)
100 {return (x<0.025) ? 1.-x*x/6.*(1.-x*x/20.): sin(x)/x;}
101
102 // valeurs dans l'espace reel
103 long Nx_,Ny_,Nz_; vector<long> N_;
104 long NCz_,NTz_;
105 int_8 NRtot_;
106
107 double Dx_,Dy_,Dz_; vector<double> D_;
108
109 // valeurs dans l'espace des K
110 double Dkx_,Dky_,Dkz_; vector<double> Dk_;
111 double Knyqx_,Knyqy_,Knyqz_; vector<double> Knyq_;
112 double Dk3_;
113 double dVol_, Vol_;
114
115 // la gestion de la FFT
116 fftw_plan pf_,pb_;
117 unsigned short nthread_;
118 int lp_;
119
120 // le stockage du Cube de donnees et les pointeurs
121 bool array_allocated_; // true if array has been allocated
122 TArray< complex<r_8> >& T_;
123 fftw_complex *fdata_;
124 TArray<r_8> R_;
125 double *data_;
126
127 // l'observateur
128 double redshref_,kredshref_;
129};
130
131} // Fin du namespace
132
133#endif
Note: See TracBrowser for help on using the repository browser.