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

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

les AGN selon C.Jackson, une premiere approche simplifiee, recodage from Jim Rich. cmv 03/04/2007

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