source: Sophya/trunk/Cosmo/SimLSS/genefluct3d.cc@ 3155

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

PrtTim deplaces de genefluct3d.cc -> cmvobserv3d.cc cmv 19/01/2007

File size: 24.9 KB
Line 
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 "integfunc.h"
23#include "geneutils.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//-------------------------------------------------------
32GeneFluct3D::GeneFluct3D(TArray< complex<r_8 > >& T)
33 : T_(T) , Nx_(0) , Ny_(0) , Nz_(0) , array_allocated_(false) , lp_(0)
34 , redshref_(-999.) , kredshref_(-999.)
35{
36 SetNThread();
37}
38
39GeneFluct3D::~GeneFluct3D(void)
40{
41 fftw_destroy_plan(pf_);
42 fftw_destroy_plan(pb_);
43#ifdef FFTW_THREAD
44 if(nthread_>0) fftw_cleanup_threads();
45#endif
46}
47
48//-------------------------------------------------------
49void GeneFluct3D::SetSize(long nx,long ny,long nz,double dx,double dy,double dz)
50{
51 setsize(nx,ny,nz,dx,dy,dz);
52 setalloc();
53 setpointers(false);
54 init_fftw();
55}
56
57void GeneFluct3D::SetObservator(double redshref,double kredshref)
58// L'observateur est au redshift z=0
59// est situe sur la "perpendiculaire" a la face x,y
60// issue du centre de cette face
61// Il faut positionner le cube sur l'axe des z cad des redshifts:
62// redshref = redshift de reference
63// Si redshref<0 alors redshref=0
64// kredshref = indice (en double) correspondant a ce redshift
65// Si kredshref<0 alors kredshref=0
66{
67 if(redshref<0.) redshref = 0.;
68 if(kredshref<0.) kredshref = 0.;
69 redshref_ = redshref;
70 kredshref_ = kredshref;
71}
72
73void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz)
74{
75 if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<<nx<<","<<ny<<","<<nz
76 <<" D="<<dx<<","<<dy<<","<<dz<<endl;
77 if(nx<=0 || dx<=0.) {
78 cout<<"GeneFluct3D::setsize_Error: bad value(s)"<<endl;
79 throw ParmError("GeneFluct3D::setsize_Error: bad value(s)");
80 }
81
82 // Les tailles des tableaux
83 Nx_ = nx;
84 Ny_ = ny; if(Ny_ <= 0) Ny_ = Nx_;
85 Nz_ = nz; if(Nz_ <= 0) Nz_ = Nx_;
86 N_.resize(0); N_.push_back(Nx_); N_.push_back(Ny_); N_.push_back(Nz_);
87 NRtot_ = Nx_*Ny_*Nz_; // nombre de pixels dans le survey
88 NCz_ = Nz_/2 +1;
89 NTz_ = 2*NCz_;
90
91 // le pas dans l'espace (Mpc)
92 Dx_ = dx;
93 Dy_ = dy; if(Dy_ <= 0.) Dy_ = Dx_;
94 Dz_ = dz; if(Dz_ <= 0.) Dz_ = Dx_;
95 D_.resize(0); D_.push_back(Dx_); D_.push_back(Dy_); D_.push_back(Dz_);
96 dVol_ = Dx_*Dy_*Dz_;
97 Vol_ = (Nx_*Dx_)*(Ny_*Dy_)*(Nz_*Dz_);
98
99 // Le pas dans l'espace de Fourier (Mpc^-1)
100 Dkx_ = 2.*M_PI/(Nx_*Dx_);
101 Dky_ = 2.*M_PI/(Ny_*Dy_);
102 Dkz_ = 2.*M_PI/(Nz_*Dz_);
103 Dk_.resize(0); Dk_.push_back(Dkx_); Dk_.push_back(Dky_); Dk_.push_back(Dkz_);
104 Dk3_ = Dkx_*Dky_*Dkz_;
105
106 // La frequence de Nyquist en k (Mpc^-1)
107 Knyqx_ = M_PI/Dx_;
108 Knyqy_ = M_PI/Dy_;
109 Knyqz_ = M_PI/Dz_;
110 Knyq_.resize(0); Knyq_.push_back(Knyqx_); Knyq_.push_back(Knyqy_); Knyq_.push_back(Knyqz_);
111}
112
113void GeneFluct3D::setalloc(void)
114{
115 if(lp_>1) cout<<"--- GeneFluct3D::setalloc ---"<<endl;
116 // Dimensionnement du tableau complex<r_8>
117 // ATTENTION: TArray adresse en memoire a l'envers du C
118 // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
119 sa_size_t SzK_[3] = {NCz_,Ny_,Nx_}; // a l'envers
120 try {
121 T_.ReSize(3,SzK_);
122 array_allocated_ = true;
123 } catch (...) {
124 cout<<"GeneFluct3D::setalloc_Error: Problem allocating T_"<<endl;
125 }
126 T_.SetMemoryMapping(BaseArray::CMemoryMapping);
127}
128
129void GeneFluct3D::setpointers(bool from_real)
130{
131 if(lp_>1) cout<<"--- GeneFluct3D::setpointers ---"<<endl;
132 if(from_real) T_ = ArrCastR2C(R_);
133 else R_ = ArrCastC2R(T_);
134 // On remplit les pointeurs
135 fdata_ = (fftw_complex *) (&T_(0,0,0));
136 data_ = (double *) (&R_(0,0,0));
137}
138
139void GeneFluct3D::check_array_alloc(void)
140// Pour tester si le tableau T_ est alloue
141{
142 if(array_allocated_) return;
143 char bla[90];
144 sprintf(bla,"GeneFluct3D::check_array_alloc_Error: array is not allocated");
145 cout<<bla<<endl;
146 throw ParmError(bla);
147}
148
149void GeneFluct3D::init_fftw(void)
150{
151 // --- Initialisation de fftw3 (attention data est sur-ecrit a l'init)
152 if(lp_>1) cout<<"--- GeneFluct3D::init_fftw ---"<<endl;
153#ifdef FFTW_THREAD
154 if(nthread_>0) {
155 cout<<"...Computing with "<<nthread_<<" threads"<<endl;
156 fftw_init_threads();
157 fftw_plan_with_nthreads(nthread_);
158 }
159#endif
160 if(lp_>1) cout<<"...forward plan"<<endl;
161 pf_ = fftw_plan_dft_r2c_3d(Nx_,Ny_,Nz_,data_,fdata_,FFTW_ESTIMATE);
162 if(lp_>1) cout<<"...backward plan"<<endl;
163 pb_ = fftw_plan_dft_c2r_3d(Nx_,Ny_,Nz_,fdata_,data_,FFTW_ESTIMATE);
164}
165
166
167//-------------------------------------------------------
168void GeneFluct3D::WriteFits(string cfname,int bitpix)
169{
170 cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
171 try {
172 FitsImg3DWriter fwrt(cfname.c_str(),bitpix,5);
173 fwrt.WriteKey("NX",Nx_," axe transverse 1");
174 fwrt.WriteKey("NY",Ny_," axe transverse 2");
175 fwrt.WriteKey("NZ",Nz_," axe longitudinal (redshift)");
176 fwrt.WriteKey("DX",Dx_," Mpc");
177 fwrt.WriteKey("DY",Dy_," Mpc");
178 fwrt.WriteKey("DZ",Dz_," Mpc");
179 fwrt.WriteKey("DKX",Dkx_," Mpc^-1");
180 fwrt.WriteKey("DKY",Dky_," Mpc^-1");
181 fwrt.WriteKey("DKZ",Dkz_," Mpc^-1");
182 fwrt.WriteKey("ZREF",redshref_," reference redshift");
183 fwrt.WriteKey("KZREF",kredshref_," reference redshift on z axe");
184 fwrt.Write(R_);
185 } catch (PThrowable & exc) {
186 cout<<"Exception : "<<(string)typeid(exc).name()
187 <<" - Msg= "<<exc.Msg()<<endl;
188 return;
189 } catch (...) {
190 cout<<" some other exception was caught !"<<endl;
191 return;
192 }
193}
194
195void GeneFluct3D::ReadFits(string cfname)
196{
197 cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
198 try {
199 FitsImg3DRead fimg(cfname.c_str(),0,5);
200 fimg.Read(R_);
201 long nx = fimg.ReadKeyL("NX");
202 long ny = fimg.ReadKeyL("NY");
203 long nz = fimg.ReadKeyL("NZ");
204 double dx = fimg.ReadKey("DX");
205 double dy = fimg.ReadKey("DY");
206 double dz = fimg.ReadKey("DZ");
207 double zref = fimg.ReadKey("ZREF");
208 double kzref = fimg.ReadKey("KZREF");
209 setsize(nx,ny,nz,dx,dy,dz);
210 setpointers(true);
211 init_fftw();
212 SetObservator(zref,kzref);
213 } catch (PThrowable & exc) {
214 cout<<"Exception : "<<(string)typeid(exc).name()
215 <<" - Msg= "<<exc.Msg()<<endl;
216 return;
217 } catch (...) {
218 cout<<" some other exception was caught !"<<endl;
219 return;
220 }
221}
222
223void GeneFluct3D::WritePPF(string cfname,bool write_real)
224// On ecrit soit le TArray<r_8> ou le TArray<complex <r_8> >
225{
226 cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
227 try {
228 R_.Info()["NX"] = (int_8)Nx_;
229 R_.Info()["NY"] = (int_8)Ny_;
230 R_.Info()["NZ"] = (int_8)Nz_;
231 R_.Info()["DX"] = (r_8)Dx_;
232 R_.Info()["DY"] = (r_8)Dy_;
233 R_.Info()["DZ"] = (r_8)Dz_;
234 R_.Info()["ZREF"] = (r_8)redshref_;
235 R_.Info()["KZREF"] = (r_8)kredshref_;
236 POutPersist pos(cfname.c_str());
237 if(write_real) pos << PPFNameTag("rgen") << R_;
238 else pos << PPFNameTag("pkgen") << T_;
239 } catch (PThrowable & exc) {
240 cout<<"Exception : "<<(string)typeid(exc).name()
241 <<" - Msg= "<<exc.Msg()<<endl;
242 return;
243 } catch (...) {
244 cout<<" some other exception was caught !"<<endl;
245 return;
246 }
247}
248
249void GeneFluct3D::ReadPPF(string cfname)
250{
251 cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
252 try {
253 bool from_real = true;
254 PInPersist pis(cfname.c_str());
255 string name_tag_k = "pkgen";
256 bool found_tag_k = pis.GotoNameTag("pkgen");
257 if(found_tag_k) {
258 cout<<" ...reading spectrun into TArray<complex <r_8> >"<<endl;
259 pis >> PPFNameTag("pkgen") >> T_;
260 from_real = false;
261 } else {
262 cout<<" ...reading space into TArray<r_8>"<<endl;
263 pis >> PPFNameTag("rgen") >> R_;
264 }
265 setpointers(from_real); // a mettre ici pour relire les DVInfo
266 int_8 nx = R_.Info()["NX"];
267 int_8 ny = R_.Info()["NY"];
268 int_8 nz = R_.Info()["NZ"];
269 r_8 dx = R_.Info()["DX"];
270 r_8 dy = R_.Info()["DY"];
271 r_8 dz = R_.Info()["DZ"];
272 r_8 zref = R_.Info()["ZREF"];
273 r_8 kzref = R_.Info()["KZREF"];
274 setsize(nx,ny,nz,dx,dy,dz);
275 init_fftw();
276 SetObservator(zref,kzref);
277 } catch (PThrowable & exc) {
278 cout<<"Exception : "<<(string)typeid(exc).name()
279 <<" - Msg= "<<exc.Msg()<<endl;
280 return;
281 } catch (...) {
282 cout<<" some other exception was caught !"<<endl;
283 return;
284 }
285}
286
287//-------------------------------------------------------
288void GeneFluct3D::Print(void)
289{
290 cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
291 cout<<"Space Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<Nz_<<" ("<<NTz_<<") size="
292 <<NRtot_<<endl;
293 cout<<" Resol: dx="<<Dx_<<" dy="<<Dy_<<" dz="<<Dz_<<" Mpc"
294 <<", dVol="<<dVol_<<", Vol="<<Vol_<<" Mpc^3"<<endl;
295 cout<<"Fourier Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<NCz_<<endl;
296 cout<<" Resol: dkx="<<Dkx_<<" dky="<<Dky_<<" dkz="<<Dkz_<<" Mpc^-1"
297 <<", Dk3="<<Dk3_<<" Mpc^-3"<<endl;
298 cout<<" (2Pi/k: "<<2.*M_PI/Dkx_<<" "<<2.*M_PI/Dky_<<" "<<2.*M_PI/Dkz_<<" Mpc)"<<endl;
299 cout<<" Nyquist: kx="<<Knyqx_<<" ky="<<Knyqy_<<" kz="<<Knyqz_<<" Mpc^-1"
300 <<", Kmax="<<GetKmax()<<" Mpc^-1"<<endl;
301 cout<<" (2Pi/k: "<<2.*M_PI/Knyqx_<<" "<<2.*M_PI/Knyqy_<<" "<<2.*M_PI/Knyqz_<<" Mpc)"<<endl;
302 cout<<"Redshift "<<redshref_<<" for z axe at k="<<kredshref_<<endl;
303}
304
305//-------------------------------------------------------
306void GeneFluct3D::ComputeFourier0(GenericFunc& pk_at_z)
307// cf ComputeFourier() mais avec autre methode de realisation du spectre
308// (attention on fait une fft pour realiser le spectre)
309{
310
311 // --- realisation d'un tableau de tirage gaussiens
312 if(lp_>0) cout<<"--- ComputeFourier0: before gaussian filling ---"<<endl;
313 // On tient compte du pb de normalisation de FFTW3
314 double sntot = sqrt((double)NRtot_);
315 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
316 int_8 ip = IndexR(i,j,l);
317 data_[ip] = NorRand()/sntot;
318 }
319
320 // --- realisation d'un tableau de tirage gaussiens
321 if(lp_>0) cout<<"...before fft real ---"<<endl;
322 fftw_execute(pf_);
323
324 // --- On remplit avec une realisation
325 if(lp_>0) cout<<"...before Fourier realization filling ---"<<endl;
326 T_(0,0,0) = complex<r_8>(0.); // on coupe le continue et on l'initialise
327 long lmod = Nx_/10; if(lmod<1) lmod=1;
328 for(long i=0;i<Nx_;i++) {
329 long ii = (i>Nx_/2) ? Nx_-i : i;
330 double kx = ii*Dkx_; kx *= kx;
331 if(lp_>0 && i%lmod==0) cout<<"i="<<i<<" ii="<<ii<<endl;
332 for(long j=0;j<Ny_;j++) {
333 long jj = (j>Ny_/2) ? Ny_-j : j;
334 double ky = jj*Dky_; ky *= ky;
335 for(long l=0;l<NCz_;l++) {
336 double kz = l*Dkz_; kz *= kz;
337 if(i==0 && j==0 && l==0) continue; // Suppression du continu
338 double k = sqrt(kx+ky+kz);
339 // cf normalisation: Peacock, Cosmology, formule 16.38 p504
340 double pk = pk_at_z(k)/Vol_;
341 // ici pas de "/2" a cause de la remarque ci-dessus
342 T_(l,j,i) *= sqrt(pk);
343 }
344 }
345 }
346
347 if(lp_>0) cout<<"...computing power"<<endl;
348 double p = compute_power_carte();
349 if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
350
351}
352
353//-------------------------------------------------------
354void GeneFluct3D::ComputeFourier(GenericFunc& pk_at_z)
355// Calcule une realisation du spectre "pk_at_z"
356// Attention: dans TArray le premier indice varie le + vite
357// Explication normalisation: see Coles & Lucchin, Cosmology, p264-265
358// FFTW3: on note N=Nx*Ny*Nz
359// f --(FFT)--> F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f))
360// sum(f(x_i)^2) = S
361// sum(F(nu_i)^2) = S*N
362// sum(fb(x_i)^2) = S*N^2
363{
364 // --- RaZ du tableau
365 T_ = complex<r_8>(0.);
366
367 // --- On remplit avec une realisation
368 if(lp_>0) cout<<"--- ComputeFourier ---"<<endl;
369 long lmod = Nx_/10; if(lmod<1) lmod=1;
370 for(long i=0;i<Nx_;i++) {
371 long ii = (i>Nx_/2) ? Nx_-i : i;
372 double kx = ii*Dkx_; kx *= kx;
373 if(lp_>0 && i%lmod==0) cout<<"i="<<i<<" ii="<<ii<<endl;
374 for(long j=0;j<Ny_;j++) {
375 long jj = (j>Ny_/2) ? Ny_-j : j;
376 double ky = jj*Dky_; ky *= ky;
377 for(long l=0;l<NCz_;l++) {
378 double kz = l*Dkz_; kz *= kz;
379 if(i==0 && j==0 && l==0) continue; // Suppression du continu
380 double k = sqrt(kx+ky+kz);
381 // cf normalisation: Peacock, Cosmology, formule 16.38 p504
382 double pk = pk_at_z(k)/Vol_;
383 // Explication de la division par 2: voir perandom.cc
384 // ou egalement Coles & Lucchin, Cosmology formula 13.7.2 p279
385 T_(l,j,i) = ComplexGaussRan(sqrt(pk/2.));
386 }
387 }
388 }
389
390 manage_coefficients(); // gros effet pour les spectres que l'on utilise !
391
392 if(lp_>0) cout<<"...computing power"<<endl;
393 double p = compute_power_carte();
394 if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
395
396}
397
398long GeneFluct3D::manage_coefficients(void)
399// Take into account the real and complexe conjugate coefficients
400// because we want a realization of a real data in real space
401{
402 if(lp_>1) cout<<"...managing coefficients"<<endl;
403 check_array_alloc();
404
405 // 1./ Le Continu et Nyquist sont reels
406 long nreal = 0;
407 for(long kk=0;kk<2;kk++) {
408 long k=0; // continu
409 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
410 for(long jj=0;jj<2;jj++) {
411 long j=0;
412 if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
413 for(long ii=0;ii<2;ii++) {
414 long i=0;
415 if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
416 int_8 ip = IndexC(i,j,k);
417 //cout<<"i="<<i<<" j="<<j<<" k="<<k<<" = ("<<fdata_[ip][0]<<","<<fdata_[ip][1]<<")"<<endl;
418 fdata_[ip][1] = 0.; fdata_[ip][0] *= M_SQRT2;
419 nreal++;
420 }
421 }
422 }
423 if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
424
425 // 2./ Les elements complexe conjugues (tous dans le plan k=0,Nyquist)
426
427 // a./ les lignes et colonnes du continu et de nyquist
428 long nconj1 = 0;
429 for(long kk=0;kk<2;kk++) {
430 long k=0; // continu
431 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
432 for(long jj=0;jj<2;jj++) { // selon j
433 long j=0;
434 if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
435 for(long i=1;i<(Nx_+1)/2;i++) {
436 int_8 ip = IndexC(i,j,k);
437 int_8 ip1 = IndexC(Nx_-i,j,k);
438 fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
439 nconj1++;
440 }
441 }
442 for(long ii=0;ii<2;ii++) {
443 long i=0;
444 if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
445 for(long j=1;j<(Ny_+1)/2;j++) {
446 int_8 ip = IndexC(i,j,k);
447 int_8 ip1 = IndexC(i,Ny_-j,k);
448 fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
449 nconj1++;
450 }
451 }
452 }
453 if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
454
455 // b./ les lignes et colonnes hors continu et de nyquist
456 long nconj2 = 0;
457 for(long kk=0;kk<2;kk++) {
458 long k=0; // continu
459 if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
460 for(long j=1;j<(Ny_+1)/2;j++) {
461 if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
462 for(long i=1;i<Nx_;i++) {
463 if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
464 int_8 ip = IndexC(i,j,k);
465 int_8 ip1 = IndexC(Nx_-i,Ny_-j,k);
466 fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
467 nconj2++;
468 }
469 }
470 }
471 if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
472
473 if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(Nx_*Ny_*NCz_-nconj1-nconj2)-8<<endl;
474
475 return nreal+nconj1+nconj2;
476}
477
478double GeneFluct3D::compute_power_carte(void)
479// Calcul la puissance de la realisation du spectre Pk
480{
481 check_array_alloc();
482
483 double s2 = 0.;
484 for(long l=0;l<NCz_;l++)
485 for(long j=0;j<Ny_;j++)
486 for(long i=0;i<Nx_;i++) s2 += MODULE2(T_(l,j,i));
487
488 double s20 = 0.;
489 for(long j=0;j<Ny_;j++)
490 for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
491
492 double s2n = 0.;
493 if(Nz_%2==0)
494 for(long j=0;j<Ny_;j++)
495 for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
496
497 return 2.*s2 -s20 -s2n;
498}
499
500//-------------------------------------------------------------------
501void GeneFluct3D::FilterByPixel(void)
502// Filtrage par la fonction fenetre du pixel (parallelepipede)
503// TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
504// e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
505// = 2/(k_x*dx) * sin(k_x*dx/2) * (idem y) * (idem z)
506// Gestion divergence en 0: sin(y)/y = 1 - y^2/6*(1-y^2/20)
507// avec y = k_x*dx/2
508{
509 if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
510 check_array_alloc();
511
512 for(long i=0;i<Nx_;i++) {
513 long ii = (i>Nx_/2) ? Nx_-i : i;
514 double kx = ii*Dkx_ *Dx_/2;
515 double pk_x = pixelfilter(kx);
516 for(long j=0;j<Ny_;j++) {
517 long jj = (j>Ny_/2) ? Ny_-j : j;
518 double ky = jj*Dky_ *Dy_/2;
519 double pk_y = pixelfilter(ky);
520 for(long l=0;l<NCz_;l++) {
521 double kz = l*Dkz_ *Dz_/2;
522 double pk_z = pixelfilter(kz);
523 T_(l,j,i) *= pk_x*pk_y*pk_z;
524 }
525 }
526 }
527
528}
529
530//-------------------------------------------------------------------
531void GeneFluct3D::ComputeReal(void)
532// Calcule une realisation dans l'espace reel
533{
534 if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
535 check_array_alloc();
536
537 // On fait la FFT
538 fftw_execute(pb_);
539}
540
541//-------------------------------------------------------------------
542void GeneFluct3D::ReComputeFourier(void)
543{
544 if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
545 check_array_alloc();
546
547 // On fait la FFT
548 fftw_execute(pf_);
549 // On corrige du pb de la normalisation de FFTW3
550 double v = (double)NRtot_;
551 for(long i=0;i<Nx_;i++)
552 for(long j=0;j<Ny_;j++)
553 for(long l=0;l<NCz_;l++) T_(l,j,i) /= complex<r_8>(v);
554
555}
556
557//-------------------------------------------------------------------
558int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
559// Compute spectrum from "T" and fill HistoErr "herr"
560// T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
561// cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
562{
563 if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
564 check_array_alloc();
565
566 if(herr.NBins()<0) return -1;
567 herr.Zero();
568
569 // Attention a l'ordre
570 for(long i=0;i<Nx_;i++) {
571 long ii = (i>Nx_/2) ? Nx_-i : i;
572 double kx = ii*Dkx_; kx *= kx;
573 for(long j=0;j<Ny_;j++) {
574 long jj = (j>Ny_/2) ? Ny_-j : j;
575 double ky = jj*Dky_; ky *= ky;
576 for(long l=0;l<NCz_;l++) {
577 double kz = l*Dkz_; kz *= kz;
578 double k = sqrt(kx+ky+kz);
579 double pk = MODULE2(T_(l,j,i));
580 herr.Add(k,pk);
581 }
582 }
583 }
584 herr.ToVariance();
585
586 // renormalize to directly compare to original spectrum
587 double norm = Vol_;
588 herr *= norm;
589
590 return 0;
591}
592
593int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
594{
595 if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
596 check_array_alloc();
597
598 if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
599 herr.Zero();
600
601 // Attention a l'ordre
602 for(long i=0;i<Nx_;i++) {
603 long ii = (i>Nx_/2) ? Nx_-i : i;
604 double kx = ii*Dkx_; kx *= kx;
605 for(long j=0;j<Ny_;j++) {
606 long jj = (j>Ny_/2) ? Ny_-j : j;
607 double ky = jj*Dky_; ky *= ky;
608 double kt = sqrt(kx+ky);
609 for(long l=0;l<NCz_;l++) {
610 double kz = l*Dkz_;
611 double pk = MODULE2(T_(l,j,i));
612 herr.Add(kt,kz,pk);
613 }
614 }
615 }
616 herr.ToVariance();
617
618 // renormalize to directly compare to original spectrum
619 double norm = Vol_;
620 herr *= norm;
621
622 return 0;
623}
624
625//-------------------------------------------------------
626int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
627// Recompute MASS variance in spherical top-hat (rayon=R)
628{
629 if(lp_>0) cout<<"--- VarianceFrReal ---"<<endl;
630 check_array_alloc();
631
632 long dnx = long(R/Dx_+0.5); if(dnx<=0) dnx = 1;
633 long dny = long(R/Dy_+0.5); if(dny<=0) dny = 1;
634 long dnz = long(R/Dz_+0.5); if(dnz<=0) dnz = 1;
635 if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
636
637 double sum=0., sum2=0., r2 = R*R; int_8 nsum=0;
638
639 for(long i=dnx;i<Nx_-dnx;i+=dnx) {
640 for(long j=dny;j<Ny_-dny;j+=dny) {
641 for(long l=dnz;l<Nz_-dnz;l+=dnz) {
642 double s=0.; int_8 n=0;
643 for(long ii=i-dnx;ii<=i+dnx;ii++) {
644 double x = (ii-i)*Dx_; x *= x;
645 for(long jj=j-dny;jj<=j+dny;jj++) {
646 double y = (jj-j)*Dy_; y *= y;
647 for(long ll=l-dnz;ll<=l+dnz;ll++) {
648 double z = (ll-l)*Dz_; z *= z;
649 if(x+y+z>r2) continue;
650 int_8 ip = IndexR(ii,jj,ll);
651 s += 1.+data_[ip];
652 n++;
653 }
654 }
655 }
656 if(n>0) {sum += s; sum2 += s*s; nsum++;}
657 //cout<<i<<","<<j<<","<<l<<" n="<<n<<" s="<<s<<" sum="<<sum<<" sum2="<<sum2<<endl;
658 }
659 }
660 }
661
662 if(nsum<=1) {var=0.; return nsum;}
663
664 sum /= nsum;
665 sum2 = sum2/nsum - sum*sum;
666 if(lp_>0) cout<<"VarianceFrReal: nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
667
668 var = sum2/(sum*sum); // <dM>^2/<M>^2
669 if(lp_>0) cout<<"sigmaR^2="<<var<<" -> "<<sqrt(var)<<endl;
670
671 return nsum;
672}
673
674//-------------------------------------------------------
675int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
676// number of pixels outside of ]vmin,vmax[ extremites exclues
677// -> vmin and vmax are considered as bad
678{
679 check_array_alloc();
680
681 int_8 nbad = 0;
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 double v = data_[ip];
685 if(v<=vmin || v>=vmax) nbad++;
686 }
687
688 return nbad;
689}
690
691int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax)
692// mean,sigma^2 pour pixels avec valeurs ]vmin,vmax[ extremites exclues
693// -> mean and sigma^2 are NOT computed with pixels values vmin and vmax
694{
695 check_array_alloc();
696
697 int_8 n = 0;
698 rm = rs2 = 0.;
699
700 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
701 int_8 ip = IndexR(i,j,l);
702 double v = data_[ip];
703 if(v<=vmin || v>=vmax) continue;
704 rm += v;
705 rs2 += v*v;
706 n++;
707 }
708
709 if(n>1) {
710 rm /= (double)n;
711 rs2 = rs2/(double)n - rm*rm;
712 }
713
714 return n;
715}
716
717int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
718// set to "val0" if out of range ]vmin,vmax[ extremites exclues
719// -> vmin and vmax are set to val0
720{
721 check_array_alloc();
722
723 int_8 nbad = 0;
724 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
725 int_8 ip = IndexR(i,j,l);
726 double v = data_[ip];
727 if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
728 }
729
730 return nbad;
731}
732
733//-------------------------------------------------------
734void GeneFluct3D::TurnFluct2Mass(void)
735// d_rho/rho -> Mass (add one!)
736{
737 if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
738 check_array_alloc();
739
740
741 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
742 int_8 ip = IndexR(i,j,l);
743 data_[ip] += 1.;
744 }
745}
746
747double GeneFluct3D::TurnMass2MeanNumber(double n_by_mpc3)
748// do NOT treate negative or nul values
749{
750 if(lp_>0) cout<<"--- TurnMass2MeanNumber ---"<<endl;
751
752 double m,s2;
753 int_8 ngood = MeanSigma2(m,s2,0.,1e+200);
754 if(lp_>0) cout<<"...ngood="<<ngood
755 <<" m="<<m<<" s2="<<s2<<" -> "<<sqrt(s2)<<endl;
756
757 // On doit mettre n*Vol galaxies dans notre survey
758 // On en met uniquement dans les pixels de masse >0.
759 // On NE met PAS a zero les pixels <0
760 // On renormalise sur les pixels>0 pour qu'on ait n*Vol galaxies
761 // comme on ne prend que les pixels >0, on doit normaliser
762 // a la moyenne de <1+d_rho/rho> sur ces pixels
763 // (rappel sur tout les pixels <1+d_rho/rho>=1)
764 double dn = n_by_mpc3*Vol_/m /(double)ngood; // nb de gal a mettre ds 1 px
765 if(lp_>0) cout<<"...galaxy density move from "
766 <<n_by_mpc3*Vol_/double(NRtot_)<<" to "<<dn<<" / pixel"<<endl;
767 double sum = 0.;
768 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
769 int_8 ip = IndexR(i,j,l);
770 data_[ip] *= dn; // par coherence on multiplie aussi les <=0
771 if(data_[ip]>0.) sum += data_[ip]; // mais on ne les compte pas
772 }
773 if(lp_>0) cout<<sum<<"...galaxies put into survey / "<<n_by_mpc3*Vol_<<endl;
774
775 return sum;
776}
777
778double GeneFluct3D::ApplyPoisson(void)
779// do NOT treate negative or nul mass -> let it as it is
780{
781 if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
782 check_array_alloc();
783
784 double sum = 0.;
785 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
786 int_8 ip = IndexR(i,j,l);
787 double v = data_[ip];
788 if(v>0.) {
789 unsigned long dn = PoissRandLimit(v,10.);
790 data_[ip] = (double)dn;
791 sum += (double)dn;
792 }
793 }
794 if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
795
796 return sum;
797}
798
799double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
800// do NOT treate negative or nul mass -> let it as it is
801// INPUT:
802// massdist : distribution de masse (m*dn/dm)
803// axeslog = false : retourne la masse
804// = true : retourne le log10(mass)
805// RETURN la masse totale
806{
807 if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
808 check_array_alloc();
809
810 double sum = 0.;
811 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
812 int_8 ip = IndexR(i,j,l);
813 double v = data_[ip];
814 if(v>0.) {
815 long ngal = long(v+0.1);
816 data_[ip] = 0.;
817 for(long i=0;i<ngal;i++) {
818 double m = massdist.RandomInterp(); // massdist.Random();
819 if(axeslog) m = pow(10.,m);
820 data_[ip] += m;
821 }
822 sum += data_[ip];
823 }
824 }
825 if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
826
827 return sum;
828}
829
830void GeneFluct3D::AddNoise2Real(double snoise)
831// add noise to every pixels (meme les <=0 !)
832{
833 if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<endl;
834 check_array_alloc();
835
836 for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
837 int_8 ip = IndexR(i,j,l);
838 data_[ip] += snoise*NorRand();
839 }
840}
Note: See TracBrowser for help on using the repository browser.