1 | #include "machdefs.h"
|
---|
2 | #include <iostream>
|
---|
3 | #include <stdlib.h>
|
---|
4 | #include <stdio.h>
|
---|
5 | #include <string.h>
|
---|
6 | #include <math.h>
|
---|
7 | #include <unistd.h>
|
---|
8 |
|
---|
9 | #include "tarray.h"
|
---|
10 | #include "pexceptions.h"
|
---|
11 | #include "perandom.h"
|
---|
12 | #include "srandgen.h"
|
---|
13 |
|
---|
14 | #include "fabtcolread.h"
|
---|
15 | #include "fabtwriter.h"
|
---|
16 | #include "fioarr.h"
|
---|
17 | #include "ntuple.h"
|
---|
18 |
|
---|
19 | #include "arrctcast.h"
|
---|
20 |
|
---|
21 | #include "constcosmo.h"
|
---|
22 | #include "geneutils.h"
|
---|
23 | #include "schechter.h"
|
---|
24 |
|
---|
25 | #include "genefluct3d.h"
|
---|
26 |
|
---|
27 | #if defined(GEN3D_FLOAT)
|
---|
28 | #define GEN3D_FFTW_INIT_THREADS fftwf_init_threads
|
---|
29 | #define GEN3D_FFTW_CLEANUP_THREADS fftwf_cleanup_threads
|
---|
30 | #define GEN3D_FFTW_PLAN_WITH_NTHREADS fftwf_plan_with_nthreads
|
---|
31 | #define GEN3D_FFTW_PLAN_DFT_R2C_3D fftwf_plan_dft_r2c_3d
|
---|
32 | #define GEN3D_FFTW_PLAN_DFT_C2R_3D fftwf_plan_dft_c2r_3d
|
---|
33 | #define GEN3D_FFTW_DESTROY_PLAN fftwf_destroy_plan
|
---|
34 | #define GEN3D_FFTW_EXECUTE fftwf_execute
|
---|
35 | #else
|
---|
36 | #define GEN3D_FFTW_INIT_THREADS fftw_init_threads
|
---|
37 | #define GEN3D_FFTW_CLEANUP_THREADS fftw_cleanup_threads
|
---|
38 | #define GEN3D_FFTW_PLAN_WITH_NTHREADS fftw_plan_with_nthreads
|
---|
39 | #define GEN3D_FFTW_PLAN_DFT_R2C_3D fftw_plan_dft_r2c_3d
|
---|
40 | #define GEN3D_FFTW_PLAN_DFT_C2R_3D fftw_plan_dft_c2r_3d
|
---|
41 | #define GEN3D_FFTW_DESTROY_PLAN fftw_destroy_plan
|
---|
42 | #define GEN3D_FFTW_EXECUTE fftw_execute
|
---|
43 | #endif
|
---|
44 |
|
---|
45 | #define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag()))
|
---|
46 |
|
---|
47 | namespace SOPHYA {
|
---|
48 |
|
---|
49 | //-------------------------------------------------------
|
---|
50 | GeneFluct3D::GeneFluct3D(long nx,long ny,long nz,double dx,double dy,double dz
|
---|
51 | ,unsigned short nthread,int lp)
|
---|
52 | {
|
---|
53 | init_default();
|
---|
54 |
|
---|
55 | lp_ = lp;
|
---|
56 | nthread_ = nthread;
|
---|
57 |
|
---|
58 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
59 | setalloc();
|
---|
60 | setpointers(false);
|
---|
61 | init_fftw();
|
---|
62 | }
|
---|
63 |
|
---|
64 | GeneFluct3D::GeneFluct3D(unsigned short nthread)
|
---|
65 | {
|
---|
66 | init_default();
|
---|
67 | setsize(2,2,2,1.,1.,1.);
|
---|
68 | nthread_ = nthread;
|
---|
69 | setalloc();
|
---|
70 | setpointers(false);
|
---|
71 | init_fftw();
|
---|
72 | }
|
---|
73 |
|
---|
74 | GeneFluct3D::~GeneFluct3D(void)
|
---|
75 | {
|
---|
76 | delete_fftw();
|
---|
77 | }
|
---|
78 |
|
---|
79 | //-------------------------------------------------------
|
---|
80 | void GeneFluct3D::init_default(void)
|
---|
81 | {
|
---|
82 | Nx_ = Ny_ = Nz_ = 0;
|
---|
83 | is_set_fft_plan = false;
|
---|
84 | nthread_ = 0;
|
---|
85 | lp_ = 0;
|
---|
86 | array_allocated_ = false; array_type = 0;
|
---|
87 | cosmo_ = NULL;
|
---|
88 | growth_ = NULL;
|
---|
89 | good_dzinc_ = 0.01;
|
---|
90 | compute_pk_redsh_ref_ = -999.;
|
---|
91 | redsh_ref_ = -999.;
|
---|
92 | kredsh_ref_ = 0.;
|
---|
93 | dred_ref_ = -999.;
|
---|
94 | h_ref_ = -999.;
|
---|
95 | loscom_ref_ = -999.;
|
---|
96 | dtrc_ref_ = dlum_ref_ = dang_ref_ = -999.;
|
---|
97 | nu_ref_ = dnu_ref_ = -999.;
|
---|
98 | loscom_min_ = loscom_max_ = -999.;
|
---|
99 | loscom2zred_min_ = loscom2zred_max_ = 0.;
|
---|
100 | xobs_[0] = xobs_[1] = xobs_[2] = 0.;
|
---|
101 | zred_.resize(0);
|
---|
102 | loscom_.resize(0);
|
---|
103 | loscom2zred_.resize(0);
|
---|
104 | }
|
---|
105 |
|
---|
106 | void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz)
|
---|
107 | {
|
---|
108 | if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<<nx<<","<<ny<<","<<nz
|
---|
109 | <<" D="<<dx<<","<<dy<<","<<dz<<endl;
|
---|
110 | if(nx<=0 || dx<=0.) {
|
---|
111 | const char *bla = "GeneFluct3D::setsize_Error: bad value(s) for nn/dx";
|
---|
112 | cout<<bla<<endl; throw ParmError(bla);
|
---|
113 | }
|
---|
114 |
|
---|
115 | // Les tailles des tableaux
|
---|
116 | Nx_ = nx;
|
---|
117 | Ny_ = ny; if(Ny_ <= 0) Ny_ = Nx_;
|
---|
118 | Nz_ = nz; if(Nz_ <= 0) Nz_ = Nx_;
|
---|
119 | N_.resize(0); N_.push_back(Nx_); N_.push_back(Ny_); N_.push_back(Nz_);
|
---|
120 | NRtot_ = Nx_*Ny_*Nz_; // nombre de pixels dans le survey
|
---|
121 | NCz_ = Nz_/2 +1;
|
---|
122 | NTz_ = 2*NCz_;
|
---|
123 |
|
---|
124 | // le pas dans l'espace (Mpc)
|
---|
125 | Dx_ = dx;
|
---|
126 | Dy_ = dy; if(Dy_ <= 0.) Dy_ = Dx_;
|
---|
127 | Dz_ = dz; if(Dz_ <= 0.) Dz_ = Dx_;
|
---|
128 | D_.resize(0); D_.push_back(Dx_); D_.push_back(Dy_); D_.push_back(Dz_);
|
---|
129 | dVol_ = Dx_*Dy_*Dz_;
|
---|
130 | Vol_ = (Nx_*Dx_)*(Ny_*Dy_)*(Nz_*Dz_);
|
---|
131 |
|
---|
132 | // Le pas dans l'espace de Fourier (Mpc^-1)
|
---|
133 | Dkx_ = 2.*M_PI/(Nx_*Dx_);
|
---|
134 | Dky_ = 2.*M_PI/(Ny_*Dy_);
|
---|
135 | Dkz_ = 2.*M_PI/(Nz_*Dz_);
|
---|
136 | Dk_.resize(0); Dk_.push_back(Dkx_); Dk_.push_back(Dky_); Dk_.push_back(Dkz_);
|
---|
137 | Dk3_ = Dkx_*Dky_*Dkz_;
|
---|
138 |
|
---|
139 | // La frequence de Nyquist en k (Mpc^-1)
|
---|
140 | Knyqx_ = M_PI/Dx_;
|
---|
141 | Knyqy_ = M_PI/Dy_;
|
---|
142 | Knyqz_ = M_PI/Dz_;
|
---|
143 | Knyq_.resize(0); Knyq_.push_back(Knyqx_); Knyq_.push_back(Knyqy_); Knyq_.push_back(Knyqz_);
|
---|
144 | }
|
---|
145 |
|
---|
146 | void GeneFluct3D::setalloc(void)
|
---|
147 | {
|
---|
148 | #if defined(GEN3D_FLOAT)
|
---|
149 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc FLOAT ---"<<endl;
|
---|
150 | #else
|
---|
151 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc DOUBLE ---"<<endl;
|
---|
152 | #endif
|
---|
153 | // Dimensionnement du tableau complex<r_8>
|
---|
154 | // ATTENTION: TArray adresse en memoire a l'envers du C
|
---|
155 | // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
|
---|
156 | sa_size_t SzK_[3] = {NCz_,Ny_,Nx_}; // a l'envers
|
---|
157 | try {
|
---|
158 | T_.ReSize(3,SzK_);
|
---|
159 | array_allocated_ = true; array_type=0;
|
---|
160 | if(lp_>1) cout<<" allocating: "<<T_.Size()*sizeof(complex<GEN3D_TYPE>)/1.e6<<" Mo"<<endl;
|
---|
161 | } catch (...) {
|
---|
162 | cout<<"GeneFluct3D::setalloc_Error: Problem allocating T_"<<endl;
|
---|
163 | }
|
---|
164 | }
|
---|
165 |
|
---|
166 | void GeneFluct3D::setpointers(bool from_real)
|
---|
167 | {
|
---|
168 | if(lp_>1) cout<<"--- GeneFluct3D::setpointers ---"<<endl;
|
---|
169 | if(from_real) T_ = ArrCastR2C(R_);
|
---|
170 | else R_ = ArrCastC2R(T_);
|
---|
171 | // On remplit les pointeurs
|
---|
172 | fdata_ = (GEN3D_FFTW_COMPLEX *) (&T_(0,0,0));
|
---|
173 | data_ = (GEN3D_TYPE *) (&R_(0,0,0));
|
---|
174 | }
|
---|
175 |
|
---|
176 | void GeneFluct3D::init_fftw(void)
|
---|
177 | {
|
---|
178 | if( is_set_fft_plan ) delete_fftw();
|
---|
179 |
|
---|
180 | // --- Initialisation de fftw3 (attention data est sur-ecrit a l'init)
|
---|
181 | if(lp_>1) cout<<"--- GeneFluct3D::init_fftw ---"<<endl;
|
---|
182 | #ifdef WITH_FFTW_THREAD
|
---|
183 | if(nthread_>0) {
|
---|
184 | cout<<"...Computing with "<<nthread_<<" threads"<<endl;
|
---|
185 | GEN3D_FFTW_INIT_THREADS();
|
---|
186 | GEN3D_FFTW_PLAN_WITH_NTHREADS(nthread_);
|
---|
187 | }
|
---|
188 | #endif
|
---|
189 | if(lp_>1) cout<<"...forward plan"<<endl;
|
---|
190 | pf_ = GEN3D_FFTW_PLAN_DFT_R2C_3D(Nx_,Ny_,Nz_,data_,fdata_,FFTW_ESTIMATE);
|
---|
191 | if(lp_>1) cout<<"...backward plan"<<endl;
|
---|
192 | pb_ = GEN3D_FFTW_PLAN_DFT_C2R_3D(Nx_,Ny_,Nz_,fdata_,data_,FFTW_ESTIMATE);
|
---|
193 | is_set_fft_plan = true;
|
---|
194 | }
|
---|
195 |
|
---|
196 | void GeneFluct3D::delete_fftw(void)
|
---|
197 | {
|
---|
198 | if( !is_set_fft_plan ) return;
|
---|
199 | GEN3D_FFTW_DESTROY_PLAN(pf_);
|
---|
200 | GEN3D_FFTW_DESTROY_PLAN(pb_);
|
---|
201 | #ifdef WITH_FFTW_THREAD
|
---|
202 | if(nthread_>0) GEN3D_FFTW_CLEANUP_THREADS();
|
---|
203 | #endif
|
---|
204 | is_set_fft_plan = false;
|
---|
205 | }
|
---|
206 |
|
---|
207 | void GeneFluct3D::check_array_alloc(void)
|
---|
208 | // Pour tester si le tableau T_ est alloue
|
---|
209 | {
|
---|
210 | if(array_allocated_) return;
|
---|
211 | char bla[90];
|
---|
212 | sprintf(bla,"GeneFluct3D::check_array_alloc_Error: array is not allocated");
|
---|
213 | cout<<bla<<endl; throw ParmError(bla);
|
---|
214 | }
|
---|
215 |
|
---|
216 | //-------------------------------------------------------
|
---|
217 | void GeneFluct3D::SetObservator(double redshref,double kredshref)
|
---|
218 | // L'observateur est au redshift z=0
|
---|
219 | // est situe sur la "perpendiculaire" a la face x,y
|
---|
220 | // issue du centre de cette face
|
---|
221 | // Il faut positionner le cube sur l'axe des z cad des redshifts:
|
---|
222 | // redshref = redshift de reference
|
---|
223 | // Si redshref<0 alors redshref=0
|
---|
224 | // kredshref = indice (en double) correspondant a ce redshift
|
---|
225 | // Si kredshref<0 alors kredshref=nz/2 (milieu du cube)
|
---|
226 | // Exemple: redshref=1.5 kredshref=250.75
|
---|
227 | // -> Le pixel i=nx/2 j=ny/2 k=250.75 est au redshift 1.5
|
---|
228 | {
|
---|
229 | if(redshref<0.) redshref = 0.;
|
---|
230 | if(kredshref<0.) {
|
---|
231 | if(Nz_<=0) {
|
---|
232 | const char *bla = "GeneFluct3D::SetObservator_Error: for kredsh_ref<0 define cube geometry first";
|
---|
233 | cout<<bla<<endl; throw ParmError(bla);
|
---|
234 | }
|
---|
235 | kredshref = Nz_/2.;
|
---|
236 | }
|
---|
237 | redsh_ref_ = redshref;
|
---|
238 | kredsh_ref_ = kredshref;
|
---|
239 | if(lp_>0)
|
---|
240 | cout<<"--- GeneFluct3D::SetObservator zref="<<redsh_ref_<<" kref="<<kredsh_ref_<<endl;
|
---|
241 | }
|
---|
242 |
|
---|
243 | void GeneFluct3D::SetCosmology(CosmoCalc& cosmo)
|
---|
244 | {
|
---|
245 | cosmo_ = &cosmo;
|
---|
246 | if(lp_>1) cosmo_->Print();
|
---|
247 | }
|
---|
248 |
|
---|
249 | void GeneFluct3D::SetGrowthFactor(GrowthFactor& growth)
|
---|
250 | {
|
---|
251 | growth_ = &growth;
|
---|
252 | }
|
---|
253 |
|
---|
254 | long GeneFluct3D::LosComRedshift(double zinc,long npoints)
|
---|
255 | // Given a position of the cube relative to the observer
|
---|
256 | // and a cosmology
|
---|
257 | // (SetObservator() and SetCosmology() should have been called !)
|
---|
258 | // This routine filled:
|
---|
259 | // the vector "zred_" of scanned redshift (by zinc increments)
|
---|
260 | // the vector "loscom_" of corresponding los comoving distance
|
---|
261 | // -- Input:
|
---|
262 | // zinc : redshift increment for computation
|
---|
263 | // npoints : number of points required for inverting loscom -> zred
|
---|
264 | //
|
---|
265 | {
|
---|
266 | if(zinc<=0.) zinc = good_dzinc_;
|
---|
267 | if(lp_>0) cout<<"--- LosComRedshift: zinc="<<zinc<<" , npoints="<<npoints<<endl;
|
---|
268 |
|
---|
269 | if(cosmo_ == NULL || growth_==NULL || redsh_ref_<0.) {
|
---|
270 | const char *bla = "GeneFluct3D::LosComRedshift_Error: set Observator, Cosmology and Growth first";
|
---|
271 | cout<<bla<<endl; throw ParmError(bla);
|
---|
272 | }
|
---|
273 |
|
---|
274 | // La distance angulaire/luminosite/Dnu au pixel de reference
|
---|
275 | dred_ref_ = Dz_/(cosmo_->Dhubble()/cosmo_->E(redsh_ref_));
|
---|
276 | loscom_ref_ = cosmo_->Dloscom(redsh_ref_);
|
---|
277 | dtrc_ref_ = cosmo_->Dtrcom(redsh_ref_);
|
---|
278 | dlum_ref_ = cosmo_->Dlum(redsh_ref_);
|
---|
279 | dang_ref_ = cosmo_->Dang(redsh_ref_);
|
---|
280 | h_ref_ = cosmo_->H(redsh_ref_);
|
---|
281 | growth_ref_ = (*growth_)(redsh_ref_);
|
---|
282 | dsdz_growth_ref_ = growth_->DsDz(redsh_ref_,zinc);
|
---|
283 | nu_ref_ = Fr_HyperFin_Par/(1.+redsh_ref_); // GHz
|
---|
284 | dnu_ref_ = Fr_HyperFin_Par *dred_ref_/pow(1.+redsh_ref_,2.); // GHz
|
---|
285 | if(lp_>0) {
|
---|
286 | cout<<"...reference pixel redshref="<<redsh_ref_
|
---|
287 | <<", dredref="<<dred_ref_
|
---|
288 | <<", nuref="<<nu_ref_ <<" GHz"
|
---|
289 | <<", dnuref="<<dnu_ref_ <<" GHz"<<endl
|
---|
290 | <<" H="<<h_ref_<<" km/s/Mpc"
|
---|
291 | <<", D="<<growth_ref_
|
---|
292 | <<", dD/dz="<<dsdz_growth_ref_<<endl
|
---|
293 | <<" dlosc="<<loscom_ref_<<" Mpc com"
|
---|
294 | <<", dtrc="<<dtrc_ref_<<" Mpc com"
|
---|
295 | <<", dlum="<<dlum_ref_<<" Mpc"
|
---|
296 | <<", dang="<<dang_ref_<<" Mpc"<<endl;
|
---|
297 | }
|
---|
298 |
|
---|
299 | // On calcule les coordonnees de l'observateur dans le repere du cube
|
---|
300 | // cad dans le repere ou l'origine est au centre du pixel i=j=l=0.
|
---|
301 | // L'observateur est sur un axe centre sur le milieu de la face Oxy
|
---|
302 | xobs_[0] = Nx_/2.*Dx_;
|
---|
303 | xobs_[1] = Ny_/2.*Dy_;
|
---|
304 | xobs_[2] = kredsh_ref_*Dz_ - loscom_ref_;
|
---|
305 |
|
---|
306 | // L'observateur est-il dans le cube?
|
---|
307 | bool obs_in_cube = false;
|
---|
308 | if(xobs_[2]>=0. && xobs_[2]<=Nz_*Dz_) obs_in_cube = true;
|
---|
309 |
|
---|
310 | // Find MINIMUM los com distance to the observer:
|
---|
311 | // c'est le centre de la face a k=0
|
---|
312 | // (ou zero si l'observateur est dans le cube)
|
---|
313 | loscom_min_ = 0.;
|
---|
314 | if(!obs_in_cube) loscom_min_ = -xobs_[2];
|
---|
315 |
|
---|
316 | // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED
|
---|
317 | if(loscom_min_<=1.e-50)
|
---|
318 | for(int i=0;i<10;i++)
|
---|
319 | cout<<"ATTENTION TOUTES LES PARTIES DU CODE NE MARCHENT PAS POUR UN OBSERVATEUR DANS LE CUBE"<<endl;
|
---|
320 | // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED
|
---|
321 |
|
---|
322 |
|
---|
323 | // Find MAXIMUM los com distance to the observer:
|
---|
324 | // ou que soit positionne l'observateur, la distance
|
---|
325 | // maximal est sur un des coins du cube
|
---|
326 | loscom_max_ = 0.;
|
---|
327 | for(long i=0;i<=1;i++) {
|
---|
328 | double dx2 = DXcom(i*(Nx_-1)); dx2 *= dx2;
|
---|
329 | for(long j=0;j<=1;j++) {
|
---|
330 | double dy2 = DYcom(j*(Ny_-1)); dy2 *= dy2;
|
---|
331 | for(long k=0;k<=1;k++) {
|
---|
332 | double dz2 = DZcom(k*(Nz_-1)); dz2 *= dz2;
|
---|
333 | dz2 = sqrt(dx2+dy2+dz2);
|
---|
334 | if(dz2>loscom_max_) loscom_max_ = dz2;
|
---|
335 | }
|
---|
336 | }
|
---|
337 | }
|
---|
338 | if(lp_>0) {
|
---|
339 | cout<<"...zref="<<redsh_ref_<<" kzref="<<kredsh_ref_<<" losref="<<loscom_ref_<<" Mpc\n"
|
---|
340 | <<" xobs="<<xobs_[0]<<" , "<<xobs_[1]<<" , "<<xobs_[2]<<" Mpc "
|
---|
341 | <<" in_cube="<<obs_in_cube
|
---|
342 | <<" loscom_min="<<loscom_min_<<" loscom_max="<<loscom_max_<<" Mpc (com)"<<endl;
|
---|
343 | }
|
---|
344 |
|
---|
345 | // Fill the corresponding vectors for loscom and zred
|
---|
346 | // Be shure to have one dlc<loscom_min and one dlc>loscom_max
|
---|
347 | double zmin = 0., dlcmin=0.;
|
---|
348 | while(1) {
|
---|
349 | if(lp_>0)
|
---|
350 | cout<<"...Filling zred starting at zmin="<<zmin<<" with zinc="<<zinc
|
---|
351 | <<", loscom_min-max=["<<loscom_min_<<","<<loscom_max_<<"]"<<endl;
|
---|
352 | zred_.resize(0); loscom_.resize(0);
|
---|
353 | for(double z=zmin; ; z+=zinc) {
|
---|
354 | double dlc = cosmo_->Dloscom(z);
|
---|
355 | if(dlc<loscom_min_) {
|
---|
356 | dlcmin = dlc;
|
---|
357 | zmin = z;
|
---|
358 | zred_.resize(0); loscom_.resize(0);
|
---|
359 | }
|
---|
360 | zred_.push_back(z);
|
---|
361 | loscom_.push_back(dlc);
|
---|
362 | z += zinc;
|
---|
363 | if(dlc>loscom_max_) {
|
---|
364 | if(lp_>0)
|
---|
365 | cout<<" Min: z="<<zmin<<" dlc="<<dlcmin<<", Max: z="<<z<<" dlc="<<dlc<<endl;
|
---|
366 | break; // on sort apres avoir stoque un dlc>dlcmax
|
---|
367 | }
|
---|
368 | }
|
---|
369 | if(zred_.size()>=10) break;
|
---|
370 | zinc /= 10.;
|
---|
371 | cout<<" not enough points ("<<zred_.size()<<") for zref, retry with zinc="<<zinc<<endl;
|
---|
372 | }
|
---|
373 |
|
---|
374 | if(lp_>0) {
|
---|
375 | long n = zred_.size();
|
---|
376 | cout<<"...zred/loscom tables[zinc="<<zinc<<"]: n="<<n;
|
---|
377 | if(n>0) cout<<" z="<<zred_[0]<<" -> d="<<loscom_[0];
|
---|
378 | if(n>1) cout<<" , z="<<zred_[n-1]<<" -> d="<<loscom_[n-1];
|
---|
379 | cout<<endl;
|
---|
380 | }
|
---|
381 |
|
---|
382 | // Compute the parameters and tables needed for inversion loscom->zred
|
---|
383 | if(npoints<3) npoints = zred_.size();
|
---|
384 | InverseFunc invfun(zred_,loscom_);
|
---|
385 | invfun.ComputeParab(npoints,loscom2zred_);
|
---|
386 | loscom2zred_min_ = invfun.YMin();
|
---|
387 | loscom2zred_max_ = invfun.YMax();
|
---|
388 |
|
---|
389 | if(lp_>0) {
|
---|
390 | long n = loscom2zred_.size();
|
---|
391 | cout<<"...loscom -> zred[npoints="<<npoints<<"]: n="<<n
|
---|
392 | <<" los_min="<<loscom2zred_min_
|
---|
393 | <<" los_max="<<loscom2zred_max_
|
---|
394 | <<" -> zred=[";
|
---|
395 | if(n>0) cout<<loscom2zred_[0];
|
---|
396 | cout<<",";
|
---|
397 | if(n>1) cout<<loscom2zred_[n-1];
|
---|
398 | cout<<"]"<<endl;
|
---|
399 | if(lp_>1 && n>0)
|
---|
400 | for(int i=0;i<n;i++)
|
---|
401 | if(i<2 || abs(i-n/2)<2 || i>=n-2)
|
---|
402 | cout<<" i="<<i
|
---|
403 | <<" d="<<loscom2zred_min_+i*(loscom2zred_max_-loscom2zred_min_)/(n-1.)
|
---|
404 | <<" Mpc z="<<loscom2zred_[i]<<endl;
|
---|
405 | }
|
---|
406 |
|
---|
407 |
|
---|
408 | // Compute the table for D'(z)/D(z) where D'(z)=dD/dEta (Eta is conformal time)
|
---|
409 | zredmin_dpsd_ = zred_[0];
|
---|
410 | zredmax_dpsd_ = zred_[zred_.size()-1];
|
---|
411 | long nptd = long(sqrt(Nx_*Nx_ + Ny_*Ny_ +Nz_*Nz_));
|
---|
412 | if(nptd<10) nptd = 10;
|
---|
413 | good_dzinc_ = (zredmax_dpsd_ - zredmin_dpsd_)/nptd;
|
---|
414 | if(lp_>0) cout<<"...good_dzinc changed to "<<good_dzinc_<<endl;
|
---|
415 | dpsdfrzred_.resize(0);
|
---|
416 | double dz = (zredmax_dpsd_ - zredmin_dpsd_)/(nptd-1);
|
---|
417 | int nmod = nptd/5; if(nmod==0) nmod = 1;
|
---|
418 | if(lp_>0) cout<<"...Compute table D'/D on "<<nptd<<" pts for z["
|
---|
419 | <<zredmin_dpsd_<<","<<zredmax_dpsd_<<"]"<<endl;
|
---|
420 | for(int i=0;i<nptd;i++) {
|
---|
421 | double z = zredmin_dpsd_ + i*dz;
|
---|
422 | double v = -cosmo_->H(z) * growth_->DsDz(z,good_dzinc_) / (*growth_)(z);
|
---|
423 | dpsdfrzred_.push_back(v);
|
---|
424 | if(lp_ && (i%nmod==0 || i==nptd-1)) cout<<" z="<<z<<" D'/D="<<v<<" km/s/Mpc"<<endl;
|
---|
425 | }
|
---|
426 |
|
---|
427 | return zred_.size();
|
---|
428 | }
|
---|
429 |
|
---|
430 | //-------------------------------------------------------
|
---|
431 | void GeneFluct3D::WriteFits(string cfname,int bitpix)
|
---|
432 | {
|
---|
433 | cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
|
---|
434 | try {
|
---|
435 | FitsImg3DWriter fwrt(cfname.c_str(),bitpix,5);
|
---|
436 | fwrt.WriteKey("NX",Nx_," axe transverse 1");
|
---|
437 | fwrt.WriteKey("NY",Ny_," axe transverse 2");
|
---|
438 | fwrt.WriteKey("NZ",Nz_," axe longitudinal (redshift)");
|
---|
439 | fwrt.WriteKey("DX",Dx_," Mpc");
|
---|
440 | fwrt.WriteKey("DY",Dy_," Mpc");
|
---|
441 | fwrt.WriteKey("DZ",Dz_," Mpc");
|
---|
442 | fwrt.WriteKey("DKX",Dkx_," Mpc^-1");
|
---|
443 | fwrt.WriteKey("DKY",Dky_," Mpc^-1");
|
---|
444 | fwrt.WriteKey("DKZ",Dkz_," Mpc^-1");
|
---|
445 | fwrt.WriteKey("ZREFPK",compute_pk_redsh_ref_," Pk computed redshift");
|
---|
446 | fwrt.WriteKey("ZREF",redsh_ref_," reference redshift");
|
---|
447 | fwrt.WriteKey("KZREF",kredsh_ref_," reference redshift on z axe");
|
---|
448 | fwrt.WriteKey("HREF",h_ref_," ref H(z)");
|
---|
449 | fwrt.WriteKey("Growth",Dref()," growth at reference redshift");
|
---|
450 | fwrt.WriteKey("GrowthPk",DrefPk()," growth at Pk computed redshift");
|
---|
451 | fwrt.Write(R_);
|
---|
452 | } catch (PThrowable & exc) {
|
---|
453 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
454 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
455 | return;
|
---|
456 | } catch (...) {
|
---|
457 | cout<<" some other exception was caught !"<<endl;
|
---|
458 | return;
|
---|
459 | }
|
---|
460 | }
|
---|
461 |
|
---|
462 | void GeneFluct3D::ReadFits(string cfname)
|
---|
463 | {
|
---|
464 | cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
|
---|
465 | try {
|
---|
466 | FitsImg3DRead fimg(cfname.c_str(),0,5);
|
---|
467 | fimg.Read(R_);
|
---|
468 | long nx = fimg.ReadKeyL("NX");
|
---|
469 | long ny = fimg.ReadKeyL("NY");
|
---|
470 | long nz = fimg.ReadKeyL("NZ");
|
---|
471 | double dx = fimg.ReadKey("DX");
|
---|
472 | double dy = fimg.ReadKey("DY");
|
---|
473 | double dz = fimg.ReadKey("DZ");
|
---|
474 | double pkzref = fimg.ReadKey("ZREFPK");
|
---|
475 | double zref = fimg.ReadKey("ZREF");
|
---|
476 | double kzref = fimg.ReadKey("KZREF");
|
---|
477 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
478 | setpointers(true);
|
---|
479 | init_fftw();
|
---|
480 | SetObservator(zref,kzref);
|
---|
481 | array_allocated_ = true;
|
---|
482 | compute_pk_redsh_ref_ = pkzref;
|
---|
483 | } catch (PThrowable & exc) {
|
---|
484 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
485 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
486 | return;
|
---|
487 | } catch (...) {
|
---|
488 | cout<<" some other exception was caught !"<<endl;
|
---|
489 | return;
|
---|
490 | }
|
---|
491 | }
|
---|
492 |
|
---|
493 | void GeneFluct3D::WritePPF(string cfname,bool write_real)
|
---|
494 | // On ecrit soit le TArray<r_8> ou le TArray<complex <r_8> >
|
---|
495 | {
|
---|
496 | cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
|
---|
497 | try {
|
---|
498 | R_.Info()["NX"] = (int_8)Nx_;
|
---|
499 | R_.Info()["NY"] = (int_8)Ny_;
|
---|
500 | R_.Info()["NZ"] = (int_8)Nz_;
|
---|
501 | R_.Info()["DX"] = (r_8)Dx_;
|
---|
502 | R_.Info()["DY"] = (r_8)Dy_;
|
---|
503 | R_.Info()["DZ"] = (r_8)Dz_;
|
---|
504 | R_.Info()["ZREFPK"] = (r_8)compute_pk_redsh_ref_;
|
---|
505 | R_.Info()["ZREF"] = (r_8)redsh_ref_;
|
---|
506 | R_.Info()["KZREF"] = (r_8)kredsh_ref_;
|
---|
507 | R_.Info()["HREF"] = (r_8)h_ref_;
|
---|
508 | R_.Info()["Growth"] = (r_8)Dref();
|
---|
509 | R_.Info()["GrowthPk"] = (r_8)DrefPk();
|
---|
510 | POutPersist pos(cfname.c_str());
|
---|
511 | if(write_real) pos << PPFNameTag("rgen") << R_;
|
---|
512 | else pos << PPFNameTag("pkgen") << T_;
|
---|
513 | } catch (PThrowable & exc) {
|
---|
514 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
515 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
516 | return;
|
---|
517 | } catch (...) {
|
---|
518 | cout<<" some other exception was caught !"<<endl;
|
---|
519 | return;
|
---|
520 | }
|
---|
521 | }
|
---|
522 |
|
---|
523 | void GeneFluct3D::ReadPPF(string cfname)
|
---|
524 | {
|
---|
525 | cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
|
---|
526 | try {
|
---|
527 | bool from_real = true;
|
---|
528 | PInPersist pis(cfname.c_str());
|
---|
529 | string name_tag_k = "pkgen";
|
---|
530 | bool found_tag_k = pis.GotoNameTag("pkgen");
|
---|
531 | if(found_tag_k) {
|
---|
532 | cout<<" ...reading spectrum into TArray<complex <r_8> >"<<endl;
|
---|
533 | pis >> PPFNameTag("pkgen") >> T_;
|
---|
534 | from_real = false;
|
---|
535 | } else {
|
---|
536 | cout<<" ...reading space into TArray<r_8>"<<endl;
|
---|
537 | pis >> PPFNameTag("rgen") >> R_;
|
---|
538 | }
|
---|
539 | setpointers(from_real); // a mettre ici pour relire les DVInfo
|
---|
540 | int_8 nx = R_.Info()["NX"];
|
---|
541 | int_8 ny = R_.Info()["NY"];
|
---|
542 | int_8 nz = R_.Info()["NZ"];
|
---|
543 | r_8 dx = R_.Info()["DX"];
|
---|
544 | r_8 dy = R_.Info()["DY"];
|
---|
545 | r_8 dz = R_.Info()["DZ"];
|
---|
546 | r_8 pkzref = R_.Info()["ZREFPK"];
|
---|
547 | r_8 zref = R_.Info()["ZREF"];
|
---|
548 | r_8 kzref = R_.Info()["KZREF"];
|
---|
549 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
550 | init_fftw();
|
---|
551 | SetObservator(zref,kzref);
|
---|
552 | array_allocated_ = true;
|
---|
553 | compute_pk_redsh_ref_ = pkzref;
|
---|
554 | } catch (PThrowable & exc) {
|
---|
555 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
556 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
557 | return;
|
---|
558 | } catch (...) {
|
---|
559 | cout<<" some other exception was caught !"<<endl;
|
---|
560 | return;
|
---|
561 | }
|
---|
562 | }
|
---|
563 |
|
---|
564 | void GeneFluct3D::WriteSlicePPF(string cfname)
|
---|
565 | // On ecrit 3 tranches du cube selon chaque axe
|
---|
566 | {
|
---|
567 | cout<<"--- GeneFluct3D::WriteSlicePPF: Writing Cube Slices "<<cfname<<endl;
|
---|
568 | try {
|
---|
569 |
|
---|
570 | POutPersist pos(cfname.c_str());
|
---|
571 | TMatrix<r_4> S;
|
---|
572 | char str[16];
|
---|
573 | long i,j,l;
|
---|
574 |
|
---|
575 | // Tranches en Z
|
---|
576 | for(int s=0;s<3;s++) {
|
---|
577 | S.ReSize(Nx_,Ny_);
|
---|
578 | if(s==0) l=0; else if(s==1) l=(Nz_+1)/2; else l=Nz_-1;
|
---|
579 | sprintf(str,"z%ld",l);
|
---|
580 | for(i=0;i<Nx_;i++) for(j=0;j<Ny_;j++) S(i,j)=data_[IndexR(i,j,l)];
|
---|
581 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
582 | }
|
---|
583 |
|
---|
584 | // Tranches en Y
|
---|
585 | for(int s=0;s<3;s++) {
|
---|
586 | S.ReSize(Nz_,Nx_);
|
---|
587 | if(s==0) j=0; else if(s==1) j=(Ny_+1)/2; else j=Ny_-1;
|
---|
588 | sprintf(str,"y%ld",j);
|
---|
589 | for(i=0;i<Nx_;i++) for(l=0;l<Nz_;l++) S(l,i)=data_[IndexR(i,j,l)];
|
---|
590 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
591 | }
|
---|
592 |
|
---|
593 | // Tranches en X
|
---|
594 | for(int s=0;s<3;s++) {
|
---|
595 | S.ReSize(Nz_,Ny_);
|
---|
596 | if(s==0) i=0; else if(s==1) i=(Nx_+1)/2; else i=Nx_-1;
|
---|
597 | sprintf(str,"x%ld",i);
|
---|
598 | for(j=0;j<Ny_;j++) for(l=0;l<Nz_;l++) S(l,j)=data_[IndexR(i,j,l)];
|
---|
599 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
600 | }
|
---|
601 |
|
---|
602 | } catch (PThrowable & exc) {
|
---|
603 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
604 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
605 | return;
|
---|
606 | } catch (...) {
|
---|
607 | cout<<" some other exception was caught !"<<endl;
|
---|
608 | return;
|
---|
609 | }
|
---|
610 | }
|
---|
611 |
|
---|
612 | //-------------------------------------------------------
|
---|
613 | void GeneFluct3D::NTupleCheck(POutPersist &pos,string ntname,unsigned long nent)
|
---|
614 | // Remplit le NTuple "ntname" avec "nent" valeurs du cube (reel ou complex) et l'ecrit dans "pos"
|
---|
615 | {
|
---|
616 | if(ntname.size()<=0 || nent==0) return;
|
---|
617 | int nvar = 0;
|
---|
618 | if(array_type==1) nvar = 3;
|
---|
619 | else if(array_type==2) nvar = 4;
|
---|
620 | else return;
|
---|
621 | const char *vname[4] = {"t","z","re","im"};
|
---|
622 | float xnt[4];
|
---|
623 | NTuple nt(nvar,vname);
|
---|
624 |
|
---|
625 | if(array_type==1) {
|
---|
626 | unsigned long nmod = Nx_*Ny_*Nz_/nent; if(nmod==0) nmod=1;
|
---|
627 | unsigned long n=0;
|
---|
628 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
629 | if(n==nmod) {
|
---|
630 | int_8 ip = IndexR(i,j,l);
|
---|
631 | xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=data_[ip];
|
---|
632 | nt.Fill(xnt);
|
---|
633 | n=0;
|
---|
634 | }
|
---|
635 | n++;
|
---|
636 | }
|
---|
637 | } else {
|
---|
638 | unsigned long nmod = Nx_*Ny_*NCz_/nent; if(nmod==0) nmod=1;
|
---|
639 | unsigned long n=0;
|
---|
640 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<NCz_;l++) {
|
---|
641 | if(n==nmod) {
|
---|
642 | xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=T_(l,j,i).real(); xnt[3]=T_(l,j,i).imag();
|
---|
643 | nt.Fill(xnt);
|
---|
644 | n=0;
|
---|
645 | }
|
---|
646 | n++;
|
---|
647 | }
|
---|
648 | }
|
---|
649 |
|
---|
650 | pos.PutObject(nt,ntname);
|
---|
651 | }
|
---|
652 |
|
---|
653 | //-------------------------------------------------------
|
---|
654 | void GeneFluct3D::Print(void)
|
---|
655 | {
|
---|
656 | cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
|
---|
657 | cout<<"Space Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<Nz_<<" ("<<NTz_<<") size="
|
---|
658 | <<NRtot_<<endl;
|
---|
659 | cout<<" Resol: dx="<<Dx_<<" dy="<<Dy_<<" dz="<<Dz_<<" Mpc"
|
---|
660 | <<", dVol="<<dVol_<<", Vol="<<Vol_<<" Mpc^3"<<endl;
|
---|
661 | cout<<"Fourier Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<NCz_<<endl;
|
---|
662 | cout<<" Resol: dkx="<<Dkx_<<" dky="<<Dky_<<" dkz="<<Dkz_<<" Mpc^-1"
|
---|
663 | <<", Dk3="<<Dk3_<<" Mpc^-3"<<endl;
|
---|
664 | cout<<" (2Pi/k: "<<2.*M_PI/Dkx_<<" "<<2.*M_PI/Dky_<<" "<<2.*M_PI/Dkz_<<" Mpc)"<<endl;
|
---|
665 | cout<<" Nyquist: kx="<<Knyqx_<<" ky="<<Knyqy_<<" kz="<<Knyqz_<<" Mpc^-1"
|
---|
666 | <<", Kmax="<<GetKmax()<<" Mpc^-1"<<endl;
|
---|
667 | cout<<" (2Pi/k: "<<2.*M_PI/Knyqx_<<" "<<2.*M_PI/Knyqy_<<" "<<2.*M_PI/Knyqz_<<" Mpc)"<<endl;
|
---|
668 | cout<<"Redshift "<<redsh_ref_<<" for z axe at k="<<kredsh_ref_<<endl;
|
---|
669 | }
|
---|
670 |
|
---|
671 | //-------------------------------------------------------
|
---|
672 | void GeneFluct3D::ComputeFourier0(PkSpectrumZ& pk_at_z)
|
---|
673 | // cf ComputeFourier() mais avec autre methode de realisation du spectre
|
---|
674 | // (attention on fait une fft pour realiser le spectre)
|
---|
675 | {
|
---|
676 | compute_pk_redsh_ref_ = pk_at_z.GetZ();
|
---|
677 | // --- realisation d'un tableau de tirage gaussiens
|
---|
678 | if(lp_>0) cout<<"--- ComputeFourier0 at z="<<compute_pk_redsh_ref_
|
---|
679 | <<": before gaussian filling ---"<<endl;
|
---|
680 | // On tient compte du pb de normalisation de FFTW3
|
---|
681 | double sntot = sqrt((double)NRtot_);
|
---|
682 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
683 | int_8 ip = IndexR(i,j,l);
|
---|
684 | data_[ip] = NorRand()/sntot;
|
---|
685 | }
|
---|
686 |
|
---|
687 | // --- realisation d'un tableau de tirage gaussiens
|
---|
688 | if(lp_>0) cout<<"...before fft real ---"<<endl;
|
---|
689 | GEN3D_FFTW_EXECUTE(pf_);
|
---|
690 |
|
---|
691 | // --- On remplit avec une realisation
|
---|
692 | if(lp_>0) cout<<"...before Fourier realization filling"<<endl;
|
---|
693 | T_(0,0,0) = complex<GEN3D_TYPE>(0.); // on coupe le continue et on l'initialise
|
---|
694 | long lmod = Nx_/20; if(lmod<1) lmod=1;
|
---|
695 | for(long i=0;i<Nx_;i++) {
|
---|
696 | double kx = AbsKx(i); kx *= kx;
|
---|
697 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<"\t nx-i="<<Nx_-i<<"\t kx="<<Kx(i)<<"\t pk="<<pk_at_z(kx)<<endl;
|
---|
698 | for(long j=0;j<Ny_;j++) {
|
---|
699 | double ky = AbsKy(j); ky *= ky;
|
---|
700 | for(long l=0;l<NCz_;l++) {
|
---|
701 | double kz = AbsKz(l); kz *= kz;
|
---|
702 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
703 | double k = sqrt(kx+ky+kz);
|
---|
704 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
705 | double pk = pk_at_z(k)/Vol_;
|
---|
706 | // ici pas de "/2" a cause de la remarque ci-dessus
|
---|
707 | T_(l,j,i) *= sqrt(pk);
|
---|
708 | }
|
---|
709 | }
|
---|
710 | }
|
---|
711 |
|
---|
712 | array_type = 2;
|
---|
713 |
|
---|
714 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
715 | double p = compute_power_carte();
|
---|
716 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
717 |
|
---|
718 | }
|
---|
719 |
|
---|
720 | //-------------------------------------------------------
|
---|
721 | void GeneFluct3D::ComputeFourier(PkSpectrumZ& pk_at_z)
|
---|
722 | // Calcule une realisation du spectre "pk_at_z"
|
---|
723 | // Attention: dans TArray le premier indice varie le + vite
|
---|
724 | // Explication normalisation: see Coles & Lucchin, Cosmology, p264-265
|
---|
725 | // FFTW3: on note N=Nx*Ny*Nz
|
---|
726 | // f --(FFT)--> F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f))
|
---|
727 | // sum(f(x_i)^2) = S
|
---|
728 | // sum(F(nu_i)^2) = S*N
|
---|
729 | // sum(fb(x_i)^2) = S*N^2
|
---|
730 | {
|
---|
731 | // --- RaZ du tableau
|
---|
732 | T_ = complex<GEN3D_TYPE>(0.);
|
---|
733 | compute_pk_redsh_ref_ = pk_at_z.GetZ();
|
---|
734 |
|
---|
735 | // --- On remplit avec une realisation
|
---|
736 | if(lp_>0) cout<<"--- ComputeFourier at z="<<compute_pk_redsh_ref_<<" ---"<<endl;
|
---|
737 | long lmod = Nx_/20; if(lmod<1) lmod=1;
|
---|
738 | for(long i=0;i<Nx_;i++) {
|
---|
739 | double kx = AbsKx(i); kx *= kx;
|
---|
740 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<"\t nx-i="<<Nx_-i<<"\t kx="<<Kx(i)<<"\t pk="<<pk_at_z(kx)<<endl;
|
---|
741 | for(long j=0;j<Ny_;j++) {
|
---|
742 | double ky = AbsKy(j); ky *= ky;
|
---|
743 | for(long l=0;l<NCz_;l++) {
|
---|
744 | double kz = AbsKz(l); kz *= kz;
|
---|
745 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
746 | double k = sqrt(kx+ky+kz);
|
---|
747 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
748 | double pk = pk_at_z(k)/Vol_;
|
---|
749 | // Explication de la division par 2: voir perandom.cc
|
---|
750 | // ou egalement Coles & Lucchin, Cosmology formula 13.7.2 p279
|
---|
751 | T_(l,j,i) = ComplexGaussianRand(sqrt(pk/2.));
|
---|
752 | }
|
---|
753 | }
|
---|
754 | }
|
---|
755 |
|
---|
756 | array_type = 2;
|
---|
757 | manage_coefficients(); // gros effet pour les spectres que l'on utilise !
|
---|
758 |
|
---|
759 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
760 | double p = compute_power_carte();
|
---|
761 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
762 |
|
---|
763 | }
|
---|
764 |
|
---|
765 | long GeneFluct3D::manage_coefficients(void)
|
---|
766 | // Take into account the real and complexe conjugate coefficients
|
---|
767 | // because we want a realization of a real data in real space
|
---|
768 | // On ecrit que: P(k_x,k_y,k_z) = P(-k_x,-k_y,-k_z)
|
---|
769 | // avec k_x = i, -k_x = N_x - i etc...
|
---|
770 | {
|
---|
771 | if(lp_>1) cout<<"...managing coefficients"<<endl;
|
---|
772 | check_array_alloc();
|
---|
773 |
|
---|
774 | // 1./ Le Continu et Nyquist sont reels
|
---|
775 | long nreal = 0;
|
---|
776 | for(long kk=0;kk<2;kk++) {
|
---|
777 | long k=0; // continu
|
---|
778 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
779 | for(long jj=0;jj<2;jj++) {
|
---|
780 | long j=0;
|
---|
781 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
782 | for(long ii=0;ii<2;ii++) {
|
---|
783 | long i=0;
|
---|
784 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
785 | int_8 ip = IndexC(i,j,k);
|
---|
786 | //cout<<"i="<<i<<" j="<<j<<" k="<<k<<" = ("<<fdata_[ip][0]<<","<<fdata_[ip][1]<<")"<<endl;
|
---|
787 | fdata_[ip][1] = 0.; fdata_[ip][0] *= M_SQRT2;
|
---|
788 | nreal++;
|
---|
789 | }
|
---|
790 | }
|
---|
791 | }
|
---|
792 | if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
|
---|
793 |
|
---|
794 | // 2./ Les elements complexe conjugues (tous dans le plan k=0,Nyquist)
|
---|
795 |
|
---|
796 | // a./ les lignes et colonnes du continu et de nyquist
|
---|
797 | long nconj1 = 0;
|
---|
798 | for(long kk=0;kk<2;kk++) {
|
---|
799 | long k=0; // continu
|
---|
800 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
801 | for(long jj=0;jj<2;jj++) { // selon j
|
---|
802 | long j=0;
|
---|
803 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
804 | for(long i=1;i<(Nx_+1)/2;i++) {
|
---|
805 | int_8 ip = IndexC(i,j,k);
|
---|
806 | int_8 ip1 = IndexC(Nx_-i,j,k);
|
---|
807 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
808 | nconj1++;
|
---|
809 | }
|
---|
810 | }
|
---|
811 | for(long ii=0;ii<2;ii++) {
|
---|
812 | long i=0;
|
---|
813 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
814 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
815 | int_8 ip = IndexC(i,j,k);
|
---|
816 | int_8 ip1 = IndexC(i,Ny_-j,k);
|
---|
817 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
818 | nconj1++;
|
---|
819 | }
|
---|
820 | }
|
---|
821 | }
|
---|
822 | if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
|
---|
823 |
|
---|
824 | // b./ les lignes et colonnes hors continu et de nyquist
|
---|
825 | long nconj2 = 0;
|
---|
826 | for(long kk=0;kk<2;kk++) {
|
---|
827 | long k=0; // continu
|
---|
828 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
829 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
830 | if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
|
---|
831 | for(long i=1;i<Nx_;i++) {
|
---|
832 | if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
|
---|
833 | int_8 ip = IndexC(i,j,k);
|
---|
834 | int_8 ip1 = IndexC(Nx_-i,Ny_-j,k);
|
---|
835 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
836 | nconj2++;
|
---|
837 | }
|
---|
838 | }
|
---|
839 | }
|
---|
840 | if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
|
---|
841 |
|
---|
842 | if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(Nx_*Ny_*NCz_-nconj1-nconj2)-nreal<<endl;
|
---|
843 |
|
---|
844 | return nreal+nconj1+nconj2;
|
---|
845 | }
|
---|
846 |
|
---|
847 | double GeneFluct3D::compute_power_carte(void)
|
---|
848 | // Calcul la puissance de la realisation du spectre Pk
|
---|
849 | {
|
---|
850 | check_array_alloc();
|
---|
851 |
|
---|
852 | double s2 = 0.;
|
---|
853 | for(long l=0;l<NCz_;l++)
|
---|
854 | for(long j=0;j<Ny_;j++)
|
---|
855 | for(long i=0;i<Nx_;i++) s2 += MODULE2(T_(l,j,i));
|
---|
856 |
|
---|
857 | double s20 = 0.;
|
---|
858 | for(long j=0;j<Ny_;j++)
|
---|
859 | for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
|
---|
860 |
|
---|
861 | double s2n = 0.;
|
---|
862 | if(Nz_%2==0)
|
---|
863 | for(long j=0;j<Ny_;j++)
|
---|
864 | for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
|
---|
865 |
|
---|
866 | return 2.*s2 -s20 -s2n;
|
---|
867 | }
|
---|
868 |
|
---|
869 | //-------------------------------------------------------------------
|
---|
870 | void GeneFluct3D::FilterByPixel(void)
|
---|
871 | // Filtrage par la fonction fenetre du pixel (parallelepipede)
|
---|
872 | // TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
|
---|
873 | // e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
|
---|
874 | // = 2/(k_x*dx) * sin(k_x*dx/2) * (idem y) * (idem z)
|
---|
875 | // Gestion divergence en 0: sin(y)/y = 1 - y^2/6*(1-y^2/20)
|
---|
876 | // avec y = k_x*dx/2
|
---|
877 | {
|
---|
878 | if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
|
---|
879 | check_array_alloc();
|
---|
880 |
|
---|
881 | for(long i=0;i<Nx_;i++) {
|
---|
882 | double kx = AbsKx(i) *Dx_/2;
|
---|
883 | double pk_x = pixelfilter(kx);
|
---|
884 | for(long j=0;j<Ny_;j++) {
|
---|
885 | double ky = AbsKy(j) *Dy_/2;
|
---|
886 | double pk_y = pixelfilter(ky);
|
---|
887 | for(long l=0;l<NCz_;l++) {
|
---|
888 | double kz = AbsKz(l) *Dz_/2;
|
---|
889 | double pk_z = pixelfilter(kz);
|
---|
890 | T_(l,j,i) *= pk_x*pk_y*pk_z;
|
---|
891 | }
|
---|
892 | }
|
---|
893 | }
|
---|
894 |
|
---|
895 | }
|
---|
896 |
|
---|
897 | //-------------------------------------------------------------------
|
---|
898 | void GeneFluct3D::ToVelLoS(void)
|
---|
899 | {
|
---|
900 | double zpk = compute_pk_redsh_ref_;
|
---|
901 | double dpsd = -cosmo_->H(zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
902 | if(lp_>0) cout<<"--- ToVelLoS --- at z="<<zpk<<", D'/D="<<dpsd
|
---|
903 | <<" (km/s)/Mpc (comp. for dz="<<good_dzinc_<<")"<<endl;
|
---|
904 | check_array_alloc();
|
---|
905 |
|
---|
906 | for(long i=0;i<Nx_;i++) {
|
---|
907 | double kx = Kx(i);
|
---|
908 | for(long j=0;j<Ny_;j++) {
|
---|
909 | double ky = Ky(j);
|
---|
910 | double kt2 = kx*kx + ky*ky;
|
---|
911 | for(long l=0;l<NCz_;l++) {
|
---|
912 | double kz = Kz(l);
|
---|
913 | double k2 = kt2 + kz*kz;
|
---|
914 | if(k2<=0.) continue;
|
---|
915 | T_(l,j,i) *= complex<double>(0.,dpsd*kz/k2);
|
---|
916 | }
|
---|
917 | }
|
---|
918 | }
|
---|
919 | }
|
---|
920 |
|
---|
921 | //-------------------------------------------------------------------
|
---|
922 | void GeneFluct3D::ApplyGrowthFactor(int type_evol)
|
---|
923 | // Apply Growth to real space
|
---|
924 | // Using the correspondance between redshift and los comoving distance
|
---|
925 | // describe in vector "zred_" "loscom_"
|
---|
926 | // type_evol = 1 : evolution avec la distance a l'observateur
|
---|
927 | // 2 : evolution avec la distance du plan Z
|
---|
928 | // (tous les pixels d'un plan Z sont mis au meme redshift z que celui du milieu)
|
---|
929 | {
|
---|
930 | if(lp_>0) cout<<"--- ApplyGrowthFactor: evol="<<type_evol<<endl;
|
---|
931 | check_array_alloc();
|
---|
932 |
|
---|
933 | if(growth_ == NULL) {
|
---|
934 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
|
---|
935 | cout<<bla<<endl; throw ParmError(bla);
|
---|
936 | }
|
---|
937 | if(type_evol<1 || type_evol>2) {
|
---|
938 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
|
---|
939 | cout<<bla<<endl; throw ParmError(bla);
|
---|
940 | }
|
---|
941 |
|
---|
942 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
943 |
|
---|
944 | for(long i=0;i<Nx_;i++) {
|
---|
945 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
946 | for(long j=0;j<Ny_;j++) {
|
---|
947 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
948 | for(long l=0;l<Nz_;l++) {
|
---|
949 | double dz = DZcom(l);
|
---|
950 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
951 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
952 | double z = interpinv(dz); // interpolation par morceau
|
---|
953 | double dzgr = (*growth_)(z);
|
---|
954 | int_8 ip = IndexR(i,j,l);
|
---|
955 | data_[ip] *= dzgr;
|
---|
956 | }
|
---|
957 | }
|
---|
958 | }
|
---|
959 |
|
---|
960 | }
|
---|
961 |
|
---|
962 | //-------------------------------------------------------------------
|
---|
963 | void GeneFluct3D::ApplyDerGrowthFactor(int type_evol)
|
---|
964 | // Apply Conformal derivative of Growth to real space for transverse velocity cube
|
---|
965 | {
|
---|
966 | if(lp_>0) cout<<"--- ApplyDerGrowthFactor: evol="<<type_evol<<endl;
|
---|
967 | check_array_alloc();
|
---|
968 |
|
---|
969 | if(growth_ == NULL) {
|
---|
970 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
|
---|
971 | cout<<bla<<endl; throw ParmError(bla);
|
---|
972 | }
|
---|
973 | if(type_evol<1 || type_evol>2) {
|
---|
974 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
|
---|
975 | cout<<bla<<endl; throw ParmError(bla);
|
---|
976 | }
|
---|
977 |
|
---|
978 | double zpk = compute_pk_redsh_ref_;
|
---|
979 | double dpsd_orig = - cosmo_->H(zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
980 |
|
---|
981 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
982 | InterpFunc interdpd(zredmin_dpsd_,zredmax_dpsd_,dpsdfrzred_);
|
---|
983 |
|
---|
984 | for(long i=0;i<Nx_;i++) {
|
---|
985 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
986 | for(long j=0;j<Ny_;j++) {
|
---|
987 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
988 | for(long l=0;l<Nz_;l++) {
|
---|
989 | double dz = DZcom(l);
|
---|
990 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
991 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
992 | double z = interpinv(dz); // interpolation par morceau
|
---|
993 | double dpsd = interdpd(z);
|
---|
994 | int_8 ip = IndexR(i,j,l);
|
---|
995 | data_[ip] *= dpsd / dpsd_orig;
|
---|
996 | }
|
---|
997 | }
|
---|
998 | }
|
---|
999 |
|
---|
1000 | }
|
---|
1001 |
|
---|
1002 | //-------------------------------------------------------------------
|
---|
1003 | void GeneFluct3D::ComputeReal(void)
|
---|
1004 | // Calcule une realisation dans l'espace reel
|
---|
1005 | {
|
---|
1006 | if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
|
---|
1007 | check_array_alloc();
|
---|
1008 |
|
---|
1009 | // On fait la FFT
|
---|
1010 | GEN3D_FFTW_EXECUTE(pb_);
|
---|
1011 | array_type = 1;
|
---|
1012 | }
|
---|
1013 |
|
---|
1014 | //-------------------------------------------------------------------
|
---|
1015 | void GeneFluct3D::ReComputeFourier(void)
|
---|
1016 | {
|
---|
1017 | if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
|
---|
1018 | check_array_alloc();
|
---|
1019 |
|
---|
1020 | // On fait la FFT
|
---|
1021 | GEN3D_FFTW_EXECUTE(pf_);
|
---|
1022 | array_type = 2;
|
---|
1023 |
|
---|
1024 | // On corrige du pb de la normalisation de FFTW3
|
---|
1025 | complex<r_8> v((r_8)NRtot_,0.);
|
---|
1026 | for(long i=0;i<Nx_;i++)
|
---|
1027 | for(long j=0;j<Ny_;j++)
|
---|
1028 | for(long l=0;l<NCz_;l++) T_(l,j,i) /= v;
|
---|
1029 | }
|
---|
1030 |
|
---|
1031 | //-------------------------------------------------------------------
|
---|
1032 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
|
---|
1033 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
1034 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
1035 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
1036 | {
|
---|
1037 | if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
|
---|
1038 | check_array_alloc();
|
---|
1039 |
|
---|
1040 | if(herr.NBins()<0) return -1;
|
---|
1041 | herr.Zero();
|
---|
1042 |
|
---|
1043 | // Attention a l'ordre
|
---|
1044 | for(long i=0;i<Nx_;i++) {
|
---|
1045 | double kx = AbsKx(i); kx *= kx;
|
---|
1046 | for(long j=0;j<Ny_;j++) {
|
---|
1047 | double ky = AbsKy(j); ky *= ky;
|
---|
1048 | for(long l=0;l<NCz_;l++) {
|
---|
1049 | double kz = AbsKz(l);
|
---|
1050 | double k = sqrt(kx+ky+kz*kz);
|
---|
1051 | double pk = MODULE2(T_(l,j,i));
|
---|
1052 | herr.Add(k,pk);
|
---|
1053 | }
|
---|
1054 | }
|
---|
1055 | }
|
---|
1056 | herr.ToVariance();
|
---|
1057 |
|
---|
1058 | // renormalize to directly compare to original spectrum
|
---|
1059 | double norm = Vol_;
|
---|
1060 | herr *= norm;
|
---|
1061 |
|
---|
1062 | return 0;
|
---|
1063 | }
|
---|
1064 |
|
---|
1065 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
|
---|
1066 | {
|
---|
1067 | if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
|
---|
1068 | check_array_alloc();
|
---|
1069 |
|
---|
1070 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
1071 | herr.Zero();
|
---|
1072 |
|
---|
1073 | // Attention a l'ordre
|
---|
1074 | for(long i=0;i<Nx_;i++) {
|
---|
1075 | double kx = AbsKx(i); kx *= kx;
|
---|
1076 | for(long j=0;j<Ny_;j++) {
|
---|
1077 | double ky = AbsKy(j); ky *= ky;
|
---|
1078 | double kt = sqrt(kx+ky);
|
---|
1079 | for(long l=0;l<NCz_;l++) {
|
---|
1080 | double kz = AbsKz(l);
|
---|
1081 | double pk = MODULE2(T_(l,j,i));
|
---|
1082 | herr.Add(kt,kz,pk);
|
---|
1083 | }
|
---|
1084 | }
|
---|
1085 | }
|
---|
1086 | herr.ToVariance();
|
---|
1087 |
|
---|
1088 | // renormalize to directly compare to original spectrum
|
---|
1089 | double norm = Vol_;
|
---|
1090 | herr *= norm;
|
---|
1091 |
|
---|
1092 | return 0;
|
---|
1093 | }
|
---|
1094 |
|
---|
1095 | //-------------------------------------------------------------------
|
---|
1096 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor)
|
---|
1097 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
1098 | // AVEC la soustraction du niveau de bruit et la correction par filterpixel()
|
---|
1099 | // Si on ne fait pas ca, alors on obtient un spectre non-isotrope!
|
---|
1100 | //
|
---|
1101 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
1102 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
1103 | {
|
---|
1104 | if(lp_>0) cout<<"--- ComputeSpectrum: sigma="<<sigma<<endl;
|
---|
1105 | check_array_alloc();
|
---|
1106 |
|
---|
1107 | if(sigma<=0.) sigma = 0.;
|
---|
1108 | double sigma2 = sigma*sigma / (double)NRtot_;
|
---|
1109 |
|
---|
1110 | if(herr.NBins()<0) return -1;
|
---|
1111 | herr.Zero();
|
---|
1112 |
|
---|
1113 | TVector<r_8> vfz(NCz_);
|
---|
1114 | if(pixcor) // kz = l*Dkz_
|
---|
1115 | for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
|
---|
1116 |
|
---|
1117 | // Attention a l'ordre
|
---|
1118 | for(long i=0;i<Nx_;i++) {
|
---|
1119 | double kx = AbsKx(i);
|
---|
1120 | double fx = (pixcor) ? pixelfilter(kx*Dx_/2): 1.;
|
---|
1121 | kx *= kx; fx *= fx;
|
---|
1122 | for(long j=0;j<Ny_;j++) {
|
---|
1123 | double ky = AbsKy(j);
|
---|
1124 | double fy = (pixcor) ? pixelfilter(ky*Dy_/2): 1.;
|
---|
1125 | ky *= ky; fy *= fy;
|
---|
1126 | for(long l=0;l<NCz_;l++) {
|
---|
1127 | double kz = AbsKz(l);
|
---|
1128 | double k = sqrt(kx+ky+kz*kz);
|
---|
1129 | double pk = MODULE2(T_(l,j,i)) - sigma2;
|
---|
1130 | double fz = (pixcor) ? vfz(l): 1.;
|
---|
1131 | double f = fx*fy*fz;
|
---|
1132 | if(f>0.) herr.Add(k,pk/f);
|
---|
1133 | }
|
---|
1134 | }
|
---|
1135 | }
|
---|
1136 | herr.ToVariance();
|
---|
1137 | for(int i=0;i<herr.NBins();i++) herr(i) += sigma2;
|
---|
1138 |
|
---|
1139 | // renormalize to directly compare to original spectrum
|
---|
1140 | double norm = Vol_;
|
---|
1141 | herr *= norm;
|
---|
1142 |
|
---|
1143 | return 0;
|
---|
1144 | }
|
---|
1145 |
|
---|
1146 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor)
|
---|
1147 | // AVEC la soustraction du niveau de bruit et la correction par filterpixel()
|
---|
1148 | {
|
---|
1149 | if(lp_>0) cout<<"--- ComputeSpectrum2D: sigma="<<sigma<<endl;
|
---|
1150 | check_array_alloc();
|
---|
1151 |
|
---|
1152 | if(sigma<=0.) sigma = 0.;
|
---|
1153 | double sigma2 = sigma*sigma / (double)NRtot_;
|
---|
1154 |
|
---|
1155 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
1156 | herr.Zero();
|
---|
1157 |
|
---|
1158 | TVector<r_8> vfz(NCz_);
|
---|
1159 | if(pixcor) // kz = l*Dkz_
|
---|
1160 | for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
|
---|
1161 |
|
---|
1162 | // Attention a l'ordre
|
---|
1163 | for(long i=0;i<Nx_;i++) {
|
---|
1164 | double kx = AbsKx(i);
|
---|
1165 | double fx = (pixcor) ? pixelfilter(kx*Dx_/2) : 1.;
|
---|
1166 | kx *= kx; fx *= fx;
|
---|
1167 | for(long j=0;j<Ny_;j++) {
|
---|
1168 | double ky = AbsKy(j);
|
---|
1169 | double fy = (pixcor) ? pixelfilter(ky*Dy_/2) : 1.;
|
---|
1170 | ky *= ky; fy *= fy;
|
---|
1171 | double kt = sqrt(kx+ky);
|
---|
1172 | for(long l=0;l<NCz_;l++) {
|
---|
1173 | double kz = AbsKz(l);
|
---|
1174 | double pk = MODULE2(T_(l,j,i)) - sigma2;
|
---|
1175 | double fz = (pixcor) ? vfz(l): 1.;
|
---|
1176 | double f = fx*fy*fz;
|
---|
1177 | if(f>0.) herr.Add(kt,kz,pk/f);
|
---|
1178 | }
|
---|
1179 | }
|
---|
1180 | }
|
---|
1181 | herr.ToVariance();
|
---|
1182 | for(int i=0;i<herr.NBinX();i++)
|
---|
1183 | for(int j=0;j<herr.NBinY();j++) herr(i,j) += sigma2;
|
---|
1184 |
|
---|
1185 | // renormalize to directly compare to original spectrum
|
---|
1186 | double norm = Vol_;
|
---|
1187 | herr *= norm;
|
---|
1188 |
|
---|
1189 | return 0;
|
---|
1190 | }
|
---|
1191 |
|
---|
1192 | //-------------------------------------------------------
|
---|
1193 | int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
|
---|
1194 | // Recompute MASS variance in spherical top-hat (rayon=R)
|
---|
1195 | // Par definition: SigmaR^2 = <(M-<M>)^2>/<M>^2
|
---|
1196 | // ou M = masse dans sphere de rayon R
|
---|
1197 | // --- ATTENTION: la variance calculee a une tres grande dispersion
|
---|
1198 | // (surtout si le volume du cube est petit). Pour verifier
|
---|
1199 | // que le sigmaR calcule par cette methode est en accord avec
|
---|
1200 | // le sigmaR en input, il faut faire plusieurs simulations (~100)
|
---|
1201 | // et regarder la moyenne des sigmaR reconstruits
|
---|
1202 | {
|
---|
1203 | if(lp_>0) cout<<"--- VarianceFrReal R="<<R<<endl;
|
---|
1204 | check_array_alloc();
|
---|
1205 |
|
---|
1206 | long dnx = long(R/Dx_)+1; if(dnx<=0) dnx = 1;
|
---|
1207 | long dny = long(R/Dy_)+1; if(dny<=0) dny = 1;
|
---|
1208 | long dnz = long(R/Dz_)+1; if(dnz<=0) dnz = 1;
|
---|
1209 | if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
|
---|
1210 |
|
---|
1211 | double sum=0., sum2=0., sn=0., r2 = R*R;
|
---|
1212 | int_8 nsum=0;
|
---|
1213 |
|
---|
1214 | for(long i=dnx;i<Nx_-dnx;i+=2*dnx) {
|
---|
1215 | for(long j=dny;j<Ny_-dny;j+=2*dny) {
|
---|
1216 | for(long l=dnz;l<Nz_-dnz;l+=2*dnz) {
|
---|
1217 | double m=0.; int_8 n=0;
|
---|
1218 | for(long ii=i-dnx;ii<=i+dnx;ii++) {
|
---|
1219 | double x = (ii-i)*Dx_; x *= x;
|
---|
1220 | for(long jj=j-dny;jj<=j+dny;jj++) {
|
---|
1221 | double y = (jj-j)*Dy_; y *= y;
|
---|
1222 | for(long ll=l-dnz;ll<=l+dnz;ll++) {
|
---|
1223 | double z = (ll-l)*Dz_; z *= z;
|
---|
1224 | if(x+y+z>r2) continue;
|
---|
1225 | int_8 ip = IndexR(ii,jj,ll);
|
---|
1226 | m += 1.+data_[ip]; // 1+drho/rho
|
---|
1227 | n++;
|
---|
1228 | }
|
---|
1229 | }
|
---|
1230 | }
|
---|
1231 | if(n>0) {sum += m; sum2 += m*m; nsum++; sn += n;}
|
---|
1232 | //cout<<i<<","<<j<<","<<l<<" n="<<n<<" m="<<m<<" sum="<<sum<<" sum2="<<sum2<<endl;
|
---|
1233 | }
|
---|
1234 | }
|
---|
1235 | }
|
---|
1236 |
|
---|
1237 | if(nsum<=1) {var=0.; return nsum;}
|
---|
1238 | sum /= nsum;
|
---|
1239 | sum2 = sum2/nsum - sum*sum;
|
---|
1240 | sn /= nsum;
|
---|
1241 | if(lp_>0) cout<<"...<n>="<<sn<<", nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
|
---|
1242 | var = sum2/(sum*sum); // <dM^2>/<M>^2
|
---|
1243 | if(lp_>0) cout<<"...sigmaR^2 = <(M-<M>)^2>/<M>^2 = "<<var
|
---|
1244 | <<" -> sigmaR = "<<sqrt(var)<<endl;
|
---|
1245 |
|
---|
1246 | return nsum;
|
---|
1247 | }
|
---|
1248 |
|
---|
1249 | //-------------------------------------------------------
|
---|
1250 | int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
|
---|
1251 | // number of pixels outside of ]vmin,vmax[ extremites exclues
|
---|
1252 | // -> vmin and vmax are considered as bad
|
---|
1253 | {
|
---|
1254 | check_array_alloc();
|
---|
1255 |
|
---|
1256 | int_8 nbad = 0;
|
---|
1257 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1258 | int_8 ip = IndexR(i,j,l);
|
---|
1259 | double v = data_[ip];
|
---|
1260 | if(v<=vmin || v>=vmax) nbad++;
|
---|
1261 | }
|
---|
1262 |
|
---|
1263 | if(lp_>0) cout<<"--- NumberOfBad "<<nbad<<" px out of ]"<<vmin<<","<<vmax
|
---|
1264 | <<"[ i.e. frac="<<nbad/(double)NRtot_<<endl;
|
---|
1265 | return nbad;
|
---|
1266 | }
|
---|
1267 |
|
---|
1268 | int_8 GeneFluct3D::MinMax(double& xmin,double& xmax,double vmin,double vmax)
|
---|
1269 | // Calcul des valeurs xmin et xmax dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
1270 | {
|
---|
1271 | bool tstval = (vmax>vmin)? true: false;
|
---|
1272 | if(lp_>0) {
|
---|
1273 | cout<<"--- MinMax";
|
---|
1274 | if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
|
---|
1275 | cout<<endl;
|
---|
1276 | }
|
---|
1277 | check_array_alloc();
|
---|
1278 |
|
---|
1279 | int_8 n = 0;
|
---|
1280 | xmin = xmax = data_[0];
|
---|
1281 |
|
---|
1282 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1283 | int_8 ip = IndexR(i,j,l);
|
---|
1284 | double x = data_[ip];
|
---|
1285 | if(tstval && (x<=vmin || x>=vmax)) continue;
|
---|
1286 | if(x<xmin) xmin = x;
|
---|
1287 | if(x>xmax) xmax = x;
|
---|
1288 | n++;
|
---|
1289 | }
|
---|
1290 |
|
---|
1291 | if(lp_>0) cout<<" n="<<n<<" min="<<xmin<<" max="<<xmax<<endl;
|
---|
1292 |
|
---|
1293 | return n;
|
---|
1294 | }
|
---|
1295 |
|
---|
1296 | int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax
|
---|
1297 | ,bool useout,double vout)
|
---|
1298 | // Calcul de mean,sigma2 dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
1299 | // useout = false: ne pas utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
1300 | // true : utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
1301 | // en remplacant leurs valeurs par "vout"
|
---|
1302 | {
|
---|
1303 | bool tstval = (vmax>vmin)? true: false;
|
---|
1304 | if(lp_>0) {
|
---|
1305 | cout<<"--- MeanSigma2";
|
---|
1306 | if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
|
---|
1307 | if(useout) cout<<", useout="<<useout<<" vout="<<vout;
|
---|
1308 | cout<<endl;
|
---|
1309 | }
|
---|
1310 | check_array_alloc();
|
---|
1311 |
|
---|
1312 | int_8 n = 0;
|
---|
1313 | rm = rs2 = 0.;
|
---|
1314 |
|
---|
1315 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1316 | int_8 ip = IndexR(i,j,l);
|
---|
1317 | double v = data_[ip];
|
---|
1318 | if(tstval) {
|
---|
1319 | if(v<=vmin || v>=vmax) {if(useout) v=vout; else continue;}
|
---|
1320 | }
|
---|
1321 | rm += v;
|
---|
1322 | rs2 += v*v;
|
---|
1323 | n++;
|
---|
1324 | }
|
---|
1325 |
|
---|
1326 | if(n>1) {
|
---|
1327 | rm /= (double)n;
|
---|
1328 | rs2 = rs2/(double)n - rm*rm;
|
---|
1329 | }
|
---|
1330 |
|
---|
1331 | if(lp_>0) cout<<" n="<<n<<" m="<<rm<<" s2="<<rs2<<" s="<<sqrt(fabs(rs2))<<endl;
|
---|
1332 |
|
---|
1333 | return n;
|
---|
1334 | }
|
---|
1335 |
|
---|
1336 | int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
|
---|
1337 | // set to "val0" if out of range ]vmin,vmax[ extremites exclues
|
---|
1338 | // cad set to "val0" if in [vmin,vmax] -> vmin and vmax are set to val0
|
---|
1339 | {
|
---|
1340 | check_array_alloc();
|
---|
1341 |
|
---|
1342 | int_8 nbad = 0;
|
---|
1343 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1344 | int_8 ip = IndexR(i,j,l);
|
---|
1345 | double v = data_[ip];
|
---|
1346 | if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
|
---|
1347 | }
|
---|
1348 |
|
---|
1349 | if(lp_>0) cout<<"--- SetToVal "<<nbad<<" px set to="<<val0
|
---|
1350 | <<" because out of range=]"<<vmin<<","<<vmax<<"["<<endl;
|
---|
1351 | return nbad;
|
---|
1352 | }
|
---|
1353 |
|
---|
1354 | void GeneFluct3D::ScaleOffset(double scalecube,double offsetcube)
|
---|
1355 | // Replace "V" by "scalecube * ( V + offsetcube )"
|
---|
1356 | {
|
---|
1357 | if(lp_>0) cout<<"--- ScaleCube scale="<<scalecube<<" offset="<<offsetcube<<endl;
|
---|
1358 |
|
---|
1359 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1360 | int_8 ip = IndexR(i,j,l);
|
---|
1361 | data_[ip] = scalecube * ( data_[ip] + offsetcube );
|
---|
1362 | }
|
---|
1363 |
|
---|
1364 | return;
|
---|
1365 | }
|
---|
1366 |
|
---|
1367 | //-------------------------------------------------------
|
---|
1368 | void GeneFluct3D::TurnFluct2Mass(void)
|
---|
1369 | // d_rho/rho -> Mass (add one!)
|
---|
1370 | {
|
---|
1371 | if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
|
---|
1372 | check_array_alloc();
|
---|
1373 |
|
---|
1374 |
|
---|
1375 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1376 | int_8 ip = IndexR(i,j,l);
|
---|
1377 | data_[ip] += 1.;
|
---|
1378 | }
|
---|
1379 | }
|
---|
1380 |
|
---|
1381 | double GeneFluct3D::TurnFluct2MeanNumber(double val_by_mpc3)
|
---|
1382 | // ATTENTION: la gestion des pixels<0 proposee ici induit une perte de variance
|
---|
1383 | // sur la carte, le spectre Pk reconstruit sera plus faible!
|
---|
1384 | // L'effet sera d'autant plus grand que le nombre de pixels<0 sera grand.
|
---|
1385 | {
|
---|
1386 | if(lp_>0) cout<<"--- TurnFluct2MeanNumber : "<<val_by_mpc3<<" quantity (gal or mass)/Mpc^3"<<endl;
|
---|
1387 |
|
---|
1388 | // First convert dRho/Rho into 1+dRho/Rho
|
---|
1389 | int_8 nball = 0; double sumall = 0., sumall2 = 0.;
|
---|
1390 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1391 | int_8 ip = IndexR(i,j,l);
|
---|
1392 | data_[ip] += 1.;
|
---|
1393 | nball++; sumall += data_[ip]; sumall2 += data_[ip]*data_[ip];
|
---|
1394 | }
|
---|
1395 | if(nball>2) {
|
---|
1396 | sumall /= (double)nball;
|
---|
1397 | sumall2 = sumall2/(double)nball - sumall*sumall;
|
---|
1398 | if(lp_>0) cout<<"1+dRho/Rho: mean="<<sumall<<" variance="<<sumall2
|
---|
1399 | <<" -> "<<sqrt(fabs(sumall2))<<endl;
|
---|
1400 | }
|
---|
1401 |
|
---|
1402 | // Find contribution for positive pixels
|
---|
1403 | int_8 nbpos = 0; double sumpos = 0. , sumpos2 = 0.;
|
---|
1404 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1405 | int_8 ip = IndexR(i,j,l);
|
---|
1406 | double v = data_[ip];
|
---|
1407 | if(data_[ip]>0.) {nbpos++; sumpos += v; sumpos2 += v*v;}
|
---|
1408 | }
|
---|
1409 | if(nbpos<1) {
|
---|
1410 | cout<<"TurnFluct2MeanNumber_Error: nbpos<1"<<endl;
|
---|
1411 | throw RangeCheckError("TurnFluct2MeanNumber_Error: nbpos<1");
|
---|
1412 | }
|
---|
1413 | sumpos2 = sumpos2/nball - sumpos*sumpos/(nball*nball);
|
---|
1414 | if(lp_>0)
|
---|
1415 | cout<<"1+dRho/Rho with v<0 set to zero: mean="<<sumpos/nball
|
---|
1416 | <<" variance="<<sumpos2<<" -> "<<sqrt(fabs(sumpos2))<<endl;
|
---|
1417 | cout<<"Sum of positive values: sumpos="<<sumpos
|
---|
1418 | <<" (n(v>0) = "<<nbpos<<" frac(v>0)="<<nbpos/(double)NRtot_<<")"<<endl;
|
---|
1419 |
|
---|
1420 | // - Mettre exactement val_by_mpc3*Vol galaxies (ou Msol) dans notre survey
|
---|
1421 | // - Uniquement dans les pixels de masse >0.
|
---|
1422 | // - Mise a zero des pixels <0
|
---|
1423 | double dn = val_by_mpc3 * Vol_ / sumpos;
|
---|
1424 | if(lp_>0) cout<<"...density move from "
|
---|
1425 | <<val_by_mpc3*dVol_<<" to "<<dn<<" / pixel"<<endl;
|
---|
1426 |
|
---|
1427 | double sum = 0.;
|
---|
1428 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1429 | int_8 ip = IndexR(i,j,l);
|
---|
1430 | if(data_[ip]<=0.) data_[ip] = 0.;
|
---|
1431 | else {
|
---|
1432 | data_[ip] *= dn;
|
---|
1433 | sum += data_[ip];
|
---|
1434 | }
|
---|
1435 | }
|
---|
1436 |
|
---|
1437 | if(lp_>0) cout<<"...quantity put into survey "<<sum<<" / "<<val_by_mpc3*Vol_<<endl;
|
---|
1438 |
|
---|
1439 | return sum;
|
---|
1440 | }
|
---|
1441 |
|
---|
1442 | double GeneFluct3D::ApplyPoisson(void)
|
---|
1443 | // do NOT treate negative or nul mass -> let it as it is
|
---|
1444 | {
|
---|
1445 | if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
|
---|
1446 | check_array_alloc();
|
---|
1447 |
|
---|
1448 | double sum = 0.;
|
---|
1449 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1450 | int_8 ip = IndexR(i,j,l);
|
---|
1451 | double v = data_[ip];
|
---|
1452 | if(v>0.) {
|
---|
1453 | uint_8 dn = PoissonRand(v,10.);
|
---|
1454 | data_[ip] = (double)dn;
|
---|
1455 | sum += (double)dn;
|
---|
1456 | }
|
---|
1457 | }
|
---|
1458 | if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
|
---|
1459 |
|
---|
1460 | return sum;
|
---|
1461 | }
|
---|
1462 |
|
---|
1463 | double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
|
---|
1464 | // do NOT treate negative or nul mass -> let it as it is
|
---|
1465 | // INPUT:
|
---|
1466 | // massdist : distribution de masse (m*dn/dm)
|
---|
1467 | // axeslog = false : retourne la masse
|
---|
1468 | // = true : retourne le log10(mass)
|
---|
1469 | // RETURN la masse totale
|
---|
1470 | {
|
---|
1471 | if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
|
---|
1472 | check_array_alloc();
|
---|
1473 |
|
---|
1474 | double sum = 0.;
|
---|
1475 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1476 | int_8 ip = IndexR(i,j,l);
|
---|
1477 | double v = data_[ip];
|
---|
1478 | if(v>0.) {
|
---|
1479 | long ngal = long(v+0.1);
|
---|
1480 | data_[ip] = 0.;
|
---|
1481 | for(long i=0;i<ngal;i++) {
|
---|
1482 | double m = massdist.RandomInterp(); // massdist.Random();
|
---|
1483 | if(axeslog) m = pow(10.,m);
|
---|
1484 | data_[ip] += m;
|
---|
1485 | }
|
---|
1486 | sum += data_[ip];
|
---|
1487 | }
|
---|
1488 | }
|
---|
1489 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
1490 |
|
---|
1491 | return sum;
|
---|
1492 | }
|
---|
1493 |
|
---|
1494 | double GeneFluct3D::TurnNGal2MassQuick(SchechterMassDist& schmdist)
|
---|
1495 | // idem TurnNGal2Mass mais beaucoup plus rapide
|
---|
1496 | {
|
---|
1497 | if(lp_>0) cout<<"--- TurnNGal2MassQuick ---"<<endl;
|
---|
1498 | check_array_alloc();
|
---|
1499 |
|
---|
1500 | double sum = 0.;
|
---|
1501 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1502 | int_8 ip = IndexR(i,j,l);
|
---|
1503 | double v = data_[ip];
|
---|
1504 | if(v>0.) {
|
---|
1505 | long ngal = long(v+0.1);
|
---|
1506 | data_[ip] = schmdist.TirMass(ngal);
|
---|
1507 | sum += data_[ip];
|
---|
1508 | }
|
---|
1509 | }
|
---|
1510 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
1511 |
|
---|
1512 | return sum;
|
---|
1513 | }
|
---|
1514 |
|
---|
1515 | void GeneFluct3D::AddNoise2Real(double snoise,int type_evol)
|
---|
1516 | // add noise to every pixels (meme les <=0 !)
|
---|
1517 | // type_evol = 0 : pas d'evolution de la puissance du bruit
|
---|
1518 | // 1 : evolution de la puissance du bruit avec la distance a l'observateur
|
---|
1519 | // 2 : evolution de la puissance du bruit avec la distance du plan Z
|
---|
1520 | // (tous les plans Z sont mis au meme redshift z de leur milieu)
|
---|
1521 | {
|
---|
1522 | if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<" evol="<<type_evol<<endl;
|
---|
1523 | check_array_alloc();
|
---|
1524 |
|
---|
1525 | if(type_evol<0) type_evol = 0;
|
---|
1526 | if(type_evol>2) {
|
---|
1527 | const char *bla = "GeneFluct3D::AddNoise2Real_Error: bad type_evol value";
|
---|
1528 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1529 | }
|
---|
1530 |
|
---|
1531 | vector<double> correction;
|
---|
1532 | InterpFunc *intercor = NULL;
|
---|
1533 |
|
---|
1534 | if(type_evol>0) {
|
---|
1535 | // Sigma_Noise(en mass) :
|
---|
1536 | // Slim ~ 1/sqrt(DNu) * sqrt(nlobe) en W/m^2Hz
|
---|
1537 | // Flim ~ sqrt(DNu) * sqrt(nlobe) en W/m^2
|
---|
1538 | // Mlim ~ sqrt(DNu) * (Dlum)^2 * sqrt(nlobe) en Msol
|
---|
1539 | // nlobe ~ 1/Dtrcom^2
|
---|
1540 | // Mlim ~ sqrt(DNu) * (Dlum)^2 / Dtrcom
|
---|
1541 | if(cosmo_ == NULL || redsh_ref_<0. || loscom2zred_.size()<1) {
|
---|
1542 | const char *bla = "GeneFluct3D::AddNoise2Real_Error: set Observator and Cosmology first";
|
---|
1543 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1544 | }
|
---|
1545 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
1546 | long nsz = loscom2zred_.size(), nszmod=((nsz>10)? nsz/10: 1);
|
---|
1547 | for(long i=0;i<nsz;i++) {
|
---|
1548 | double d = interpinv.X(i);
|
---|
1549 | double zred = interpinv(d);
|
---|
1550 | double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
|
---|
1551 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
1552 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
1553 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
1554 | double corr = sqrt(dnu/dnu_ref_) * pow(dlum/dlum_ref_,2.) * dtrc_ref_/dtrc;
|
---|
1555 | if(lp_>0 && (i==0 || i==nsz-1 || i%nszmod==0))
|
---|
1556 | cout<<"i="<<i<<" d="<<d<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
|
---|
1557 | <<" dtrc="<<dtrc<<" dlum="<<dlum<<" -> cor="<<corr<<endl;
|
---|
1558 | correction.push_back(corr);
|
---|
1559 | }
|
---|
1560 | intercor = new InterpFunc(loscom2zred_min_,loscom2zred_max_,correction);
|
---|
1561 | }
|
---|
1562 |
|
---|
1563 | double corrlim[2] = {1.,1.};
|
---|
1564 | for(long i=0;i<Nx_;i++) {
|
---|
1565 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
1566 | for(long j=0;j<Ny_;j++) {
|
---|
1567 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
1568 | for(long l=0;l<Nz_;l++) {
|
---|
1569 | double corr = 1.;
|
---|
1570 | if(type_evol>0) {
|
---|
1571 | double dz = DZcom(l);
|
---|
1572 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
1573 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
1574 | corr = (*intercor)(dz);
|
---|
1575 | if(corr<corrlim[0]) corrlim[0]=corr; else if(corr>corrlim[1]) corrlim[1]=corr;
|
---|
1576 | }
|
---|
1577 | int_8 ip = IndexR(i,j,l);
|
---|
1578 | data_[ip] += snoise*corr*NorRand();
|
---|
1579 | }
|
---|
1580 | }
|
---|
1581 | }
|
---|
1582 | if(type_evol>0)
|
---|
1583 | cout<<"correction factor range: ["<<corrlim[0]<<","<<corrlim[1]<<"]"<<endl;
|
---|
1584 |
|
---|
1585 | if(intercor!=NULL) delete intercor;
|
---|
1586 | }
|
---|
1587 |
|
---|
1588 | } // Fin namespace SOPHYA
|
---|
1589 |
|
---|
1590 |
|
---|
1591 |
|
---|
1592 |
|
---|
1593 | /*********************************************************************
|
---|
1594 | void GeneFluct3D::AddAGN(double lfjy,double lsigma,double powlaw)
|
---|
1595 | // Add AGN flux into simulation:
|
---|
1596 | // --- Procedure:
|
---|
1597 | // 1. lancer "cmvdefsurv" avec les parametres du survey
|
---|
1598 | // (au redshift de reference du survey)
|
---|
1599 | // et recuperer l'angle solide "angsol sr" du pixel elementaire
|
---|
1600 | // au centre du cube.
|
---|
1601 | // 2. lancer "cmvtstagn" pour cet angle solide -> cmvtstagn.ppf
|
---|
1602 | // 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne
|
---|
1603 | // cad une moyenne <log10(S)> et un sigma "sig"
|
---|
1604 | // Attention: la distribution n'est pas gaussienne les "mean,sigma"
|
---|
1605 | // de l'histo ne sont pas vraiment ce que l'on veut
|
---|
1606 | // --- Limitations actuelle du code:
|
---|
1607 | // . les AGN sont supposes evoluer avec la meme loi de puissance pour tout theta,phi
|
---|
1608 | // . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee
|
---|
1609 | // . la distribution est approximee a une gaussienne
|
---|
1610 | // ... C'est une approximation pour un observateur loin du centre du cube
|
---|
1611 | // et pour un cube peu epais / distance observateur
|
---|
1612 | // --- Parametres de la routine:
|
---|
1613 | // llfy : c'est le <log10(S)> du flux depose par les AGN
|
---|
1614 | // dans l'angle solide du pixel elementaire de reference du cube
|
---|
1615 | // lsigma : c'est le sigma de la distribution des log10(S)
|
---|
1616 | // powlaw : c'est la pente de la distribution cad que le flux "lmsol"
|
---|
1617 | // et considere comme le flux a 1.4GHz et qu'on suppose une loi
|
---|
1618 | // F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz)
|
---|
1619 | // - Comme on est en echelle log10():
|
---|
1620 | // on tire log10(Msol) + X
|
---|
1621 | // ou X est une realisation sur une gaussienne de variance "sig^2"
|
---|
1622 | // La masse realisee est donc: Msol*10^X
|
---|
1623 | // - Pas de probleme de pixel negatif car on a une multiplication!
|
---|
1624 | {
|
---|
1625 | if(lp_>0) cout<<"--- AddAGN: <log10(S Jy)> = "<<lfjy<<" , sigma = "<<lsigma<<endl;
|
---|
1626 | check_array_alloc();
|
---|
1627 |
|
---|
1628 | if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) {
|
---|
1629 | char *bla = "GeneFluct3D::AddAGN_Error: set Observator and Cosmology first";
|
---|
1630 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1631 | }
|
---|
1632 |
|
---|
1633 | // Le flux des AGN en Jy et en mass solaire
|
---|
1634 | double fagnref = pow(10.,lfjy)*(dnu_ref_*1.e9); // Jy.Hz = W/m^2
|
---|
1635 | double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlum_ref_); // Msol
|
---|
1636 | if(lp_>0)
|
---|
1637 | cout<<"Au pixel de ref: fagnref="<<fagnref
|
---|
1638 | <<" Jy.Hz (a 1.4GHz), magnref="<<magnref<<" Msol"<<endl;
|
---|
1639 |
|
---|
1640 | if(powlaw!=0.) {
|
---|
1641 | // F(nu) = F(1.4GHz)*(nu GHz/1.4 Ghz)^p = F(1.4GHz)*(1/(1+z))^p , car nu = 1.4 GHz/(1+z)
|
---|
1642 | magnref *= pow(1/(1.+redsh_ref_),powlaw);
|
---|
1643 | if(lp_>0) cout<<" powlaw="<<powlaw<<" -> change magnref to "<<magnref<<" Msol"<<endl;
|
---|
1644 | }
|
---|
1645 |
|
---|
1646 | // Les infos en fonction de l'indice "l" selon Oz
|
---|
1647 | vector<double> correction;
|
---|
1648 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
1649 | long nzmod = ((Nz_>10)?Nz_/10:1);
|
---|
1650 | for(long l=0;l<Nz_;l++) {
|
---|
1651 | double z = fabs(DZcom(l));
|
---|
1652 | double zred = interpinv(z);
|
---|
1653 | double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
|
---|
1654 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
1655 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
1656 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
1657 | // on a: Mass ~ DNu * Dlum^2 / Dtrcom^2
|
---|
1658 | double corr = dnu/dnu_ref_*pow(dtrc_ref_/dtrc*dlum/dlum_ref_,2.);
|
---|
1659 | // F(nu) = F(1.4GHz)*(nu GHz/1.4 Ghz)^p = F(1.4GHz)*(1/(1+z))^p , car nu = 1.4 GHz/(1+z)
|
---|
1660 | if(powlaw!=0.) corr *= pow((1.+redsh_ref_)/(1.+zred),powlaw);
|
---|
1661 | correction.push_back(corr);
|
---|
1662 | if(lp_>0 && (l==0 || l==Nz_-1 || l%nzmod==0)) {
|
---|
1663 | cout<<"l="<<l<<" z="<<z<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
|
---|
1664 | <<" dtrc="<<dtrc<<" dlum="<<dlum
|
---|
1665 | <<" -> cor="<<corr<<endl;
|
---|
1666 | }
|
---|
1667 | }
|
---|
1668 |
|
---|
1669 | double sum=0., sum2=0., nsum=0.;
|
---|
1670 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) {
|
---|
1671 | double a = lsigma*NorRand();
|
---|
1672 | a = magnref*pow(10.,a);
|
---|
1673 | // On met le meme tirage le long de Oz (indice k)
|
---|
1674 | for(long l=0;l<Nz_;l++) {
|
---|
1675 | int_8 ip = IndexR(i,j,l);
|
---|
1676 | data_[ip] += a*correction[l];
|
---|
1677 | }
|
---|
1678 | sum += a; sum2 += a*a; nsum += 1.;
|
---|
1679 | }
|
---|
1680 |
|
---|
1681 | if(lp_>0 && nsum>1.) {
|
---|
1682 | sum /= nsum;
|
---|
1683 | sum2 = sum2/nsum - sum*sum;
|
---|
1684 | cout<<"...Mean mass="<<sum<<" Msol , s^2="<<sum2<<" s="<<sqrt(fabs(sum2))<<endl;
|
---|
1685 | }
|
---|
1686 |
|
---|
1687 | }
|
---|
1688 | *********************************************************************/
|
---|