#include "machdefs.h" #include #include #include #include #include #include #include "tarray.h" #include "pexceptions.h" #include "perandom.h" #include "srandgen.h" #include "fabtcolread.h" #include "fabtwriter.h" #include "fioarr.h" #include "arrctcast.h" #include "constcosmo.h" #include "geneutils.h" #include "schechter.h" #include "genefluct3d.h" #define WITH_FFTW_THREAD #define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag())) namespace SOPHYA { //------------------------------------------------------- GeneFluct3D::GeneFluct3D(long nx,long ny,long nz,double dx,double dy,double dz ,unsigned short nthread,int lp) { init_default(); lp_ = lp; nthread_ = nthread; setsize(nx,ny,nz,dx,dy,dz); setalloc(); setpointers(false); init_fftw(); } GeneFluct3D::GeneFluct3D(unsigned short nthread) { init_default(); setsize(2,2,2,1.,1.,1.); nthread_ = nthread; setalloc(); setpointers(false); init_fftw(); } GeneFluct3D::~GeneFluct3D(void) { delete_fftw(); } //------------------------------------------------------- void GeneFluct3D::init_default(void) { Nx_ = Ny_ = Nz_ = 0; is_set_fftw_plan = false; nthread_ = 0; lp_ = 0; array_allocated_ = false; cosmo_ = NULL; growth_ = NULL; redsh_ref_ = -999.; kredsh_ref_ = 0.; dred_ref_ = -999.; loscom_ref_ = -999.; dtrc_ref_ = dlum_ref_ = dang_ref_ = -999.; nu_ref_ = dnu_ref_ = -999.; loscom_min_ = loscom_max_ = -999.; loscom2zred_min_ = loscom2zred_max_ = 0.; xobs_[0] = xobs_[1] = xobs_[2] = 0.; zred_.resize(0); loscom_.resize(0); loscom2zred_.resize(0); } void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz) { if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<1) cout<<"--- GeneFluct3D::setalloc ---"< // ATTENTION: TArray adresse en memoire a l'envers du C // Tarray(n1,n2,n3) == Carray[n3][n2][n1] sa_size_t SzK_[3] = {NCz_,Ny_,Nx_}; // a l'envers try { T_.ReSize(3,SzK_); array_allocated_ = true; if(lp_>1) cout<<" allocating: "<)/1.e6<<" Mo"<1) cout<<"--- GeneFluct3D::setpointers ---"<1) cout<<"--- GeneFluct3D::init_fftw ---"<0) { cout<<"...Computing with "<1) cout<<"...forward plan"<1) cout<<"...backward plan"<0) fftw_cleanup_threads(); #endif is_set_fftw_plan = false; } void GeneFluct3D::check_array_alloc(void) // Pour tester si le tableau T_ est alloue { if(array_allocated_) return; char bla[90]; sprintf(bla,"GeneFluct3D::check_array_alloc_Error: array is not allocated"); cout< Le pixel i=nx/2 j=ny/2 k=250.75 est au redshift 1.5 { if(redshref<0.) redshref = 0.; if(kredshref<0.) { if(Nz_<=0) { char *bla = "GeneFluct3D::SetObservator_Error: for kredsh_ref<0 define cube geometry first"; cout<0) cout<<"--- GeneFluct3D::SetObservator zref="< zred // { if(lp_>0) cout<<"--- LosComRedshift: zinc="<Dhubble()/cosmo_->E(redsh_ref_)); loscom_ref_ = cosmo_->Dloscom(redsh_ref_); dtrc_ref_ = cosmo_->Dtrcom(redsh_ref_); dlum_ref_ = cosmo_->Dlum(redsh_ref_); dang_ref_ = cosmo_->Dang(redsh_ref_); nu_ref_ = Fr_HyperFin_Par/(1.+redsh_ref_); // GHz dnu_ref_ = Fr_HyperFin_Par *dred_ref_/pow(1.+redsh_ref_,2.); // GHz if(lp_>0) { cout<<"...reference pixel redshref="<=0. && xobs_[2]<=Nz_*Dz_) obs_in_cube = true; // Find MINIMUM los com distance to the observer: // c'est le centre de la face a k=0 // (ou zero si l'observateur est dans le cube) loscom_min_ = 0.; if(!obs_in_cube) loscom_min_ = -xobs_[2]; // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED if(loscom_min_<=1.e-50) for(int i=0;i<50;i++) cout<<"ATTENTION TOUTES LES PARTIES DU CODE NE MARCHENT PAS POUR UN OBSERVATEUR DANS LE CUBE"<loscom_max_) loscom_max_ = dz2; } } } if(lp_>0) { cout<<"...zref="<loscom_max if(zinc<=0.) zinc = 0.01; for(double z=0.; ; z+=zinc) { double dlc = cosmo_->Dloscom(z); if(dlcloscom_max_) break; // on sort apres avoir stoque un dlc>dlcmax } if(lp_>0) { long n = zred_.size(); cout<<"...zred/loscom tables[zinc="< d="<1) cout<<" , z="< d="<zred if(npoints<3) npoints = zred_.size(); InverseFunc invfun(zred_,loscom_); invfun.ComputeParab(npoints,loscom2zred_); loscom2zred_min_ = invfun.YMin(); loscom2zred_max_ = invfun.YMax(); if(lp_>0) { long n = loscom2zred_.size(); cout<<"...loscom -> zred[npoints="< zred=["; if(n>0) cout<1) cout<1 && n>0) for(int i=0;i=n-2) cout<<" i="< ou le TArray > { cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="< >"<> PPFNameTag("pkgen") >> T_; from_real = false; } else { cout<<" ...reading space into TArray"<> PPFNameTag("rgen") >> R_; } setpointers(from_real); // a mettre ici pour relire les DVInfo int_8 nx = R_.Info()["NX"]; int_8 ny = R_.Info()["NY"]; int_8 nz = R_.Info()["NZ"]; r_8 dx = R_.Info()["DX"]; r_8 dy = R_.Info()["DY"]; r_8 dz = R_.Info()["DZ"]; r_8 zref = R_.Info()["ZREF"]; r_8 kzref = R_.Info()["KZREF"]; setsize(nx,ny,nz,dx,dy,dz); init_fftw(); SetObservator(zref,kzref); array_allocated_ = true; } catch (PThrowable & exc) { cout<<"Exception : "<<(string)typeid(exc).name() <<" - Msg= "< S; char str[16]; long i,j,l; // Tranches en Z for(int s=0;s<3;s++) { S.ReSize(Nx_,Ny_); if(s==0) l=0; else if(s==1) l=(Nz_+1)/2; else l=Nz_-1; sprintf(str,"z%ld",l); for(i=0;i0) cout<<"...before fft real ---"<0) cout<<"...before Fourier realization filling"<(0.); // on coupe le continue et on l'initialise long lmod = Nx_/10; if(lmod<1) lmod=1; for(long i=0;iNx_/2) ? Nx_-i : i; double kx = ii*Dkx_; kx *= kx; if(lp_>0 && i%lmod==0) cout<<"i="<0) cout<<"...computing power"<0) cout<<"Puissance dans la realisation: "< F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f)) // sum(f(x_i)^2) = S // sum(F(nu_i)^2) = S*N // sum(fb(x_i)^2) = S*N^2 { // --- RaZ du tableau T_ = complex(0.); // --- On remplit avec une realisation if(lp_>0) cout<<"--- ComputeFourier ---"<Nx_/2) ? Nx_-i : i; double kx = ii*Dkx_; kx *= kx; if(lp_>0 && i%lmod==0) cout<<"i="<0) cout<<"Puissance dans la realisation: "<1) cout<<"...managing coefficients"<1) cout<<"Number of forced real number ="<1) cout<<"Number of forced conjugate on cont+nyq ="<1) cout<<"Number of forced conjugate hors cont+nyq ="<1) cout<<"Check: ddl= "<Nx_/2) ? Nx_-i : i; double kx = ii*Dkx_ *Dx_/2; double pk_x = pixelfilter(kx); for(long j=0;jNy_/2) ? Ny_-j : j; double ky = jj*Dky_ *Dy_/2; double pk_y = pixelfilter(ky); for(long l=0;l0) cout<<"--- ApplyGrowthFactor: evol="<2) { char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value"; cout<Linear(z,ok); // interpolation lineaire //double dzgr = growth_->Parab(z,ok); // interpolation parabolique int_8 ip = IndexR(i,j,l); data_[ip] *= dzgr; } } } //CHECK: {POutPersist pos("applygrowth.ppf"); string tag="hgr"; pos.PutObject(hgr,tag);} } //------------------------------------------------------------------- void GeneFluct3D::ComputeReal(void) // Calcule une realisation dans l'espace reel { if(lp_>0) cout<<"--- ComputeReal ---"<0) cout<<"--- ReComputeFourier ---"<(v); } //------------------------------------------------------------------- int GeneFluct3D::ComputeSpectrum(HistoErr& herr) // Compute spectrum from "T" and fill HistoErr "herr" // T : dans le format standard de GeneFuct3D: T(nz,ny,nx) // cad T(kz,ky,kx) avec 00) cout<<"--- ComputeSpectrum ---"<Nx_/2) ? Nx_-i : i; double kx = ii*Dkx_; kx *= kx; for(long j=0;jNy_/2) ? Ny_-j : j; double ky = jj*Dky_; ky *= ky; for(long l=0;l0) cout<<"--- ComputeSpectrum2D ---"<Nx_/2) ? Nx_-i : i; double kx = ii*Dkx_; kx *= kx; for(long j=0;jNy_/2) ? Ny_-j : j; double ky = jj*Dky_; ky *= ky; double kt = sqrt(kx+ky); for(long l=0;l0) cout<<"--- ComputeSpectrum: sigma="< vfz(NCz_); if(pixcor) // kz = l*Dkz_ for(long l=0;lNx_/2) ? Nx_-i : i; double kx = ii*Dkx_; double fx = (pixcor) ? pixelfilter(kx*Dx_/2): 1.; kx *= kx; fx *= fx; for(long j=0;jNy_/2) ? Ny_-j : j; double ky = jj*Dky_; double fy = (pixcor) ? pixelfilter(ky*Dy_/2): 1.; ky *= ky; fy *= fy; for(long l=0;l0.) herr.Add(k,pk/f); } } } herr.ToVariance(); for(int i=0;i0) cout<<"--- ComputeSpectrum2D: sigma="< vfz(NCz_); if(pixcor) // kz = l*Dkz_ for(long l=0;lNx_/2) ? Nx_-i : i; double kx = ii*Dkx_; double fx = (pixcor) ? pixelfilter(kx*Dx_/2) : 1.; kx *= kx; fx *= fx; for(long j=0;jNy_/2) ? Ny_-j : j; double ky = jj*Dky_; double fy = (pixcor) ? pixelfilter(ky*Dy_/2) : 1.; ky *= ky; fy *= fy; double kt = sqrt(kx+ky); for(long l=0;l0.) herr.Add(kt,kz,pk/f); } } } herr.ToVariance(); for(int i=0;i0) cout<<"--- VarianceFrReal R="<0) cout<<"dnx="<="<)^2>="<^2/^2 if(lp_>0) cout<<"...sigmaR^2="< "< vmin and vmax are considered as bad { check_array_alloc(); int_8 nbad = 0; for(long i=0;i=vmax) nbad++; } if(lp_>0) cout<<"--- NumberOfBad "<vmin)? true: false; if(lp_>0) { cout<<"--- MinMax"; if(tstval) cout<<" range=]"<=vmax)) continue; if(xxmax) xmax = x; n++; } if(lp_>0) cout<<" n="<vmin)? true: false; if(lp_>0) { cout<<"--- MeanSigma2"; if(tstval) cout<<" range=]"< vmin and vmax are set to val0 { check_array_alloc(); int_8 nbad = 0; for(long i=0;i=vmax) {data_[ip] = val0; nbad++;} } if(lp_>0) cout<<"--- SetToVal "<0) cout<<"--- ScaleCube scale="<0) cout<<"--- TurnMass2MeanNumber : "<0.) {mgood += data_[ip]; ngood++;} } if(ngood>0) mgood /= (double)ngood; if(nall>0) mall /= (double)nall; if(lp_>0) cout<<"...ngood="<0. // On MET a zero les pixels <0 // On renormalise sur les pixels>0 pour qu'on ait n*Vol galaxies (ou Msol) // comme on ne prend que les pixels >0, on doit normaliser // a la moyenne de <1+d_rho/rho> sur ces pixels // (rappel sur tout les pixels <1+d_rho/rho>=1) // nb de gal (ou Msol) a mettre ds 1 px: double dn = val_by_mpc3*Vol_/ (mgood/mall) /(double)ngood; if(lp_>0) cout<<"...quantity density move from " <0) cout< let it as it is { if(lp_>0) cout<<"--- ApplyPoisson ---"<0.) { unsigned long dn = PoissRandLimit(v,10.); data_[ip] = (double)dn; sum += (double)dn; } } if(lp_>0) cout< let it as it is // INPUT: // massdist : distribution de masse (m*dn/dm) // axeslog = false : retourne la masse // = true : retourne le log10(mass) // RETURN la masse totale { if(lp_>0) cout<<"--- TurnNGal2Mass ---"<0.) { long ngal = long(v+0.1); data_[ip] = 0.; for(long i=0;i0) cout<0) cout<<"--- TurnNGal2MassQuick ---"<0.) { long ngal = long(v+0.1); data_[ip] = schmdist.TirMass(ngal); sum += data_[ip]; } } if(lp_>0) cout<0) cout<<"--- AddNoise2Real: snoise = "< correction; InterpFunc *intercor = NULL; if(type_evol>0) { // Sigma_Noise(en mass) : // Slim ~ 1/sqrt(DNu) * sqrt(nlobe) en W/m^2Hz // Flim ~ sqrt(DNu) * sqrt(nlobe) en W/m^2 // Mlim ~ sqrt(DNu) * (Dlum)^2 * sqrt(nlobe) en Msol // nlobe ~ 1/Dtrcom^2 // Mlim ~ sqrt(DNu) * (Dlum)^2 / Dtrcom if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) { char *bla = "GeneFluct3D::AddNoise2Real_Error: set Observator and Cosmology first"; cout<10)? nsz/10: 1); for(long i=0;iDtrcom(zred); // pour variation angle solide double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred)); double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu double corr = sqrt(dnu/dnu_ref_) * pow(dlum/dlum_ref_,2.) * dtrc_ref_/dtrc; if(lp_>0 && (i==0 || i==nsz-1 || i%nszmod==0)) cout<<"i="< cor="<0) { double dz = DZcom(l); if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz); else dz = fabs(dz); // tous les plans Z au meme redshift corr = (*intercor)(dz); if(corrcorrlim[1]) corrlim[1]=corr; } int_8 ip = IndexR(i,j,l); data_[ip] += snoise*corr*NorRand(); } } } if(type_evol>0) cout<<"correction factor range: ["< cmvtstagn.ppf // 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne // cad une moyenne et un sigma "sig" // Attention: la distribution n'est pas gaussienne les "mean,sigma" // de l'histo ne sont pas vraiment ce que l'on veut // --- Limitations actuelle du code: // . les AGN sont supposes evoluer avec la meme loi de puissance pour tout theta,phi // . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee // . la distribution est approximee a une gaussienne // ... C'est une approximation pour un observateur loin du centre du cube // et pour un cube peu epais / distance observateur // --- Parametres de la routine: // llfy : c'est le du flux depose par les AGN // dans l'angle solide du pixel elementaire de reference du cube // lsigma : c'est le sigma de la distribution des log10(S) // powlaw : c'est la pente de la distribution cad que le flux "lmsol" // et considere comme le flux a 1.4GHz et qu'on suppose une loi // F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz) // - Comme on est en echelle log10(): // on tire log10(Msol) + X // ou X est une realisation sur une gaussienne de variance "sig^2" // La masse realisee est donc: Msol*10^X // - Pas de probleme de pixel negatif car on a une multiplication! { if(lp_>0) cout<<"--- AddAGN: = "<0) cout<<"Au pixel de ref: fagnref="<0) cout<<" powlaw="< change magnref to "< correction; InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_); long nzmod = ((Nz_>10)?Nz_/10:1); for(long l=0;lDtrcom(zred); // pour variation angle solide double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred)); double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu // on a: Mass ~ DNu * Dlum^2 / Dtrcom^2 double corr = dnu/dnu_ref_*pow(dtrc_ref_/dtrc*dlum/dlum_ref_,2.); // F(nu) = F(1.4GHz)*(nu GHz/1.4 Ghz)^p = F(1.4GHz)*(1/(1+z))^p , car nu = 1.4 GHz/(1+z) if(powlaw!=0.) corr *= pow((1.+redsh_ref_)/(1.+zred),powlaw); correction.push_back(corr); if(lp_>0 && (l==0 || l==Nz_-1 || l%nzmod==0)) { cout<<"l="< cor="<0 && nsum>1.) { sum /= nsum; sum2 = sum2/nsum - sum*sum; cout<<"...Mean mass="<