#include "sopnamsp.h" #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 FFTW_THREAD #define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag())) //------------------------------------------------------- GeneFluct3D::GeneFluct3D(TArray< complex >& T) : T_(T) , Nx_(0) , Ny_(0) , Nz_(0) , array_allocated_(false) , lp_(0) , redshref_(-999.) , kredshref_(0.) , cosmo_(NULL) , growth_(NULL) , loscom_ref_(-999.), loscom_min_(-999.), loscom_max_(-999.) , loscom2zred_min_(0.) , loscom2zred_max_(0.) { xobs_[0] = xobs_[1] = xobs_[2] = 0.; zred_.resize(0); loscom_.resize(0); loscom2zred_.resize(0); SetNThread(); } GeneFluct3D::~GeneFluct3D(void) { fftw_destroy_plan(pf_); fftw_destroy_plan(pb_); #ifdef FFTW_THREAD if(nthread_>0) fftw_cleanup_threads(); #endif } //------------------------------------------------------- void GeneFluct3D::SetSize(long nx,long ny,long nz,double dx,double dy,double dz) { setsize(nx,ny,nz,dx,dy,dz); setalloc(); setpointers(false); init_fftw(); } void GeneFluct3D::SetObservator(double redshref,double kredshref) // L'observateur est au redshift z=0 // est situe sur la "perpendiculaire" a la face x,y // issue du centre de cette face // Il faut positionner le cube sur l'axe des z cad des redshifts: // redshref = redshift de reference // Si redshref<0 alors redshref=0 // kredshref = indice (en double) correspondant a ce redshift // Si kredshref<0 alors kredshref=0 // Exemple: redshref=1.5 kredshref=250.75 // -> Le pixel i=nx/2 j=ny/2 k=250.75 est au redshift 1.5 { if(redshref<0.) redshref = 0.; if(kredshref<0.) kredshref = 0.; redshref_ = redshref; kredshref_ = kredshref; if(lp_>0) cout<<"--- GeneFluct3D::SetObservator zref="<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"< zred // { if(lp_>0) cout<<"--- LosComRedshift: zinc="<Dloscom(redshref_); xobs_[0] = Nx_/2.*Dx_; xobs_[1] = Ny_/2.*Dy_; xobs_[2] = kredshref_*Dz_ - loscom_ref_; // L'observateur est-il dans le cube? bool obs_in_cube = false; if(xobs_[2]>=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]; // Find MAXIMUM los com distance to the observer: // ou que soit positionne l'observateur, la distance // maximal est sur un des coins du cube loscom_max_ = 0.; for(long i=0;i<=1;i++) { double dx2 = xobs_[0] - i*Nx_*Dx_; dx2 *= dx2; for(long j=0;j<=1;j++) { double dy2 = xobs_[1] - j*Ny_*Dy_; dy2 *= dy2; for(long k=0;k<=1;k++) { double dz2 = xobs_[2] - k*Nz_*Dz_; dz2 *= dz2; dz2 = sqrt(dx2+dy2+dz2); if(dz2>loscom_max_) loscom_max_ = dz2; } } } if(lp_>0) { cout<<"...zref="<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 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); } catch (PThrowable & exc) { cout<<"Exception : "<<(string)typeid(exc).name() <<" - Msg= "<0) 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 ---"<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<<"--- VarianceFrReal ---"<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++; } return nbad; } int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax) // mean,sigma^2 pour pixels avec valeurs ]vmin,vmax[ extremites exclues // -> mean and sigma^2 are NOT computed with pixels values vmin and vmax { check_array_alloc(); int_8 n = 0; rm = rs2 = 0.; for(long i=0;i=vmax) continue; rm += v; rs2 += v*v; n++; } if(n>1) { rm /= (double)n; rs2 = rs2/(double)n - rm*rm; } return n; } int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0) // set to "val0" if out of range ]vmin,vmax[ extremites exclues // -> vmin and vmax are set to val0 { check_array_alloc(); int_8 nbad = 0; for(long i=0;i=vmax) {data_[ip] = val0; nbad++;} } return nbad; } //------------------------------------------------------- void GeneFluct3D::TurnFluct2Mass(void) // d_rho/rho -> Mass (add one!) { if(lp_>0) cout<<"--- TurnFluct2Mass ---"<0) cout<<"--- TurnMass2MeanNumber ---"<0) cout<<"...ngood="< "<0. // On NE met PAS a zero les pixels <0 // On renormalise sur les pixels>0 pour qu'on ait n*Vol galaxies // 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) double dn = n_by_mpc3*Vol_/m /(double)ngood; // nb de gal a mettre ds 1 px if(lp_>0) cout<<"...galaxy density move from " <0.) sum += data_[ip]; // mais on ne les compte pas } if(lp_>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< 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 en loi de puissance IDENTIQUE 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 dans un pixel par les AGN // lsigma : c'est le sigma de la distribution // powlaw : c'est la pente de ls 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: = "<Dang(redshref_); double dlumref = cosmo_->Dlum(redshref_); double dredref = Dz_/(cosmo_->Dhubble()/cosmo_->E(redshref_)); double dnuref = Fr_HyperFin_Par *dredref/pow(1.+redshref_,2.); // GHz double fagnref = pow(10.,lfjy)*(dnuref*1.e9); // Jy.Hz double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlumref); // Msol if(lp_>0) { cout<<"Au centre: z="<0) cout<<" powlaw="< change magnref to "< correction; InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_); for(long l=0;lDang(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 = dnu/dnuref*pow(dangref/dang*dlum/dlumref,2.); if(powlaw!=0.) corr *= pow((1.+redshref_)/(1.+zred),powlaw); correction.push_back(corr); if(lp_>0 && (l==0 || abs(l-Nz_/2)<2 || l==Nz_-1)) { cout<<"l="< corr="<0 && nsum>1.) { sum /= nsum; sum2 = sum2/nsum - sum*sum; cout<<"...Mean mass="<0) cout<<"--- ComputeSpectrum_bricolo ---"<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=1;l0) cout<<"--- ComputeSpectrum2D_bricolo ---"<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=1;l