1 | #include "sopnamsp.h"
|
---|
2 | #include "machdefs.h"
|
---|
3 | #include <iostream>
|
---|
4 | #include <stdlib.h>
|
---|
5 | #include <stdio.h>
|
---|
6 | #include <string.h>
|
---|
7 | #include <math.h>
|
---|
8 | #include <unistd.h>
|
---|
9 |
|
---|
10 | #include "tarray.h"
|
---|
11 | #include "pexceptions.h"
|
---|
12 | #include "perandom.h"
|
---|
13 | #include "srandgen.h"
|
---|
14 |
|
---|
15 | #include "fabtcolread.h"
|
---|
16 | #include "fabtwriter.h"
|
---|
17 | #include "fioarr.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 | //#define FFTW_THREAD
|
---|
28 |
|
---|
29 | #define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag()))
|
---|
30 |
|
---|
31 | //-------------------------------------------------------
|
---|
32 | GeneFluct3D::GeneFluct3D(TArray< complex<r_8 > >& T)
|
---|
33 | : T_(T) , Nx_(0) , Ny_(0) , Nz_(0) , array_allocated_(false) , lp_(0)
|
---|
34 | , redshref_(-999.) , kredshref_(0.) , cosmo_(NULL) , growth_(NULL)
|
---|
35 | , loscom_ref_(-999.), loscom_min_(-999.), loscom_max_(-999.)
|
---|
36 | , loscom2zred_min_(0.) , loscom2zred_max_(0.)
|
---|
37 | {
|
---|
38 | xobs_[0] = xobs_[1] = xobs_[2] = 0.;
|
---|
39 | zred_.resize(0);
|
---|
40 | loscom_.resize(0);
|
---|
41 | loscom2zred_.resize(0);
|
---|
42 | SetNThread();
|
---|
43 | }
|
---|
44 |
|
---|
45 | GeneFluct3D::~GeneFluct3D(void)
|
---|
46 | {
|
---|
47 | fftw_destroy_plan(pf_);
|
---|
48 | fftw_destroy_plan(pb_);
|
---|
49 | #ifdef FFTW_THREAD
|
---|
50 | if(nthread_>0) fftw_cleanup_threads();
|
---|
51 | #endif
|
---|
52 | }
|
---|
53 |
|
---|
54 | //-------------------------------------------------------
|
---|
55 | void GeneFluct3D::SetSize(long nx,long ny,long nz,double dx,double dy,double dz)
|
---|
56 | {
|
---|
57 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
58 | setalloc();
|
---|
59 | setpointers(false);
|
---|
60 | init_fftw();
|
---|
61 | }
|
---|
62 |
|
---|
63 | void GeneFluct3D::SetObservator(double redshref,double kredshref)
|
---|
64 | // L'observateur est au redshift z=0
|
---|
65 | // est situe sur la "perpendiculaire" a la face x,y
|
---|
66 | // issue du centre de cette face
|
---|
67 | // Il faut positionner le cube sur l'axe des z cad des redshifts:
|
---|
68 | // redshref = redshift de reference
|
---|
69 | // Si redshref<0 alors redshref=0
|
---|
70 | // kredshref = indice (en double) correspondant a ce redshift
|
---|
71 | // Si kredshref<0 alors kredshref=0
|
---|
72 | // Exemple: redshref=1.5 kredshref=250.75
|
---|
73 | // -> Le pixel i=nx/2 j=ny/2 k=250.75 est au redshift 1.5
|
---|
74 | {
|
---|
75 | if(redshref<0.) redshref = 0.;
|
---|
76 | if(kredshref<0.) kredshref = 0.;
|
---|
77 | redshref_ = redshref;
|
---|
78 | kredshref_ = kredshref;
|
---|
79 | if(lp_>0)
|
---|
80 | cout<<"--- GeneFluct3D::SetObservator zref="<<redshref_<<" kref="<<kredshref_<<endl;
|
---|
81 | }
|
---|
82 |
|
---|
83 | void GeneFluct3D::SetCosmology(CosmoCalc& cosmo)
|
---|
84 | {
|
---|
85 | cosmo_ = &cosmo;
|
---|
86 | if(lp_>1) cosmo_->Print();
|
---|
87 | }
|
---|
88 |
|
---|
89 | void GeneFluct3D::SetGrowthFactor(GrowthFactor& growth)
|
---|
90 | {
|
---|
91 | growth_ = &growth;
|
---|
92 | }
|
---|
93 |
|
---|
94 | void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz)
|
---|
95 | {
|
---|
96 | if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<<nx<<","<<ny<<","<<nz
|
---|
97 | <<" D="<<dx<<","<<dy<<","<<dz<<endl;
|
---|
98 | if(nx<=0 || dx<=0.) {
|
---|
99 | char *bla = "GeneFluct3D::setsize_Error: bad value(s";
|
---|
100 | cout<<bla<<endl; throw ParmError(bla);
|
---|
101 | }
|
---|
102 |
|
---|
103 | // Les tailles des tableaux
|
---|
104 | Nx_ = nx;
|
---|
105 | Ny_ = ny; if(Ny_ <= 0) Ny_ = Nx_;
|
---|
106 | Nz_ = nz; if(Nz_ <= 0) Nz_ = Nx_;
|
---|
107 | N_.resize(0); N_.push_back(Nx_); N_.push_back(Ny_); N_.push_back(Nz_);
|
---|
108 | NRtot_ = Nx_*Ny_*Nz_; // nombre de pixels dans le survey
|
---|
109 | NCz_ = Nz_/2 +1;
|
---|
110 | NTz_ = 2*NCz_;
|
---|
111 |
|
---|
112 | // le pas dans l'espace (Mpc)
|
---|
113 | Dx_ = dx;
|
---|
114 | Dy_ = dy; if(Dy_ <= 0.) Dy_ = Dx_;
|
---|
115 | Dz_ = dz; if(Dz_ <= 0.) Dz_ = Dx_;
|
---|
116 | D_.resize(0); D_.push_back(Dx_); D_.push_back(Dy_); D_.push_back(Dz_);
|
---|
117 | dVol_ = Dx_*Dy_*Dz_;
|
---|
118 | Vol_ = (Nx_*Dx_)*(Ny_*Dy_)*(Nz_*Dz_);
|
---|
119 |
|
---|
120 | // Le pas dans l'espace de Fourier (Mpc^-1)
|
---|
121 | Dkx_ = 2.*M_PI/(Nx_*Dx_);
|
---|
122 | Dky_ = 2.*M_PI/(Ny_*Dy_);
|
---|
123 | Dkz_ = 2.*M_PI/(Nz_*Dz_);
|
---|
124 | Dk_.resize(0); Dk_.push_back(Dkx_); Dk_.push_back(Dky_); Dk_.push_back(Dkz_);
|
---|
125 | Dk3_ = Dkx_*Dky_*Dkz_;
|
---|
126 |
|
---|
127 | // La frequence de Nyquist en k (Mpc^-1)
|
---|
128 | Knyqx_ = M_PI/Dx_;
|
---|
129 | Knyqy_ = M_PI/Dy_;
|
---|
130 | Knyqz_ = M_PI/Dz_;
|
---|
131 | Knyq_.resize(0); Knyq_.push_back(Knyqx_); Knyq_.push_back(Knyqy_); Knyq_.push_back(Knyqz_);
|
---|
132 | }
|
---|
133 |
|
---|
134 | void GeneFluct3D::setalloc(void)
|
---|
135 | {
|
---|
136 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc ---"<<endl;
|
---|
137 | // Dimensionnement du tableau complex<r_8>
|
---|
138 | // ATTENTION: TArray adresse en memoire a l'envers du C
|
---|
139 | // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
|
---|
140 | sa_size_t SzK_[3] = {NCz_,Ny_,Nx_}; // a l'envers
|
---|
141 | try {
|
---|
142 | T_.ReSize(3,SzK_);
|
---|
143 | array_allocated_ = true;
|
---|
144 | if(lp_>1) cout<<" allocating: "<<T_.Size()*sizeof(complex<r_8>)/1.e6<<" Mo"<<endl;
|
---|
145 | } catch (...) {
|
---|
146 | cout<<"GeneFluct3D::setalloc_Error: Problem allocating T_"<<endl;
|
---|
147 | }
|
---|
148 | T_.SetMemoryMapping(BaseArray::CMemoryMapping);
|
---|
149 | }
|
---|
150 |
|
---|
151 | void GeneFluct3D::setpointers(bool from_real)
|
---|
152 | {
|
---|
153 | if(lp_>1) cout<<"--- GeneFluct3D::setpointers ---"<<endl;
|
---|
154 | if(from_real) T_ = ArrCastR2C(R_);
|
---|
155 | else R_ = ArrCastC2R(T_);
|
---|
156 | // On remplit les pointeurs
|
---|
157 | fdata_ = (fftw_complex *) (&T_(0,0,0));
|
---|
158 | data_ = (double *) (&R_(0,0,0));
|
---|
159 | }
|
---|
160 |
|
---|
161 | void GeneFluct3D::check_array_alloc(void)
|
---|
162 | // Pour tester si le tableau T_ est alloue
|
---|
163 | {
|
---|
164 | if(array_allocated_) return;
|
---|
165 | char bla[90];
|
---|
166 | sprintf(bla,"GeneFluct3D::check_array_alloc_Error: array is not allocated");
|
---|
167 | cout<<bla<<endl; throw ParmError(bla);
|
---|
168 | }
|
---|
169 |
|
---|
170 | void GeneFluct3D::init_fftw(void)
|
---|
171 | {
|
---|
172 | // --- Initialisation de fftw3 (attention data est sur-ecrit a l'init)
|
---|
173 | if(lp_>1) cout<<"--- GeneFluct3D::init_fftw ---"<<endl;
|
---|
174 | #ifdef FFTW_THREAD
|
---|
175 | if(nthread_>0) {
|
---|
176 | cout<<"...Computing with "<<nthread_<<" threads"<<endl;
|
---|
177 | fftw_init_threads();
|
---|
178 | fftw_plan_with_nthreads(nthread_);
|
---|
179 | }
|
---|
180 | #endif
|
---|
181 | if(lp_>1) cout<<"...forward plan"<<endl;
|
---|
182 | pf_ = fftw_plan_dft_r2c_3d(Nx_,Ny_,Nz_,data_,fdata_,FFTW_ESTIMATE);
|
---|
183 | if(lp_>1) cout<<"...backward plan"<<endl;
|
---|
184 | pb_ = fftw_plan_dft_c2r_3d(Nx_,Ny_,Nz_,fdata_,data_,FFTW_ESTIMATE);
|
---|
185 | }
|
---|
186 |
|
---|
187 | //-------------------------------------------------------
|
---|
188 | long GeneFluct3D::LosComRedshift(double zinc,long npoints)
|
---|
189 | // Given a position of the cube relative to the observer
|
---|
190 | // and a cosmology
|
---|
191 | // (SetObservator() and SetCosmology() should have been called !)
|
---|
192 | // This routine filled:
|
---|
193 | // the vector "zred_" of scanned redshift (by zinc increments)
|
---|
194 | // the vector "loscom_" of corresponding los comoving distance
|
---|
195 | // -- Input:
|
---|
196 | // zinc : redshift increment for computation
|
---|
197 | // npoints : number of points required for inverting loscom -> zred
|
---|
198 | //
|
---|
199 | {
|
---|
200 | if(lp_>0) cout<<"--- LosComRedshift: zinc="<<zinc<<" , npoints="<<npoints<<endl;
|
---|
201 |
|
---|
202 | if(cosmo_ == NULL || redshref_<0.) {
|
---|
203 | char *bla = "GeneFluct3D::LosComRedshift_Error: set Observator and Cosmology first";
|
---|
204 | cout<<bla<<endl; throw ParmError(bla);
|
---|
205 | }
|
---|
206 |
|
---|
207 | // On calcule les coordonnees de l'observateur dans le repere du cube
|
---|
208 | // cad dans le repere ou l'origine est au centre du pixel i=j=l=0.
|
---|
209 | // L'observateur est sur un axe centre sur le milieu de la face Oxy
|
---|
210 | double loscom_ref_ = cosmo_->Dloscom(redshref_);
|
---|
211 | xobs_[0] = Nx_/2.*Dx_;
|
---|
212 | xobs_[1] = Ny_/2.*Dy_;
|
---|
213 | xobs_[2] = kredshref_*Dz_ - loscom_ref_;
|
---|
214 |
|
---|
215 | // L'observateur est-il dans le cube?
|
---|
216 | bool obs_in_cube = false;
|
---|
217 | if(xobs_[2]>=0. && xobs_[2]<=Nz_*Dz_) obs_in_cube = true;
|
---|
218 |
|
---|
219 | // Find MINIMUM los com distance to the observer:
|
---|
220 | // c'est le centre de la face a k=0
|
---|
221 | // (ou zero si l'observateur est dans le cube)
|
---|
222 | loscom_min_ = 0.;
|
---|
223 | if(!obs_in_cube) loscom_min_ = -xobs_[2];
|
---|
224 |
|
---|
225 | // Find MAXIMUM los com distance to the observer:
|
---|
226 | // ou que soit positionne l'observateur, la distance
|
---|
227 | // maximal est sur un des coins du cube
|
---|
228 | loscom_max_ = 0.;
|
---|
229 | for(long i=0;i<=1;i++) {
|
---|
230 | double dx2 = xobs_[0] - i*Nx_*Dx_; dx2 *= dx2;
|
---|
231 | for(long j=0;j<=1;j++) {
|
---|
232 | double dy2 = xobs_[1] - j*Ny_*Dy_; dy2 *= dy2;
|
---|
233 | for(long k=0;k<=1;k++) {
|
---|
234 | double dz2 = xobs_[2] - k*Nz_*Dz_; dz2 *= dz2;
|
---|
235 | dz2 = sqrt(dx2+dy2+dz2);
|
---|
236 | if(dz2>loscom_max_) loscom_max_ = dz2;
|
---|
237 | }
|
---|
238 | }
|
---|
239 | }
|
---|
240 | if(lp_>0) {
|
---|
241 | cout<<"...zref="<<redshref_<<" kzref="<<kredshref_<<" losref="<<loscom_ref_<<" Mpc\n"
|
---|
242 | <<" xobs="<<xobs_[0]<<" , "<<xobs_[1]<<" , "<<xobs_[2]<<" Mpc "
|
---|
243 | <<" in_cube="<<obs_in_cube
|
---|
244 | <<" loscom_min="<<loscom_min_<<" loscom_max="<<loscom_max_<<" Mpc "<<endl;
|
---|
245 | }
|
---|
246 |
|
---|
247 | // Fill the corresponding vectors for loscom and zred
|
---|
248 | if(zinc<=0.) zinc = 0.01;
|
---|
249 | for(double z=0.; ; z+=zinc) {
|
---|
250 | double dlc = cosmo_->Dloscom(z);
|
---|
251 | if(dlc<loscom_min_) {zred_.resize(0); loscom_.resize(0);}
|
---|
252 | zred_.push_back(z);
|
---|
253 | loscom_.push_back(dlc);
|
---|
254 | z += zinc;
|
---|
255 | if(dlc>loscom_max_) break; // on sort apres avoir stoque un dlc>dlcmax
|
---|
256 | }
|
---|
257 |
|
---|
258 | if(lp_>0) {
|
---|
259 | long n = zred_.size();
|
---|
260 | cout<<"...zred/loscom tables[zinc="<<zinc<<"]: n="<<n;
|
---|
261 | if(n>0) cout<<" z="<<zred_[0]<<" -> d="<<loscom_[0];
|
---|
262 | if(n>1) cout<<" , z="<<zred_[n-1]<<" -> d="<<loscom_[n-1];
|
---|
263 | cout<<endl;
|
---|
264 | }
|
---|
265 |
|
---|
266 | // Compute the parameters and tables needed for inversion loscom->zred
|
---|
267 | if(npoints<3) npoints = zred_.size();
|
---|
268 | InverseFunc invfun(zred_,loscom_);
|
---|
269 | invfun.ComputeParab(npoints,loscom2zred_);
|
---|
270 | loscom2zred_min_ = invfun.YMin();
|
---|
271 | loscom2zred_max_ = invfun.YMax();
|
---|
272 |
|
---|
273 | if(lp_>0) {
|
---|
274 | long n = loscom2zred_.size();
|
---|
275 | cout<<"...loscom -> zred[npoints="<<npoints<<"]: n="<<n
|
---|
276 | <<" los_min="<<loscom2zred_min_
|
---|
277 | <<" los_max="<<loscom2zred_max_
|
---|
278 | <<" -> zred=[";
|
---|
279 | if(n>0) cout<<loscom2zred_[0];
|
---|
280 | cout<<",";
|
---|
281 | if(n>1) cout<<loscom2zred_[n-1];
|
---|
282 | cout<<"]"<<endl;
|
---|
283 | if(lp_>1 && n>0)
|
---|
284 | for(int i=0;i<n;i++)
|
---|
285 | if(i==0 || abs(i-n/2)<2 || i==n-1)
|
---|
286 | cout<<" "<<i<<" "<<loscom2zred_[i]<<endl;
|
---|
287 | }
|
---|
288 |
|
---|
289 | return zred_.size();
|
---|
290 | }
|
---|
291 |
|
---|
292 | //-------------------------------------------------------
|
---|
293 | void GeneFluct3D::WriteFits(string cfname,int bitpix)
|
---|
294 | {
|
---|
295 | cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
|
---|
296 | try {
|
---|
297 | FitsImg3DWriter fwrt(cfname.c_str(),bitpix,5);
|
---|
298 | fwrt.WriteKey("NX",Nx_," axe transverse 1");
|
---|
299 | fwrt.WriteKey("NY",Ny_," axe transverse 2");
|
---|
300 | fwrt.WriteKey("NZ",Nz_," axe longitudinal (redshift)");
|
---|
301 | fwrt.WriteKey("DX",Dx_," Mpc");
|
---|
302 | fwrt.WriteKey("DY",Dy_," Mpc");
|
---|
303 | fwrt.WriteKey("DZ",Dz_," Mpc");
|
---|
304 | fwrt.WriteKey("DKX",Dkx_," Mpc^-1");
|
---|
305 | fwrt.WriteKey("DKY",Dky_," Mpc^-1");
|
---|
306 | fwrt.WriteKey("DKZ",Dkz_," Mpc^-1");
|
---|
307 | fwrt.WriteKey("ZREF",redshref_," reference redshift");
|
---|
308 | fwrt.WriteKey("KZREF",kredshref_," reference redshift on z axe");
|
---|
309 | fwrt.Write(R_);
|
---|
310 | } catch (PThrowable & exc) {
|
---|
311 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
312 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
313 | return;
|
---|
314 | } catch (...) {
|
---|
315 | cout<<" some other exception was caught !"<<endl;
|
---|
316 | return;
|
---|
317 | }
|
---|
318 | }
|
---|
319 |
|
---|
320 | void GeneFluct3D::ReadFits(string cfname)
|
---|
321 | {
|
---|
322 | cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
|
---|
323 | try {
|
---|
324 | FitsImg3DRead fimg(cfname.c_str(),0,5);
|
---|
325 | fimg.Read(R_);
|
---|
326 | long nx = fimg.ReadKeyL("NX");
|
---|
327 | long ny = fimg.ReadKeyL("NY");
|
---|
328 | long nz = fimg.ReadKeyL("NZ");
|
---|
329 | double dx = fimg.ReadKey("DX");
|
---|
330 | double dy = fimg.ReadKey("DY");
|
---|
331 | double dz = fimg.ReadKey("DZ");
|
---|
332 | double zref = fimg.ReadKey("ZREF");
|
---|
333 | double kzref = fimg.ReadKey("KZREF");
|
---|
334 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
335 | setpointers(true);
|
---|
336 | init_fftw();
|
---|
337 | SetObservator(zref,kzref);
|
---|
338 | } catch (PThrowable & exc) {
|
---|
339 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
340 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
341 | return;
|
---|
342 | } catch (...) {
|
---|
343 | cout<<" some other exception was caught !"<<endl;
|
---|
344 | return;
|
---|
345 | }
|
---|
346 | }
|
---|
347 |
|
---|
348 | void GeneFluct3D::WritePPF(string cfname,bool write_real)
|
---|
349 | // On ecrit soit le TArray<r_8> ou le TArray<complex <r_8> >
|
---|
350 | {
|
---|
351 | cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
|
---|
352 | try {
|
---|
353 | R_.Info()["NX"] = (int_8)Nx_;
|
---|
354 | R_.Info()["NY"] = (int_8)Ny_;
|
---|
355 | R_.Info()["NZ"] = (int_8)Nz_;
|
---|
356 | R_.Info()["DX"] = (r_8)Dx_;
|
---|
357 | R_.Info()["DY"] = (r_8)Dy_;
|
---|
358 | R_.Info()["DZ"] = (r_8)Dz_;
|
---|
359 | R_.Info()["ZREF"] = (r_8)redshref_;
|
---|
360 | R_.Info()["KZREF"] = (r_8)kredshref_;
|
---|
361 | POutPersist pos(cfname.c_str());
|
---|
362 | if(write_real) pos << PPFNameTag("rgen") << R_;
|
---|
363 | else pos << PPFNameTag("pkgen") << T_;
|
---|
364 | } catch (PThrowable & exc) {
|
---|
365 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
366 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
367 | return;
|
---|
368 | } catch (...) {
|
---|
369 | cout<<" some other exception was caught !"<<endl;
|
---|
370 | return;
|
---|
371 | }
|
---|
372 | }
|
---|
373 |
|
---|
374 | void GeneFluct3D::ReadPPF(string cfname)
|
---|
375 | {
|
---|
376 | cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
|
---|
377 | try {
|
---|
378 | bool from_real = true;
|
---|
379 | PInPersist pis(cfname.c_str());
|
---|
380 | string name_tag_k = "pkgen";
|
---|
381 | bool found_tag_k = pis.GotoNameTag("pkgen");
|
---|
382 | if(found_tag_k) {
|
---|
383 | cout<<" ...reading spectrun into TArray<complex <r_8> >"<<endl;
|
---|
384 | pis >> PPFNameTag("pkgen") >> T_;
|
---|
385 | from_real = false;
|
---|
386 | } else {
|
---|
387 | cout<<" ...reading space into TArray<r_8>"<<endl;
|
---|
388 | pis >> PPFNameTag("rgen") >> R_;
|
---|
389 | }
|
---|
390 | setpointers(from_real); // a mettre ici pour relire les DVInfo
|
---|
391 | int_8 nx = R_.Info()["NX"];
|
---|
392 | int_8 ny = R_.Info()["NY"];
|
---|
393 | int_8 nz = R_.Info()["NZ"];
|
---|
394 | r_8 dx = R_.Info()["DX"];
|
---|
395 | r_8 dy = R_.Info()["DY"];
|
---|
396 | r_8 dz = R_.Info()["DZ"];
|
---|
397 | r_8 zref = R_.Info()["ZREF"];
|
---|
398 | r_8 kzref = R_.Info()["KZREF"];
|
---|
399 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
400 | init_fftw();
|
---|
401 | SetObservator(zref,kzref);
|
---|
402 | } catch (PThrowable & exc) {
|
---|
403 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
404 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
405 | return;
|
---|
406 | } catch (...) {
|
---|
407 | cout<<" some other exception was caught !"<<endl;
|
---|
408 | return;
|
---|
409 | }
|
---|
410 | }
|
---|
411 |
|
---|
412 | //-------------------------------------------------------
|
---|
413 | void GeneFluct3D::Print(void)
|
---|
414 | {
|
---|
415 | cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
|
---|
416 | cout<<"Space Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<Nz_<<" ("<<NTz_<<") size="
|
---|
417 | <<NRtot_<<endl;
|
---|
418 | cout<<" Resol: dx="<<Dx_<<" dy="<<Dy_<<" dz="<<Dz_<<" Mpc"
|
---|
419 | <<", dVol="<<dVol_<<", Vol="<<Vol_<<" Mpc^3"<<endl;
|
---|
420 | cout<<"Fourier Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<NCz_<<endl;
|
---|
421 | cout<<" Resol: dkx="<<Dkx_<<" dky="<<Dky_<<" dkz="<<Dkz_<<" Mpc^-1"
|
---|
422 | <<", Dk3="<<Dk3_<<" Mpc^-3"<<endl;
|
---|
423 | cout<<" (2Pi/k: "<<2.*M_PI/Dkx_<<" "<<2.*M_PI/Dky_<<" "<<2.*M_PI/Dkz_<<" Mpc)"<<endl;
|
---|
424 | cout<<" Nyquist: kx="<<Knyqx_<<" ky="<<Knyqy_<<" kz="<<Knyqz_<<" Mpc^-1"
|
---|
425 | <<", Kmax="<<GetKmax()<<" Mpc^-1"<<endl;
|
---|
426 | cout<<" (2Pi/k: "<<2.*M_PI/Knyqx_<<" "<<2.*M_PI/Knyqy_<<" "<<2.*M_PI/Knyqz_<<" Mpc)"<<endl;
|
---|
427 | cout<<"Redshift "<<redshref_<<" for z axe at k="<<kredshref_<<endl;
|
---|
428 | }
|
---|
429 |
|
---|
430 | //-------------------------------------------------------
|
---|
431 | void GeneFluct3D::ComputeFourier0(GenericFunc& pk_at_z)
|
---|
432 | // cf ComputeFourier() mais avec autre methode de realisation du spectre
|
---|
433 | // (attention on fait une fft pour realiser le spectre)
|
---|
434 | {
|
---|
435 |
|
---|
436 | // --- realisation d'un tableau de tirage gaussiens
|
---|
437 | if(lp_>0) cout<<"--- ComputeFourier0: before gaussian filling ---"<<endl;
|
---|
438 | // On tient compte du pb de normalisation de FFTW3
|
---|
439 | double sntot = sqrt((double)NRtot_);
|
---|
440 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
441 | int_8 ip = IndexR(i,j,l);
|
---|
442 | data_[ip] = NorRand()/sntot;
|
---|
443 | }
|
---|
444 |
|
---|
445 | // --- realisation d'un tableau de tirage gaussiens
|
---|
446 | if(lp_>0) cout<<"...before fft real ---"<<endl;
|
---|
447 | fftw_execute(pf_);
|
---|
448 |
|
---|
449 | // --- On remplit avec une realisation
|
---|
450 | if(lp_>0) cout<<"...before Fourier realization filling"<<endl;
|
---|
451 | T_(0,0,0) = complex<r_8>(0.); // on coupe le continue et on l'initialise
|
---|
452 | long lmod = Nx_/10; if(lmod<1) lmod=1;
|
---|
453 | for(long i=0;i<Nx_;i++) {
|
---|
454 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
455 | double kx = ii*Dkx_; kx *= kx;
|
---|
456 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<" ii="<<ii<<endl;
|
---|
457 | for(long j=0;j<Ny_;j++) {
|
---|
458 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
459 | double ky = jj*Dky_; ky *= ky;
|
---|
460 | for(long l=0;l<NCz_;l++) {
|
---|
461 | double kz = l*Dkz_; kz *= kz;
|
---|
462 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
463 | double k = sqrt(kx+ky+kz);
|
---|
464 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
465 | double pk = pk_at_z(k)/Vol_;
|
---|
466 | // ici pas de "/2" a cause de la remarque ci-dessus
|
---|
467 | T_(l,j,i) *= sqrt(pk);
|
---|
468 | }
|
---|
469 | }
|
---|
470 | }
|
---|
471 |
|
---|
472 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
473 | double p = compute_power_carte();
|
---|
474 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
475 |
|
---|
476 | }
|
---|
477 |
|
---|
478 | //-------------------------------------------------------
|
---|
479 | void GeneFluct3D::ComputeFourier(GenericFunc& pk_at_z)
|
---|
480 | // Calcule une realisation du spectre "pk_at_z"
|
---|
481 | // Attention: dans TArray le premier indice varie le + vite
|
---|
482 | // Explication normalisation: see Coles & Lucchin, Cosmology, p264-265
|
---|
483 | // FFTW3: on note N=Nx*Ny*Nz
|
---|
484 | // f --(FFT)--> F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f))
|
---|
485 | // sum(f(x_i)^2) = S
|
---|
486 | // sum(F(nu_i)^2) = S*N
|
---|
487 | // sum(fb(x_i)^2) = S*N^2
|
---|
488 | {
|
---|
489 | // --- RaZ du tableau
|
---|
490 | T_ = complex<r_8>(0.);
|
---|
491 |
|
---|
492 | // --- On remplit avec une realisation
|
---|
493 | if(lp_>0) cout<<"--- ComputeFourier ---"<<endl;
|
---|
494 | long lmod = Nx_/10; if(lmod<1) lmod=1;
|
---|
495 | for(long i=0;i<Nx_;i++) {
|
---|
496 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
497 | double kx = ii*Dkx_; kx *= kx;
|
---|
498 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<" ii="<<ii<<endl;
|
---|
499 | for(long j=0;j<Ny_;j++) {
|
---|
500 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
501 | double ky = jj*Dky_; ky *= ky;
|
---|
502 | for(long l=0;l<NCz_;l++) {
|
---|
503 | double kz = l*Dkz_; kz *= kz;
|
---|
504 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
505 | double k = sqrt(kx+ky+kz);
|
---|
506 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
507 | double pk = pk_at_z(k)/Vol_;
|
---|
508 | // Explication de la division par 2: voir perandom.cc
|
---|
509 | // ou egalement Coles & Lucchin, Cosmology formula 13.7.2 p279
|
---|
510 | T_(l,j,i) = ComplexGaussRan(sqrt(pk/2.));
|
---|
511 | }
|
---|
512 | }
|
---|
513 | }
|
---|
514 |
|
---|
515 | manage_coefficients(); // gros effet pour les spectres que l'on utilise !
|
---|
516 |
|
---|
517 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
518 | double p = compute_power_carte();
|
---|
519 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
520 |
|
---|
521 | }
|
---|
522 |
|
---|
523 | long GeneFluct3D::manage_coefficients(void)
|
---|
524 | // Take into account the real and complexe conjugate coefficients
|
---|
525 | // because we want a realization of a real data in real space
|
---|
526 | {
|
---|
527 | if(lp_>1) cout<<"...managing coefficients"<<endl;
|
---|
528 | check_array_alloc();
|
---|
529 |
|
---|
530 | // 1./ Le Continu et Nyquist sont reels
|
---|
531 | long nreal = 0;
|
---|
532 | for(long kk=0;kk<2;kk++) {
|
---|
533 | long k=0; // continu
|
---|
534 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
535 | for(long jj=0;jj<2;jj++) {
|
---|
536 | long j=0;
|
---|
537 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
538 | for(long ii=0;ii<2;ii++) {
|
---|
539 | long i=0;
|
---|
540 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
541 | int_8 ip = IndexC(i,j,k);
|
---|
542 | //cout<<"i="<<i<<" j="<<j<<" k="<<k<<" = ("<<fdata_[ip][0]<<","<<fdata_[ip][1]<<")"<<endl;
|
---|
543 | fdata_[ip][1] = 0.; fdata_[ip][0] *= M_SQRT2;
|
---|
544 | nreal++;
|
---|
545 | }
|
---|
546 | }
|
---|
547 | }
|
---|
548 | if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
|
---|
549 |
|
---|
550 | // 2./ Les elements complexe conjugues (tous dans le plan k=0,Nyquist)
|
---|
551 |
|
---|
552 | // a./ les lignes et colonnes du continu et de nyquist
|
---|
553 | long nconj1 = 0;
|
---|
554 | for(long kk=0;kk<2;kk++) {
|
---|
555 | long k=0; // continu
|
---|
556 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
557 | for(long jj=0;jj<2;jj++) { // selon j
|
---|
558 | long j=0;
|
---|
559 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
560 | for(long i=1;i<(Nx_+1)/2;i++) {
|
---|
561 | int_8 ip = IndexC(i,j,k);
|
---|
562 | int_8 ip1 = IndexC(Nx_-i,j,k);
|
---|
563 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
564 | nconj1++;
|
---|
565 | }
|
---|
566 | }
|
---|
567 | for(long ii=0;ii<2;ii++) {
|
---|
568 | long i=0;
|
---|
569 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
570 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
571 | int_8 ip = IndexC(i,j,k);
|
---|
572 | int_8 ip1 = IndexC(i,Ny_-j,k);
|
---|
573 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
574 | nconj1++;
|
---|
575 | }
|
---|
576 | }
|
---|
577 | }
|
---|
578 | if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
|
---|
579 |
|
---|
580 | // b./ les lignes et colonnes hors continu et de nyquist
|
---|
581 | long nconj2 = 0;
|
---|
582 | for(long kk=0;kk<2;kk++) {
|
---|
583 | long k=0; // continu
|
---|
584 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
585 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
586 | if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
|
---|
587 | for(long i=1;i<Nx_;i++) {
|
---|
588 | if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
|
---|
589 | int_8 ip = IndexC(i,j,k);
|
---|
590 | int_8 ip1 = IndexC(Nx_-i,Ny_-j,k);
|
---|
591 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
592 | nconj2++;
|
---|
593 | }
|
---|
594 | }
|
---|
595 | }
|
---|
596 | if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
|
---|
597 |
|
---|
598 | if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(Nx_*Ny_*NCz_-nconj1-nconj2)-8<<endl;
|
---|
599 |
|
---|
600 | return nreal+nconj1+nconj2;
|
---|
601 | }
|
---|
602 |
|
---|
603 | double GeneFluct3D::compute_power_carte(void)
|
---|
604 | // Calcul la puissance de la realisation du spectre Pk
|
---|
605 | {
|
---|
606 | check_array_alloc();
|
---|
607 |
|
---|
608 | double s2 = 0.;
|
---|
609 | for(long l=0;l<NCz_;l++)
|
---|
610 | for(long j=0;j<Ny_;j++)
|
---|
611 | for(long i=0;i<Nx_;i++) s2 += MODULE2(T_(l,j,i));
|
---|
612 |
|
---|
613 | double s20 = 0.;
|
---|
614 | for(long j=0;j<Ny_;j++)
|
---|
615 | for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
|
---|
616 |
|
---|
617 | double s2n = 0.;
|
---|
618 | if(Nz_%2==0)
|
---|
619 | for(long j=0;j<Ny_;j++)
|
---|
620 | for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
|
---|
621 |
|
---|
622 | return 2.*s2 -s20 -s2n;
|
---|
623 | }
|
---|
624 |
|
---|
625 | //-------------------------------------------------------------------
|
---|
626 | void GeneFluct3D::FilterByPixel(void)
|
---|
627 | // Filtrage par la fonction fenetre du pixel (parallelepipede)
|
---|
628 | // TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
|
---|
629 | // e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
|
---|
630 | // = 2/(k_x*dx) * sin(k_x*dx/2) * (idem y) * (idem z)
|
---|
631 | // Gestion divergence en 0: sin(y)/y = 1 - y^2/6*(1-y^2/20)
|
---|
632 | // avec y = k_x*dx/2
|
---|
633 | {
|
---|
634 | if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
|
---|
635 | check_array_alloc();
|
---|
636 |
|
---|
637 | for(long i=0;i<Nx_;i++) {
|
---|
638 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
639 | double kx = ii*Dkx_ *Dx_/2;
|
---|
640 | double pk_x = pixelfilter(kx);
|
---|
641 | for(long j=0;j<Ny_;j++) {
|
---|
642 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
643 | double ky = jj*Dky_ *Dy_/2;
|
---|
644 | double pk_y = pixelfilter(ky);
|
---|
645 | for(long l=0;l<NCz_;l++) {
|
---|
646 | double kz = l*Dkz_ *Dz_/2;
|
---|
647 | double pk_z = pixelfilter(kz);
|
---|
648 | T_(l,j,i) *= pk_x*pk_y*pk_z;
|
---|
649 | }
|
---|
650 | }
|
---|
651 | }
|
---|
652 |
|
---|
653 | }
|
---|
654 |
|
---|
655 | //-------------------------------------------------------------------
|
---|
656 | void GeneFluct3D::ApplyGrowthFactor(void)
|
---|
657 | // Apply Growth to real space
|
---|
658 | // Using the correspondance between redshift and los comoving distance
|
---|
659 | // describe in vector "zred_" "loscom_"
|
---|
660 | {
|
---|
661 | if(lp_>0) cout<<"--- ApplyGrowthFactor ---"<<endl;
|
---|
662 | check_array_alloc();
|
---|
663 |
|
---|
664 | if(growth_ == NULL) {
|
---|
665 | char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
|
---|
666 | cout<<bla<<endl; throw ParmError(bla);
|
---|
667 | }
|
---|
668 |
|
---|
669 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
670 | unsigned short ok;
|
---|
671 |
|
---|
672 | //CHECK: Histo hgr(0.9*zred_[0],1.1*zred_[n-1],1000);
|
---|
673 | for(long i=0;i<Nx_;i++) {
|
---|
674 | double dx2 = xobs_[0] - i*Dx_; dx2 *= dx2;
|
---|
675 | for(long j=0;j<Ny_;j++) {
|
---|
676 | double dy2 = xobs_[1] - j*Dy_; dy2 *= dy2;
|
---|
677 | for(long l=0;l<Nz_;l++) {
|
---|
678 | double dz2 = xobs_[2] - l*Dz_; dz2 *= dz2;
|
---|
679 | dz2 = sqrt(dx2+dy2+dz2);
|
---|
680 | double z = interpinv(dz2);
|
---|
681 | //CHECK: hgr.Add(z);
|
---|
682 | double dzgr = (*growth_)(z); // interpolation par morceau
|
---|
683 | //double dzgr = growth_->Linear(z,ok); // interpolation lineaire
|
---|
684 | //double dzgr = growth_->Parab(z,ok); // interpolation parabolique
|
---|
685 | int_8 ip = IndexR(i,j,l);
|
---|
686 | data_[ip] *= dzgr;
|
---|
687 | }
|
---|
688 | }
|
---|
689 | }
|
---|
690 |
|
---|
691 | //CHECK: {POutPersist pos("applygrowth.ppf"); string tag="hgr"; pos.PutObject(hgr,tag);}
|
---|
692 |
|
---|
693 | }
|
---|
694 |
|
---|
695 | //-------------------------------------------------------------------
|
---|
696 | void GeneFluct3D::ComputeReal(void)
|
---|
697 | // Calcule une realisation dans l'espace reel
|
---|
698 | {
|
---|
699 | if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
|
---|
700 | check_array_alloc();
|
---|
701 |
|
---|
702 | // On fait la FFT
|
---|
703 | fftw_execute(pb_);
|
---|
704 | }
|
---|
705 |
|
---|
706 | //-------------------------------------------------------------------
|
---|
707 | void GeneFluct3D::ReComputeFourier(void)
|
---|
708 | {
|
---|
709 | if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
|
---|
710 | check_array_alloc();
|
---|
711 |
|
---|
712 | // On fait la FFT
|
---|
713 | fftw_execute(pf_);
|
---|
714 | // On corrige du pb de la normalisation de FFTW3
|
---|
715 | double v = (double)NRtot_;
|
---|
716 | for(long i=0;i<Nx_;i++)
|
---|
717 | for(long j=0;j<Ny_;j++)
|
---|
718 | for(long l=0;l<NCz_;l++) T_(l,j,i) /= complex<r_8>(v);
|
---|
719 |
|
---|
720 | }
|
---|
721 |
|
---|
722 | //-------------------------------------------------------------------
|
---|
723 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
|
---|
724 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
725 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
726 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
727 | {
|
---|
728 | if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
|
---|
729 | check_array_alloc();
|
---|
730 |
|
---|
731 | if(herr.NBins()<0) return -1;
|
---|
732 | herr.Zero();
|
---|
733 |
|
---|
734 | // Attention a l'ordre
|
---|
735 | for(long i=0;i<Nx_;i++) {
|
---|
736 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
737 | double kx = ii*Dkx_; kx *= kx;
|
---|
738 | for(long j=0;j<Ny_;j++) {
|
---|
739 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
740 | double ky = jj*Dky_; ky *= ky;
|
---|
741 | for(long l=0;l<NCz_;l++) {
|
---|
742 | double kz = l*Dkz_; kz *= kz;
|
---|
743 | double k = sqrt(kx+ky+kz);
|
---|
744 | double pk = MODULE2(T_(l,j,i));
|
---|
745 | herr.Add(k,pk);
|
---|
746 | }
|
---|
747 | }
|
---|
748 | }
|
---|
749 | herr.ToVariance();
|
---|
750 |
|
---|
751 | // renormalize to directly compare to original spectrum
|
---|
752 | double norm = Vol_;
|
---|
753 | herr *= norm;
|
---|
754 |
|
---|
755 | return 0;
|
---|
756 | }
|
---|
757 |
|
---|
758 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
|
---|
759 | {
|
---|
760 | if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
|
---|
761 | check_array_alloc();
|
---|
762 |
|
---|
763 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
764 | herr.Zero();
|
---|
765 |
|
---|
766 | // Attention a l'ordre
|
---|
767 | for(long i=0;i<Nx_;i++) {
|
---|
768 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
769 | double kx = ii*Dkx_; kx *= kx;
|
---|
770 | for(long j=0;j<Ny_;j++) {
|
---|
771 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
772 | double ky = jj*Dky_; ky *= ky;
|
---|
773 | double kt = sqrt(kx+ky);
|
---|
774 | for(long l=0;l<NCz_;l++) {
|
---|
775 | double kz = l*Dkz_;
|
---|
776 | double pk = MODULE2(T_(l,j,i));
|
---|
777 | herr.Add(kt,kz,pk);
|
---|
778 | }
|
---|
779 | }
|
---|
780 | }
|
---|
781 | herr.ToVariance();
|
---|
782 |
|
---|
783 | // renormalize to directly compare to original spectrum
|
---|
784 | double norm = Vol_;
|
---|
785 | herr *= norm;
|
---|
786 |
|
---|
787 | return 0;
|
---|
788 | }
|
---|
789 |
|
---|
790 | //-------------------------------------------------------
|
---|
791 | int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
|
---|
792 | // Recompute MASS variance in spherical top-hat (rayon=R)
|
---|
793 | {
|
---|
794 | if(lp_>0) cout<<"--- VarianceFrReal ---"<<endl;
|
---|
795 | check_array_alloc();
|
---|
796 |
|
---|
797 | long dnx = long(R/Dx_+0.5); if(dnx<=0) dnx = 1;
|
---|
798 | long dny = long(R/Dy_+0.5); if(dny<=0) dny = 1;
|
---|
799 | long dnz = long(R/Dz_+0.5); if(dnz<=0) dnz = 1;
|
---|
800 | if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
|
---|
801 |
|
---|
802 | double sum=0., sum2=0., r2 = R*R; int_8 nsum=0;
|
---|
803 |
|
---|
804 | for(long i=dnx;i<Nx_-dnx;i+=dnx) {
|
---|
805 | for(long j=dny;j<Ny_-dny;j+=dny) {
|
---|
806 | for(long l=dnz;l<Nz_-dnz;l+=dnz) {
|
---|
807 | double s=0.; int_8 n=0;
|
---|
808 | for(long ii=i-dnx;ii<=i+dnx;ii++) {
|
---|
809 | double x = (ii-i)*Dx_; x *= x;
|
---|
810 | for(long jj=j-dny;jj<=j+dny;jj++) {
|
---|
811 | double y = (jj-j)*Dy_; y *= y;
|
---|
812 | for(long ll=l-dnz;ll<=l+dnz;ll++) {
|
---|
813 | double z = (ll-l)*Dz_; z *= z;
|
---|
814 | if(x+y+z>r2) continue;
|
---|
815 | int_8 ip = IndexR(ii,jj,ll);
|
---|
816 | s += 1.+data_[ip];
|
---|
817 | n++;
|
---|
818 | }
|
---|
819 | }
|
---|
820 | }
|
---|
821 | if(n>0) {sum += s; sum2 += s*s; nsum++;}
|
---|
822 | //cout<<i<<","<<j<<","<<l<<" n="<<n<<" s="<<s<<" sum="<<sum<<" sum2="<<sum2<<endl;
|
---|
823 | }
|
---|
824 | }
|
---|
825 | }
|
---|
826 |
|
---|
827 | if(nsum<=1) {var=0.; return nsum;}
|
---|
828 |
|
---|
829 | sum /= nsum;
|
---|
830 | sum2 = sum2/nsum - sum*sum;
|
---|
831 | if(lp_>0) cout<<"VarianceFrReal: nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
|
---|
832 |
|
---|
833 | var = sum2/(sum*sum); // <dM>^2/<M>^2
|
---|
834 | if(lp_>0) cout<<"sigmaR^2="<<var<<" -> "<<sqrt(var)<<endl;
|
---|
835 |
|
---|
836 | return nsum;
|
---|
837 | }
|
---|
838 |
|
---|
839 | //-------------------------------------------------------
|
---|
840 | int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
|
---|
841 | // number of pixels outside of ]vmin,vmax[ extremites exclues
|
---|
842 | // -> vmin and vmax are considered as bad
|
---|
843 | {
|
---|
844 | check_array_alloc();
|
---|
845 |
|
---|
846 | int_8 nbad = 0;
|
---|
847 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
848 | int_8 ip = IndexR(i,j,l);
|
---|
849 | double v = data_[ip];
|
---|
850 | if(v<=vmin || v>=vmax) nbad++;
|
---|
851 | }
|
---|
852 |
|
---|
853 | return nbad;
|
---|
854 | }
|
---|
855 |
|
---|
856 | int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax
|
---|
857 | ,bool useout,double vout)
|
---|
858 | // Calcul de mean,sigma2 dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
859 | // useout = false: ne pas utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
860 | // true : utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
861 | // en remplacant leurs valeurs par "vout"
|
---|
862 | {
|
---|
863 | bool tstval = (vmax>vmin)? true: false;
|
---|
864 | if(lp_>0) {
|
---|
865 | cout<<"--- MeanSigma2: ";
|
---|
866 | if(tstval) cout<<"range=]"<<vmin<<","<<vmax<<"[";
|
---|
867 | if(useout) cout<<", useout="<<useout<<" vout="<<vout;
|
---|
868 | cout<<endl;
|
---|
869 | }
|
---|
870 | check_array_alloc();
|
---|
871 |
|
---|
872 | int_8 n = 0;
|
---|
873 | rm = rs2 = 0.;
|
---|
874 |
|
---|
875 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
876 | int_8 ip = IndexR(i,j,l);
|
---|
877 | double v = data_[ip];
|
---|
878 | if(tstval) {
|
---|
879 | if(v<=vmin || v>=vmax) {if(useout) v=vout; else continue;}
|
---|
880 | }
|
---|
881 | rm += v;
|
---|
882 | rs2 += v*v;
|
---|
883 | n++;
|
---|
884 | }
|
---|
885 |
|
---|
886 | if(n>1) {
|
---|
887 | rm /= (double)n;
|
---|
888 | rs2 = rs2/(double)n - rm*rm;
|
---|
889 | }
|
---|
890 |
|
---|
891 | if(lp_>0) cout<<" n="<<n<<" m="<<rm<<" s2="<<rs2<<" s="<<sqrt(fabs(rs2))<<endl;
|
---|
892 |
|
---|
893 | return n;
|
---|
894 | }
|
---|
895 |
|
---|
896 | int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
|
---|
897 | // set to "val0" if out of range ]vmin,vmax[ extremites exclues
|
---|
898 | // cad set to "val0" if in [vmin,vmax] -> vmin and vmax are set to val0
|
---|
899 | {
|
---|
900 | if(lp_>0) cout<<"--- SetToVal set to="<<val0
|
---|
901 | <<" when in range=["<<vmin<<","<<vmax<<"]"<<endl;
|
---|
902 | check_array_alloc();
|
---|
903 |
|
---|
904 | int_8 nbad = 0;
|
---|
905 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
906 | int_8 ip = IndexR(i,j,l);
|
---|
907 | double v = data_[ip];
|
---|
908 | if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
|
---|
909 | }
|
---|
910 |
|
---|
911 | return nbad;
|
---|
912 | }
|
---|
913 |
|
---|
914 | //-------------------------------------------------------
|
---|
915 | void GeneFluct3D::TurnFluct2Mass(void)
|
---|
916 | // d_rho/rho -> Mass (add one!)
|
---|
917 | {
|
---|
918 | if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
|
---|
919 | check_array_alloc();
|
---|
920 |
|
---|
921 |
|
---|
922 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
923 | int_8 ip = IndexR(i,j,l);
|
---|
924 | data_[ip] += 1.;
|
---|
925 | }
|
---|
926 | }
|
---|
927 |
|
---|
928 | double GeneFluct3D::TurnMass2MeanNumber(double n_by_mpc3)
|
---|
929 | // do NOT treate negative or nul values
|
---|
930 | {
|
---|
931 | if(lp_>0) cout<<"--- TurnMass2MeanNumber ---"<<endl;
|
---|
932 |
|
---|
933 | double m,s2;
|
---|
934 | int_8 ngood = MeanSigma2(m,s2,0.,1e+200);
|
---|
935 | if(lp_>0) cout<<"...ngood="<<ngood
|
---|
936 | <<" m="<<m<<" s2="<<s2<<" -> "<<sqrt(s2)<<endl;
|
---|
937 |
|
---|
938 | // On doit mettre n*Vol galaxies dans notre survey
|
---|
939 | // On en met uniquement dans les pixels de masse >0.
|
---|
940 | // On NE met PAS a zero les pixels <0
|
---|
941 | // On renormalise sur les pixels>0 pour qu'on ait n*Vol galaxies
|
---|
942 | // comme on ne prend que les pixels >0, on doit normaliser
|
---|
943 | // a la moyenne de <1+d_rho/rho> sur ces pixels
|
---|
944 | // (rappel sur tout les pixels <1+d_rho/rho>=1)
|
---|
945 | double dn = n_by_mpc3*Vol_/m /(double)ngood; // nb de gal a mettre ds 1 px
|
---|
946 | if(lp_>0) cout<<"...galaxy density move from "
|
---|
947 | <<n_by_mpc3*Vol_/double(NRtot_)<<" to "<<dn<<" / pixel"<<endl;
|
---|
948 | double sum = 0.;
|
---|
949 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
950 | int_8 ip = IndexR(i,j,l);
|
---|
951 | data_[ip] *= dn; // par coherence on multiplie aussi les <=0
|
---|
952 | if(data_[ip]>0.) sum += data_[ip]; // mais on ne les compte pas
|
---|
953 | }
|
---|
954 | if(lp_>0) cout<<sum<<"...galaxies put into survey / "<<n_by_mpc3*Vol_<<endl;
|
---|
955 |
|
---|
956 | return sum;
|
---|
957 | }
|
---|
958 |
|
---|
959 | double GeneFluct3D::ApplyPoisson(void)
|
---|
960 | // do NOT treate negative or nul mass -> let it as it is
|
---|
961 | {
|
---|
962 | if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
|
---|
963 | check_array_alloc();
|
---|
964 |
|
---|
965 | double sum = 0.;
|
---|
966 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
967 | int_8 ip = IndexR(i,j,l);
|
---|
968 | double v = data_[ip];
|
---|
969 | if(v>0.) {
|
---|
970 | unsigned long dn = PoissRandLimit(v,10.);
|
---|
971 | data_[ip] = (double)dn;
|
---|
972 | sum += (double)dn;
|
---|
973 | }
|
---|
974 | }
|
---|
975 | if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
|
---|
976 |
|
---|
977 | return sum;
|
---|
978 | }
|
---|
979 |
|
---|
980 | double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
|
---|
981 | // do NOT treate negative or nul mass -> let it as it is
|
---|
982 | // INPUT:
|
---|
983 | // massdist : distribution de masse (m*dn/dm)
|
---|
984 | // axeslog = false : retourne la masse
|
---|
985 | // = true : retourne le log10(mass)
|
---|
986 | // RETURN la masse totale
|
---|
987 | {
|
---|
988 | if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
|
---|
989 | check_array_alloc();
|
---|
990 |
|
---|
991 | double sum = 0.;
|
---|
992 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
993 | int_8 ip = IndexR(i,j,l);
|
---|
994 | double v = data_[ip];
|
---|
995 | if(v>0.) {
|
---|
996 | long ngal = long(v+0.1);
|
---|
997 | data_[ip] = 0.;
|
---|
998 | for(long i=0;i<ngal;i++) {
|
---|
999 | double m = massdist.RandomInterp(); // massdist.Random();
|
---|
1000 | if(axeslog) m = pow(10.,m);
|
---|
1001 | data_[ip] += m;
|
---|
1002 | }
|
---|
1003 | sum += data_[ip];
|
---|
1004 | }
|
---|
1005 | }
|
---|
1006 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
1007 |
|
---|
1008 | return sum;
|
---|
1009 | }
|
---|
1010 |
|
---|
1011 | void GeneFluct3D::AddAGN(double lfjy,double lsigma,double powlaw)
|
---|
1012 | // Add AGN flux into simulation:
|
---|
1013 | // --- Procedure:
|
---|
1014 | // 1. lancer "cmvdefsurv" avec les parametres du survey
|
---|
1015 | // (au redshift de reference du survey)
|
---|
1016 | // et recuperer l'angle solide "angsol sr" du pixel elementaire
|
---|
1017 | // au centre du cube.
|
---|
1018 | // 2. lancer "cmvtstagn" pour cet angle solide -> cmvtstagn.ppf
|
---|
1019 | // 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne
|
---|
1020 | // cad une moyenne <log10(S)> et un sigma "sig"
|
---|
1021 | // Attention: la distribution n'est pas gaussienne les "mean,sigma"
|
---|
1022 | // de l'histo ne sont pas vraiment ce que l'on veut
|
---|
1023 | // --- Limitations actuelle du code:
|
---|
1024 | // . les AGN sont supposes en loi de puissance IDENTIQUE pour tout theta,phi
|
---|
1025 | // . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee
|
---|
1026 | // . la distribution est approximee a une gaussienne
|
---|
1027 | // ... C'est une approximation pour un observateur loin du centre du cube
|
---|
1028 | // et pour un cube peu epais / distance observateur
|
---|
1029 | // --- Parametres de la routine:
|
---|
1030 | // llfy : c'est le <log10(S)> du flux depose dans un pixel par les AGN
|
---|
1031 | // lsigma : c'est le sigma de la distribution
|
---|
1032 | // powlaw : c'est la pente de ls distribution cad que le flux "lmsol"
|
---|
1033 | // et considere comme le flux a 1.4GHz et qu'on suppose une loi
|
---|
1034 | // F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz)
|
---|
1035 | // - Comme on est en echelle log10():
|
---|
1036 | // on tire log10(Msol) + X
|
---|
1037 | // ou X est une realisation sur une gaussienne de variance "sig^2"
|
---|
1038 | // La masse realisee est donc: Msol*10^X
|
---|
1039 | // - Pas de probleme de pixel negatif car on a une multiplication!
|
---|
1040 | {
|
---|
1041 | if(lp_>0) cout<<"--- AddAGN: <log10(S Jy)> = "<<lfjy<<" , sigma = "<<lsigma<<endl;
|
---|
1042 | check_array_alloc();
|
---|
1043 |
|
---|
1044 | if(cosmo_ == NULL || redshref_<0.| loscom2zred_.size()<1) {
|
---|
1045 | char *bla = "GeneFluct3D::AddAGN_Error: set Observator and Cosmology first";
|
---|
1046 | cout<<bla<<endl; throw ParmError(bla);
|
---|
1047 | }
|
---|
1048 |
|
---|
1049 | // La distance angulaire/luminosite/Dnu au centre
|
---|
1050 | double dangref = cosmo_->Dang(redshref_);
|
---|
1051 | double dlumref = cosmo_->Dlum(redshref_);
|
---|
1052 | double dredref = Dz_/(cosmo_->Dhubble()/cosmo_->E(redshref_));
|
---|
1053 | double dnuref = Fr_HyperFin_Par *dredref/pow(1.+redshref_,2.); // GHz
|
---|
1054 | double fagnref = pow(10.,lfjy)*(dnuref*1.e9); // Jy.Hz
|
---|
1055 | double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlumref); // Msol
|
---|
1056 | if(lp_>0) {
|
---|
1057 | cout<<"Au centre: z="<<redshref_<<", dredref="<<dredref<<", dnuref="<<dnuref<<" GHz"<<endl
|
---|
1058 | <<" dang="<<dangref<<" Mpc, dlum="<<dlumref<<" Mpc,"
|
---|
1059 | <<" fagnref="<<fagnref<<" Jy.Hz (a 1.4GHz), magnref="<<magnref<<" Msol"<<endl;
|
---|
1060 | }
|
---|
1061 |
|
---|
1062 | if(powlaw!=0.) {
|
---|
1063 | // F(nu) = (nu GHz/1.4 Ghz)^p * F(1.4GHz) et nu = 1.4 GHz / (1+z)
|
---|
1064 | magnref *= pow(1/(1.+redshref_),powlaw);
|
---|
1065 | if(lp_>0) cout<<" powlaw="<<powlaw<<" -> change magnref to "<<magnref<<" Msol"<<endl;
|
---|
1066 | }
|
---|
1067 |
|
---|
1068 | // Les infos en fonction de l'indice "l" selon Oz
|
---|
1069 | vector<double> correction;
|
---|
1070 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
1071 | for(long l=0;l<Nz_;l++) {
|
---|
1072 | double z = fabs(xobs_[2] - l*Dz_);
|
---|
1073 | double zred = interpinv(z);
|
---|
1074 | double dang = cosmo_->Dang(zred); // pour variation angle solide
|
---|
1075 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
1076 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
1077 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
1078 | double corr = dnu/dnuref*pow(dangref/dang*dlum/dlumref,2.);
|
---|
1079 | if(powlaw!=0.) corr *= pow((1.+redshref_)/(1.+zred),powlaw);
|
---|
1080 | correction.push_back(corr);
|
---|
1081 | if(lp_>0 && (l==0 || abs(l-Nz_/2)<2 || l==Nz_-1)) {
|
---|
1082 | cout<<"l="<<l<<" z="<<z<<" red="<<zred
|
---|
1083 | <<" da="<<dang<<" dlu="<<dlum<<" dred="<<dred
|
---|
1084 | <<" dnu="<<dnu<<" -> corr="<<corr<<endl;
|
---|
1085 | }
|
---|
1086 | }
|
---|
1087 |
|
---|
1088 | double sum=0., sum2=0., nsum=0.;
|
---|
1089 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) {
|
---|
1090 | double a = lsigma*NorRand();
|
---|
1091 | a = magnref*pow(10.,a);
|
---|
1092 | // On met le meme tirage le long de Oz (indice k)
|
---|
1093 | for(long l=0;l<Nz_;l++) {
|
---|
1094 | int_8 ip = IndexR(i,j,l);
|
---|
1095 | data_[ip] += a*correction[l];
|
---|
1096 | }
|
---|
1097 | sum += a; sum2 += a*a; nsum += 1.;
|
---|
1098 | }
|
---|
1099 |
|
---|
1100 | if(lp_>0 && nsum>1.) {
|
---|
1101 | sum /= nsum;
|
---|
1102 | sum2 = sum2/nsum - sum*sum;
|
---|
1103 | cout<<"...Mean mass="<<sum<<" Msol , s^2="<<sum2<<" s="<<sqrt(fabs(sum2))<<endl;
|
---|
1104 | }
|
---|
1105 |
|
---|
1106 | }
|
---|
1107 |
|
---|
1108 | void GeneFluct3D::AddNoise2Real(double snoise)
|
---|
1109 | // add noise to every pixels (meme les <=0 !)
|
---|
1110 | {
|
---|
1111 | if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<endl;
|
---|
1112 | check_array_alloc();
|
---|
1113 |
|
---|
1114 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
1115 | int_8 ip = IndexR(i,j,l);
|
---|
1116 | data_[ip] += snoise*NorRand();
|
---|
1117 | }
|
---|
1118 | }
|
---|
1119 |
|
---|
1120 |
|
---|
1121 |
|
---|
1122 | //-------------------------------------------------------------------
|
---|
1123 | //-------------------------------------------------------------------
|
---|
1124 | //--------------------- BRICOLO A NE PAS GARDER ---------------------
|
---|
1125 | //-------------------------------------------------------------------
|
---|
1126 | //-------------------------------------------------------------------
|
---|
1127 | int GeneFluct3D::ComputeSpectrum_bricolo(HistoErr& herr)
|
---|
1128 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
1129 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
1130 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
1131 | {
|
---|
1132 | if(lp_>0) cout<<"--- ComputeSpectrum_bricolo ---"<<endl;
|
---|
1133 | check_array_alloc();
|
---|
1134 |
|
---|
1135 | if(herr.NBins()<0) return -1;
|
---|
1136 | herr.Zero();
|
---|
1137 |
|
---|
1138 | // Attention a l'ordre
|
---|
1139 | for(long i=0;i<Nx_;i++) {
|
---|
1140 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
1141 | double kx = ii*Dkx_; kx *= kx;
|
---|
1142 | for(long j=0;j<Ny_;j++) {
|
---|
1143 | if(i==0 && j==0) continue;
|
---|
1144 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
1145 | double ky = jj*Dky_; ky *= ky;
|
---|
1146 | for(long l=1;l<NCz_;l++) {
|
---|
1147 | double kz = l*Dkz_; kz *= kz;
|
---|
1148 | double k = sqrt(kx+ky+kz);
|
---|
1149 | double pk = MODULE2(T_(l,j,i));
|
---|
1150 | herr.Add(k,pk);
|
---|
1151 | }
|
---|
1152 | }
|
---|
1153 | }
|
---|
1154 | herr.ToVariance();
|
---|
1155 |
|
---|
1156 | // renormalize to directly compare to original spectrum
|
---|
1157 | double norm = Vol_;
|
---|
1158 | herr *= norm;
|
---|
1159 |
|
---|
1160 | return 0;
|
---|
1161 | }
|
---|
1162 |
|
---|
1163 | int GeneFluct3D::ComputeSpectrum2D_bricolo(Histo2DErr& herr)
|
---|
1164 | {
|
---|
1165 | if(lp_>0) cout<<"--- ComputeSpectrum2D_bricolo ---"<<endl;
|
---|
1166 | check_array_alloc();
|
---|
1167 |
|
---|
1168 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
1169 | herr.Zero();
|
---|
1170 |
|
---|
1171 | // Attention a l'ordre
|
---|
1172 | for(long i=0;i<Nx_;i++) {
|
---|
1173 | long ii = (i>Nx_/2) ? Nx_-i : i;
|
---|
1174 | double kx = ii*Dkx_; kx *= kx;
|
---|
1175 | for(long j=0;j<Ny_;j++) {
|
---|
1176 | if(i==0 && j==0) continue;
|
---|
1177 | long jj = (j>Ny_/2) ? Ny_-j : j;
|
---|
1178 | double ky = jj*Dky_; ky *= ky;
|
---|
1179 | double kt = sqrt(kx+ky);
|
---|
1180 | for(long l=1;l<NCz_;l++) {
|
---|
1181 | double kz = l*Dkz_;
|
---|
1182 | double pk = MODULE2(T_(l,j,i));
|
---|
1183 | herr.Add(kt,kz,pk);
|
---|
1184 | }
|
---|
1185 | }
|
---|
1186 | }
|
---|
1187 | herr.ToVariance();
|
---|
1188 |
|
---|
1189 | // renormalize to directly compare to original spectrum
|
---|
1190 | double norm = Vol_;
|
---|
1191 | herr *= norm;
|
---|
1192 |
|
---|
1193 | return 0;
|
---|
1194 | }
|
---|