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

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

Intro raffinement dans simul AGN, variation dNu Da Dlum selon Oz cmv 05/04/2007

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