[3115] | 1 | #include "machdefs.h"
|
---|
| 2 | #include <iostream>
|
---|
| 3 | #include <stdlib.h>
|
---|
| 4 | #include <stdio.h>
|
---|
| 5 | #include <string.h>
|
---|
| 6 | #include <math.h>
|
---|
| 7 | #include <unistd.h>
|
---|
| 8 |
|
---|
| 9 | #include "tarray.h"
|
---|
| 10 | #include "pexceptions.h"
|
---|
| 11 | #include "perandom.h"
|
---|
| 12 | #include "srandgen.h"
|
---|
| 13 |
|
---|
[3141] | 14 | #include "fabtcolread.h"
|
---|
| 15 | #include "fabtwriter.h"
|
---|
| 16 | #include "fioarr.h"
|
---|
[3524] | 17 | #include "ntuple.h"
|
---|
[3141] | 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 |
|
---|
[3518] | 27 | #if defined(GEN3D_FLOAT)
|
---|
| 28 | #define GEN3D_FFTW_INIT_THREADS fftwf_init_threads
|
---|
| 29 | #define GEN3D_FFTW_CLEANUP_THREADS fftwf_cleanup_threads
|
---|
| 30 | #define GEN3D_FFTW_PLAN_WITH_NTHREADS fftwf_plan_with_nthreads
|
---|
| 31 | #define GEN3D_FFTW_PLAN_DFT_R2C_3D fftwf_plan_dft_r2c_3d
|
---|
| 32 | #define GEN3D_FFTW_PLAN_DFT_C2R_3D fftwf_plan_dft_c2r_3d
|
---|
| 33 | #define GEN3D_FFTW_DESTROY_PLAN fftwf_destroy_plan
|
---|
| 34 | #define GEN3D_FFTW_EXECUTE fftwf_execute
|
---|
| 35 | #else
|
---|
| 36 | #define GEN3D_FFTW_INIT_THREADS fftw_init_threads
|
---|
| 37 | #define GEN3D_FFTW_CLEANUP_THREADS fftw_cleanup_threads
|
---|
| 38 | #define GEN3D_FFTW_PLAN_WITH_NTHREADS fftw_plan_with_nthreads
|
---|
| 39 | #define GEN3D_FFTW_PLAN_DFT_R2C_3D fftw_plan_dft_r2c_3d
|
---|
| 40 | #define GEN3D_FFTW_PLAN_DFT_C2R_3D fftw_plan_dft_c2r_3d
|
---|
| 41 | #define GEN3D_FFTW_DESTROY_PLAN fftw_destroy_plan
|
---|
| 42 | #define GEN3D_FFTW_EXECUTE fftw_execute
|
---|
| 43 | #endif
|
---|
[3115] | 44 |
|
---|
| 45 | #define MODULE2(_x_) ((double)((_x_).real()*(_x_).real() + (_x_).imag()*(_x_).imag()))
|
---|
| 46 |
|
---|
[3325] | 47 | namespace SOPHYA {
|
---|
| 48 |
|
---|
[3115] | 49 | //-------------------------------------------------------
|
---|
[3349] | 50 | GeneFluct3D::GeneFluct3D(long nx,long ny,long nz,double dx,double dy,double dz
|
---|
| 51 | ,unsigned short nthread,int lp)
|
---|
[3115] | 52 | {
|
---|
[3349] | 53 | init_default();
|
---|
[3115] | 54 |
|
---|
[3349] | 55 | lp_ = lp;
|
---|
| 56 | nthread_ = nthread;
|
---|
[3115] | 57 |
|
---|
[3349] | 58 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
| 59 | setalloc();
|
---|
| 60 | setpointers(false);
|
---|
| 61 | init_fftw();
|
---|
[3141] | 62 | }
|
---|
| 63 |
|
---|
[3349] | 64 | GeneFluct3D::GeneFluct3D(unsigned short nthread)
|
---|
[3154] | 65 | {
|
---|
[3349] | 66 | init_default();
|
---|
| 67 | setsize(2,2,2,1.,1.,1.);
|
---|
| 68 | nthread_ = nthread;
|
---|
| 69 | setalloc();
|
---|
| 70 | setpointers(false);
|
---|
| 71 | init_fftw();
|
---|
[3154] | 72 | }
|
---|
| 73 |
|
---|
[3349] | 74 | GeneFluct3D::~GeneFluct3D(void)
|
---|
[3157] | 75 | {
|
---|
[3349] | 76 | delete_fftw();
|
---|
[3157] | 77 | }
|
---|
| 78 |
|
---|
[3349] | 79 | //-------------------------------------------------------
|
---|
| 80 | void GeneFluct3D::init_default(void)
|
---|
[3157] | 81 | {
|
---|
[3349] | 82 | Nx_ = Ny_ = Nz_ = 0;
|
---|
[3518] | 83 | is_set_fft_plan = false;
|
---|
[3349] | 84 | nthread_ = 0;
|
---|
| 85 | lp_ = 0;
|
---|
[3524] | 86 | array_allocated_ = false; array_type = 0;
|
---|
[3349] | 87 | cosmo_ = NULL;
|
---|
| 88 | growth_ = NULL;
|
---|
[3768] | 89 | good_dzinc_ = 0.01;
|
---|
| 90 | compute_pk_redsh_ref_ = -999.;
|
---|
[3349] | 91 | redsh_ref_ = -999.;
|
---|
| 92 | kredsh_ref_ = 0.;
|
---|
| 93 | dred_ref_ = -999.;
|
---|
| 94 | loscom_ref_ = -999.;
|
---|
| 95 | dtrc_ref_ = dlum_ref_ = dang_ref_ = -999.;
|
---|
| 96 | nu_ref_ = dnu_ref_ = -999.;
|
---|
| 97 | loscom_min_ = loscom_max_ = -999.;
|
---|
| 98 | loscom2zred_min_ = loscom2zred_max_ = 0.;
|
---|
| 99 | xobs_[0] = xobs_[1] = xobs_[2] = 0.;
|
---|
| 100 | zred_.resize(0);
|
---|
| 101 | loscom_.resize(0);
|
---|
| 102 | loscom2zred_.resize(0);
|
---|
[3157] | 103 | }
|
---|
| 104 |
|
---|
[3141] | 105 | void GeneFluct3D::setsize(long nx,long ny,long nz,double dx,double dy,double dz)
|
---|
| 106 | {
|
---|
[3155] | 107 | if(lp_>1) cout<<"--- GeneFluct3D::setsize: N="<<nx<<","<<ny<<","<<nz
|
---|
| 108 | <<" D="<<dx<<","<<dy<<","<<dz<<endl;
|
---|
[3141] | 109 | if(nx<=0 || dx<=0.) {
|
---|
[3572] | 110 | const char *bla = "GeneFluct3D::setsize_Error: bad value(s) for nn/dx";
|
---|
[3199] | 111 | cout<<bla<<endl; throw ParmError(bla);
|
---|
[3115] | 112 | }
|
---|
| 113 |
|
---|
[3141] | 114 | // Les tailles des tableaux
|
---|
[3115] | 115 | Nx_ = nx;
|
---|
| 116 | Ny_ = ny; if(Ny_ <= 0) Ny_ = Nx_;
|
---|
| 117 | Nz_ = nz; if(Nz_ <= 0) Nz_ = Nx_;
|
---|
[3141] | 118 | N_.resize(0); N_.push_back(Nx_); N_.push_back(Ny_); N_.push_back(Nz_);
|
---|
[3115] | 119 | NRtot_ = Nx_*Ny_*Nz_; // nombre de pixels dans le survey
|
---|
| 120 | NCz_ = Nz_/2 +1;
|
---|
| 121 | NTz_ = 2*NCz_;
|
---|
| 122 |
|
---|
| 123 | // le pas dans l'espace (Mpc)
|
---|
| 124 | Dx_ = dx;
|
---|
| 125 | Dy_ = dy; if(Dy_ <= 0.) Dy_ = Dx_;
|
---|
| 126 | Dz_ = dz; if(Dz_ <= 0.) Dz_ = Dx_;
|
---|
[3141] | 127 | D_.resize(0); D_.push_back(Dx_); D_.push_back(Dy_); D_.push_back(Dz_);
|
---|
[3115] | 128 | dVol_ = Dx_*Dy_*Dz_;
|
---|
| 129 | Vol_ = (Nx_*Dx_)*(Ny_*Dy_)*(Nz_*Dz_);
|
---|
| 130 |
|
---|
| 131 | // Le pas dans l'espace de Fourier (Mpc^-1)
|
---|
| 132 | Dkx_ = 2.*M_PI/(Nx_*Dx_);
|
---|
| 133 | Dky_ = 2.*M_PI/(Ny_*Dy_);
|
---|
| 134 | Dkz_ = 2.*M_PI/(Nz_*Dz_);
|
---|
[3141] | 135 | Dk_.resize(0); Dk_.push_back(Dkx_); Dk_.push_back(Dky_); Dk_.push_back(Dkz_);
|
---|
[3115] | 136 | Dk3_ = Dkx_*Dky_*Dkz_;
|
---|
| 137 |
|
---|
| 138 | // La frequence de Nyquist en k (Mpc^-1)
|
---|
| 139 | Knyqx_ = M_PI/Dx_;
|
---|
| 140 | Knyqy_ = M_PI/Dy_;
|
---|
| 141 | Knyqz_ = M_PI/Dz_;
|
---|
[3141] | 142 | Knyq_.resize(0); Knyq_.push_back(Knyqx_); Knyq_.push_back(Knyqy_); Knyq_.push_back(Knyqz_);
|
---|
| 143 | }
|
---|
[3115] | 144 |
|
---|
[3141] | 145 | void GeneFluct3D::setalloc(void)
|
---|
| 146 | {
|
---|
[3521] | 147 | #if defined(GEN3D_FLOAT)
|
---|
| 148 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc FLOAT ---"<<endl;
|
---|
| 149 | #else
|
---|
| 150 | if(lp_>1) cout<<"--- GeneFluct3D::setalloc DOUBLE ---"<<endl;
|
---|
| 151 | #endif
|
---|
[3141] | 152 | // Dimensionnement du tableau complex<r_8>
|
---|
| 153 | // ATTENTION: TArray adresse en memoire a l'envers du C
|
---|
| 154 | // Tarray(n1,n2,n3) == Carray[n3][n2][n1]
|
---|
| 155 | sa_size_t SzK_[3] = {NCz_,Ny_,Nx_}; // a l'envers
|
---|
| 156 | try {
|
---|
[3743] | 157 | T_.ReSize(3,SzK_);
|
---|
[3524] | 158 | array_allocated_ = true; array_type=0;
|
---|
[3518] | 159 | if(lp_>1) cout<<" allocating: "<<T_.Size()*sizeof(complex<GEN3D_TYPE>)/1.e6<<" Mo"<<endl;
|
---|
[3141] | 160 | } catch (...) {
|
---|
[3155] | 161 | cout<<"GeneFluct3D::setalloc_Error: Problem allocating T_"<<endl;
|
---|
[3141] | 162 | }
|
---|
[3115] | 163 | }
|
---|
| 164 |
|
---|
[3141] | 165 | void GeneFluct3D::setpointers(bool from_real)
|
---|
| 166 | {
|
---|
[3155] | 167 | if(lp_>1) cout<<"--- GeneFluct3D::setpointers ---"<<endl;
|
---|
[3141] | 168 | if(from_real) T_ = ArrCastR2C(R_);
|
---|
| 169 | else R_ = ArrCastC2R(T_);
|
---|
| 170 | // On remplit les pointeurs
|
---|
[3518] | 171 | fdata_ = (GEN3D_FFTW_COMPLEX *) (&T_(0,0,0));
|
---|
| 172 | data_ = (GEN3D_TYPE *) (&R_(0,0,0));
|
---|
[3141] | 173 | }
|
---|
| 174 |
|
---|
[3349] | 175 | void GeneFluct3D::init_fftw(void)
|
---|
[3141] | 176 | {
|
---|
[3518] | 177 | if( is_set_fft_plan ) delete_fftw();
|
---|
[3141] | 178 |
|
---|
[3154] | 179 | // --- Initialisation de fftw3 (attention data est sur-ecrit a l'init)
|
---|
[3155] | 180 | if(lp_>1) cout<<"--- GeneFluct3D::init_fftw ---"<<endl;
|
---|
[3349] | 181 | #ifdef WITH_FFTW_THREAD
|
---|
[3154] | 182 | if(nthread_>0) {
|
---|
[3155] | 183 | cout<<"...Computing with "<<nthread_<<" threads"<<endl;
|
---|
[3518] | 184 | GEN3D_FFTW_INIT_THREADS();
|
---|
| 185 | GEN3D_FFTW_PLAN_WITH_NTHREADS(nthread_);
|
---|
[3154] | 186 | }
|
---|
| 187 | #endif
|
---|
[3155] | 188 | if(lp_>1) cout<<"...forward plan"<<endl;
|
---|
[3518] | 189 | pf_ = GEN3D_FFTW_PLAN_DFT_R2C_3D(Nx_,Ny_,Nz_,data_,fdata_,FFTW_ESTIMATE);
|
---|
[3155] | 190 | if(lp_>1) cout<<"...backward plan"<<endl;
|
---|
[3518] | 191 | pb_ = GEN3D_FFTW_PLAN_DFT_C2R_3D(Nx_,Ny_,Nz_,fdata_,data_,FFTW_ESTIMATE);
|
---|
| 192 | is_set_fft_plan = true;
|
---|
[3154] | 193 | }
|
---|
[3141] | 194 |
|
---|
[3349] | 195 | void GeneFluct3D::delete_fftw(void)
|
---|
| 196 | {
|
---|
[3518] | 197 | if( !is_set_fft_plan ) return;
|
---|
| 198 | GEN3D_FFTW_DESTROY_PLAN(pf_);
|
---|
| 199 | GEN3D_FFTW_DESTROY_PLAN(pb_);
|
---|
[3349] | 200 | #ifdef WITH_FFTW_THREAD
|
---|
[3518] | 201 | if(nthread_>0) GEN3D_FFTW_CLEANUP_THREADS();
|
---|
[3349] | 202 | #endif
|
---|
[3518] | 203 | is_set_fft_plan = false;
|
---|
[3349] | 204 | }
|
---|
| 205 |
|
---|
| 206 | void GeneFluct3D::check_array_alloc(void)
|
---|
| 207 | // Pour tester si le tableau T_ est alloue
|
---|
| 208 | {
|
---|
| 209 | if(array_allocated_) return;
|
---|
| 210 | char bla[90];
|
---|
| 211 | sprintf(bla,"GeneFluct3D::check_array_alloc_Error: array is not allocated");
|
---|
| 212 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 213 | }
|
---|
| 214 |
|
---|
[3157] | 215 | //-------------------------------------------------------
|
---|
[3349] | 216 | void GeneFluct3D::SetObservator(double redshref,double kredshref)
|
---|
| 217 | // L'observateur est au redshift z=0
|
---|
| 218 | // est situe sur la "perpendiculaire" a la face x,y
|
---|
| 219 | // issue du centre de cette face
|
---|
| 220 | // Il faut positionner le cube sur l'axe des z cad des redshifts:
|
---|
| 221 | // redshref = redshift de reference
|
---|
| 222 | // Si redshref<0 alors redshref=0
|
---|
| 223 | // kredshref = indice (en double) correspondant a ce redshift
|
---|
| 224 | // Si kredshref<0 alors kredshref=nz/2 (milieu du cube)
|
---|
| 225 | // Exemple: redshref=1.5 kredshref=250.75
|
---|
| 226 | // -> Le pixel i=nx/2 j=ny/2 k=250.75 est au redshift 1.5
|
---|
| 227 | {
|
---|
| 228 | if(redshref<0.) redshref = 0.;
|
---|
| 229 | if(kredshref<0.) {
|
---|
| 230 | if(Nz_<=0) {
|
---|
[3572] | 231 | const char *bla = "GeneFluct3D::SetObservator_Error: for kredsh_ref<0 define cube geometry first";
|
---|
[3349] | 232 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 233 | }
|
---|
| 234 | kredshref = Nz_/2.;
|
---|
| 235 | }
|
---|
| 236 | redsh_ref_ = redshref;
|
---|
| 237 | kredsh_ref_ = kredshref;
|
---|
| 238 | if(lp_>0)
|
---|
| 239 | cout<<"--- GeneFluct3D::SetObservator zref="<<redsh_ref_<<" kref="<<kredsh_ref_<<endl;
|
---|
| 240 | }
|
---|
| 241 |
|
---|
| 242 | void GeneFluct3D::SetCosmology(CosmoCalc& cosmo)
|
---|
| 243 | {
|
---|
| 244 | cosmo_ = &cosmo;
|
---|
| 245 | if(lp_>1) cosmo_->Print();
|
---|
| 246 | }
|
---|
| 247 |
|
---|
| 248 | void GeneFluct3D::SetGrowthFactor(GrowthFactor& growth)
|
---|
| 249 | {
|
---|
| 250 | growth_ = &growth;
|
---|
| 251 | }
|
---|
| 252 |
|
---|
[3199] | 253 | long GeneFluct3D::LosComRedshift(double zinc,long npoints)
|
---|
[3157] | 254 | // Given a position of the cube relative to the observer
|
---|
| 255 | // and a cosmology
|
---|
| 256 | // (SetObservator() and SetCosmology() should have been called !)
|
---|
| 257 | // This routine filled:
|
---|
| 258 | // the vector "zred_" of scanned redshift (by zinc increments)
|
---|
| 259 | // the vector "loscom_" of corresponding los comoving distance
|
---|
[3199] | 260 | // -- Input:
|
---|
| 261 | // zinc : redshift increment for computation
|
---|
| 262 | // npoints : number of points required for inverting loscom -> zred
|
---|
[3157] | 263 | //
|
---|
| 264 | {
|
---|
[3768] | 265 | if(zinc<=0.) zinc = good_dzinc_;
|
---|
[3199] | 266 | if(lp_>0) cout<<"--- LosComRedshift: zinc="<<zinc<<" , npoints="<<npoints<<endl;
|
---|
[3154] | 267 |
|
---|
[3768] | 268 | if(cosmo_ == NULL || growth_==NULL || redsh_ref_<0.) {
|
---|
| 269 | const char *bla = "GeneFluct3D::LosComRedshift_Error: set Observator, Cosmology and Growth first";
|
---|
[3199] | 270 | cout<<bla<<endl; throw ParmError(bla);
|
---|
[3157] | 271 | }
|
---|
| 272 |
|
---|
[3271] | 273 | // La distance angulaire/luminosite/Dnu au pixel de reference
|
---|
| 274 | dred_ref_ = Dz_/(cosmo_->Dhubble()/cosmo_->E(redsh_ref_));
|
---|
| 275 | loscom_ref_ = cosmo_->Dloscom(redsh_ref_);
|
---|
| 276 | dtrc_ref_ = cosmo_->Dtrcom(redsh_ref_);
|
---|
| 277 | dlum_ref_ = cosmo_->Dlum(redsh_ref_);
|
---|
| 278 | dang_ref_ = cosmo_->Dang(redsh_ref_);
|
---|
[3768] | 279 | h_ref = cosmo_->H(redsh_ref_);
|
---|
| 280 | growth_ref_ = (*growth_)(redsh_ref_);
|
---|
| 281 | dsdz_growth_ref_ = growth_->DsDz(redsh_ref_,zinc);
|
---|
[3271] | 282 | nu_ref_ = Fr_HyperFin_Par/(1.+redsh_ref_); // GHz
|
---|
| 283 | dnu_ref_ = Fr_HyperFin_Par *dred_ref_/pow(1.+redsh_ref_,2.); // GHz
|
---|
| 284 | if(lp_>0) {
|
---|
| 285 | cout<<"...reference pixel redshref="<<redsh_ref_
|
---|
| 286 | <<", dredref="<<dred_ref_
|
---|
| 287 | <<", nuref="<<nu_ref_ <<" GHz"
|
---|
| 288 | <<", dnuref="<<dnu_ref_ <<" GHz"<<endl
|
---|
[3768] | 289 | <<" H="<<h_ref<<" km/s/Mpc"
|
---|
| 290 | <<", D="<<growth_ref_
|
---|
| 291 | <<", dD/dz="<<dsdz_growth_ref_<<endl
|
---|
[3271] | 292 | <<" dlosc="<<loscom_ref_<<" Mpc com"
|
---|
| 293 | <<", dtrc="<<dtrc_ref_<<" Mpc com"
|
---|
| 294 | <<", dlum="<<dlum_ref_<<" Mpc"
|
---|
| 295 | <<", dang="<<dang_ref_<<" Mpc"<<endl;
|
---|
| 296 | }
|
---|
| 297 |
|
---|
[3199] | 298 | // On calcule les coordonnees de l'observateur dans le repere du cube
|
---|
| 299 | // cad dans le repere ou l'origine est au centre du pixel i=j=l=0.
|
---|
| 300 | // L'observateur est sur un axe centre sur le milieu de la face Oxy
|
---|
[3157] | 301 | xobs_[0] = Nx_/2.*Dx_;
|
---|
| 302 | xobs_[1] = Ny_/2.*Dy_;
|
---|
[3271] | 303 | xobs_[2] = kredsh_ref_*Dz_ - loscom_ref_;
|
---|
[3157] | 304 |
|
---|
| 305 | // L'observateur est-il dans le cube?
|
---|
| 306 | bool obs_in_cube = false;
|
---|
| 307 | if(xobs_[2]>=0. && xobs_[2]<=Nz_*Dz_) obs_in_cube = true;
|
---|
| 308 |
|
---|
| 309 | // Find MINIMUM los com distance to the observer:
|
---|
| 310 | // c'est le centre de la face a k=0
|
---|
| 311 | // (ou zero si l'observateur est dans le cube)
|
---|
| 312 | loscom_min_ = 0.;
|
---|
| 313 | if(!obs_in_cube) loscom_min_ = -xobs_[2];
|
---|
| 314 |
|
---|
[3271] | 315 | // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED
|
---|
| 316 | if(loscom_min_<=1.e-50)
|
---|
[3746] | 317 | for(int i=0;i<10;i++)
|
---|
[3271] | 318 | cout<<"ATTENTION TOUTES LES PARTIES DU CODE NE MARCHENT PAS POUR UN OBSERVATEUR DANS LE CUBE"<<endl;
|
---|
| 319 | // TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED TO BE FIXED
|
---|
| 320 |
|
---|
| 321 |
|
---|
[3157] | 322 | // Find MAXIMUM los com distance to the observer:
|
---|
| 323 | // ou que soit positionne l'observateur, la distance
|
---|
| 324 | // maximal est sur un des coins du cube
|
---|
| 325 | loscom_max_ = 0.;
|
---|
| 326 | for(long i=0;i<=1;i++) {
|
---|
[3271] | 327 | double dx2 = DXcom(i*(Nx_-1)); dx2 *= dx2;
|
---|
[3157] | 328 | for(long j=0;j<=1;j++) {
|
---|
[3271] | 329 | double dy2 = DYcom(j*(Ny_-1)); dy2 *= dy2;
|
---|
[3157] | 330 | for(long k=0;k<=1;k++) {
|
---|
[3271] | 331 | double dz2 = DZcom(k*(Nz_-1)); dz2 *= dz2;
|
---|
[3157] | 332 | dz2 = sqrt(dx2+dy2+dz2);
|
---|
| 333 | if(dz2>loscom_max_) loscom_max_ = dz2;
|
---|
| 334 | }
|
---|
| 335 | }
|
---|
| 336 | }
|
---|
| 337 | if(lp_>0) {
|
---|
[3271] | 338 | cout<<"...zref="<<redsh_ref_<<" kzref="<<kredsh_ref_<<" losref="<<loscom_ref_<<" Mpc\n"
|
---|
[3157] | 339 | <<" xobs="<<xobs_[0]<<" , "<<xobs_[1]<<" , "<<xobs_[2]<<" Mpc "
|
---|
| 340 | <<" in_cube="<<obs_in_cube
|
---|
[3271] | 341 | <<" loscom_min="<<loscom_min_<<" loscom_max="<<loscom_max_<<" Mpc (com)"<<endl;
|
---|
[3157] | 342 | }
|
---|
| 343 |
|
---|
[3199] | 344 | // Fill the corresponding vectors for loscom and zred
|
---|
[3768] | 345 | // Be shure to have one dlc<loscom_min and one dlc>loscom_max
|
---|
[3746] | 346 | double zmin = 0., dlcmin=0.;
|
---|
| 347 | while(1) {
|
---|
| 348 | if(lp_>0)
|
---|
| 349 | cout<<"...Filling zred starting at zmin="<<zmin<<" with zinc="<<zinc
|
---|
| 350 | <<", loscom_min-max=["<<loscom_min_<<","<<loscom_max_<<"]"<<endl;
|
---|
| 351 | zred_.resize(0); loscom_.resize(0);
|
---|
| 352 | for(double z=zmin; ; z+=zinc) {
|
---|
| 353 | double dlc = cosmo_->Dloscom(z);
|
---|
| 354 | if(dlc<loscom_min_) {
|
---|
| 355 | dlcmin = dlc;
|
---|
| 356 | zmin = z;
|
---|
| 357 | zred_.resize(0); loscom_.resize(0);
|
---|
| 358 | }
|
---|
| 359 | zred_.push_back(z);
|
---|
| 360 | loscom_.push_back(dlc);
|
---|
| 361 | z += zinc;
|
---|
| 362 | if(dlc>loscom_max_) {
|
---|
| 363 | if(lp_>0)
|
---|
| 364 | cout<<" Min: z="<<zmin<<" dlc="<<dlcmin<<", Max: z="<<z<<" dlc="<<dlc<<endl;
|
---|
| 365 | break; // on sort apres avoir stoque un dlc>dlcmax
|
---|
| 366 | }
|
---|
| 367 | }
|
---|
| 368 | if(zred_.size()>=10) break;
|
---|
| 369 | zinc /= 10.;
|
---|
| 370 | cout<<" not enough points ("<<zred_.size()<<") for zref, retry with zinc="<<zinc<<endl;
|
---|
[3157] | 371 | }
|
---|
| 372 |
|
---|
| 373 | if(lp_>0) {
|
---|
[3199] | 374 | long n = zred_.size();
|
---|
| 375 | cout<<"...zred/loscom tables[zinc="<<zinc<<"]: n="<<n;
|
---|
[3157] | 376 | if(n>0) cout<<" z="<<zred_[0]<<" -> d="<<loscom_[0];
|
---|
| 377 | if(n>1) cout<<" , z="<<zred_[n-1]<<" -> d="<<loscom_[n-1];
|
---|
| 378 | cout<<endl;
|
---|
| 379 | }
|
---|
| 380 |
|
---|
[3199] | 381 | // Compute the parameters and tables needed for inversion loscom->zred
|
---|
| 382 | if(npoints<3) npoints = zred_.size();
|
---|
| 383 | InverseFunc invfun(zred_,loscom_);
|
---|
| 384 | invfun.ComputeParab(npoints,loscom2zred_);
|
---|
| 385 | loscom2zred_min_ = invfun.YMin();
|
---|
| 386 | loscom2zred_max_ = invfun.YMax();
|
---|
| 387 |
|
---|
| 388 | if(lp_>0) {
|
---|
| 389 | long n = loscom2zred_.size();
|
---|
| 390 | cout<<"...loscom -> zred[npoints="<<npoints<<"]: n="<<n
|
---|
| 391 | <<" los_min="<<loscom2zred_min_
|
---|
| 392 | <<" los_max="<<loscom2zred_max_
|
---|
| 393 | <<" -> zred=[";
|
---|
| 394 | if(n>0) cout<<loscom2zred_[0];
|
---|
| 395 | cout<<",";
|
---|
| 396 | if(n>1) cout<<loscom2zred_[n-1];
|
---|
| 397 | cout<<"]"<<endl;
|
---|
| 398 | if(lp_>1 && n>0)
|
---|
| 399 | for(int i=0;i<n;i++)
|
---|
[3330] | 400 | if(i<2 || abs(i-n/2)<2 || i>=n-2)
|
---|
| 401 | cout<<" i="<<i
|
---|
| 402 | <<" d="<<loscom2zred_min_+i*(loscom2zred_max_-loscom2zred_min_)/(n-1.)
|
---|
| 403 | <<" Mpc z="<<loscom2zred_[i]<<endl;
|
---|
[3199] | 404 | }
|
---|
| 405 |
|
---|
[3768] | 406 |
|
---|
| 407 | // Compute the table for D'(z)/D(z) where D'(z)=dD/dEta (Eta is conformal time)
|
---|
| 408 | zredmin_dpsd_ = zred_[0];
|
---|
| 409 | zredmax_dpsd_ = zred_[zred_.size()-1];
|
---|
| 410 | long nptd = long(sqrt(Nx_*Nx_ + Ny_*Ny_ +Nz_*Nz_));
|
---|
| 411 | if(nptd<10) nptd = 10;
|
---|
| 412 | good_dzinc_ = (zredmax_dpsd_ - zredmin_dpsd_)/nptd;
|
---|
| 413 | if(lp_>0) cout<<"...good_dzinc changed to "<<good_dzinc_<<endl;
|
---|
| 414 | dpsdfrzred_.resize(0);
|
---|
| 415 | double dz = (zredmax_dpsd_ - zredmin_dpsd_)/(nptd-1);
|
---|
| 416 | int nmod = nptd/5; if(nmod==0) nmod = 1;
|
---|
| 417 | if(lp_>0) cout<<"...Compute table D'/D on "<<nptd<<" pts for z["
|
---|
| 418 | <<zredmin_dpsd_<<","<<zredmax_dpsd_<<"]"<<endl;
|
---|
| 419 | for(int i=0;i<nptd;i++) {
|
---|
| 420 | double z = zredmin_dpsd_ + i*dz;
|
---|
| 421 | double v = -cosmo_->H(z) * growth_->DsDz(z,good_dzinc_) / (*growth_)(z);
|
---|
| 422 | dpsdfrzred_.push_back(v);
|
---|
| 423 | if(lp_ && (i%nmod==0 || i==nptd-1)) cout<<" z="<<z<<" D'/D="<<v<<" km/s/Mpc"<<endl;
|
---|
| 424 | }
|
---|
| 425 |
|
---|
[3199] | 426 | return zred_.size();
|
---|
[3157] | 427 | }
|
---|
| 428 |
|
---|
[3115] | 429 | //-------------------------------------------------------
|
---|
[3141] | 430 | void GeneFluct3D::WriteFits(string cfname,int bitpix)
|
---|
| 431 | {
|
---|
[3155] | 432 | cout<<"--- GeneFluct3D::WriteFits: Writing Cube to "<<cfname<<endl;
|
---|
[3141] | 433 | try {
|
---|
| 434 | FitsImg3DWriter fwrt(cfname.c_str(),bitpix,5);
|
---|
| 435 | fwrt.WriteKey("NX",Nx_," axe transverse 1");
|
---|
| 436 | fwrt.WriteKey("NY",Ny_," axe transverse 2");
|
---|
| 437 | fwrt.WriteKey("NZ",Nz_," axe longitudinal (redshift)");
|
---|
| 438 | fwrt.WriteKey("DX",Dx_," Mpc");
|
---|
| 439 | fwrt.WriteKey("DY",Dy_," Mpc");
|
---|
| 440 | fwrt.WriteKey("DZ",Dz_," Mpc");
|
---|
| 441 | fwrt.WriteKey("DKX",Dkx_," Mpc^-1");
|
---|
| 442 | fwrt.WriteKey("DKY",Dky_," Mpc^-1");
|
---|
| 443 | fwrt.WriteKey("DKZ",Dkz_," Mpc^-1");
|
---|
[3768] | 444 | fwrt.WriteKey("ZREFPK",compute_pk_redsh_ref_," Pk computed redshift");
|
---|
[3271] | 445 | fwrt.WriteKey("ZREF",redsh_ref_," reference redshift");
|
---|
| 446 | fwrt.WriteKey("KZREF",kredsh_ref_," reference redshift on z axe");
|
---|
[3768] | 447 | fwrt.WriteKey("Growth",Dref()," growth at reference redshift");
|
---|
| 448 | fwrt.WriteKey("GrowthPk",DrefPk()," growth at Pk computed redshift");
|
---|
[3141] | 449 | fwrt.Write(R_);
|
---|
| 450 | } catch (PThrowable & exc) {
|
---|
| 451 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
| 452 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
| 453 | return;
|
---|
| 454 | } catch (...) {
|
---|
| 455 | cout<<" some other exception was caught !"<<endl;
|
---|
| 456 | return;
|
---|
| 457 | }
|
---|
| 458 | }
|
---|
| 459 |
|
---|
| 460 | void GeneFluct3D::ReadFits(string cfname)
|
---|
| 461 | {
|
---|
[3155] | 462 | cout<<"--- GeneFluct3D::ReadFits: Reading Cube from "<<cfname<<endl;
|
---|
[3141] | 463 | try {
|
---|
| 464 | FitsImg3DRead fimg(cfname.c_str(),0,5);
|
---|
| 465 | fimg.Read(R_);
|
---|
| 466 | long nx = fimg.ReadKeyL("NX");
|
---|
| 467 | long ny = fimg.ReadKeyL("NY");
|
---|
| 468 | long nz = fimg.ReadKeyL("NZ");
|
---|
| 469 | double dx = fimg.ReadKey("DX");
|
---|
| 470 | double dy = fimg.ReadKey("DY");
|
---|
| 471 | double dz = fimg.ReadKey("DZ");
|
---|
[3768] | 472 | double pkzref = fimg.ReadKey("ZREFPK");
|
---|
[3154] | 473 | double zref = fimg.ReadKey("ZREF");
|
---|
| 474 | double kzref = fimg.ReadKey("KZREF");
|
---|
[3141] | 475 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
| 476 | setpointers(true);
|
---|
[3154] | 477 | init_fftw();
|
---|
| 478 | SetObservator(zref,kzref);
|
---|
[3330] | 479 | array_allocated_ = true;
|
---|
[3768] | 480 | compute_pk_redsh_ref_ = pkzref;
|
---|
[3141] | 481 | } catch (PThrowable & exc) {
|
---|
| 482 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
| 483 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
| 484 | return;
|
---|
| 485 | } catch (...) {
|
---|
| 486 | cout<<" some other exception was caught !"<<endl;
|
---|
| 487 | return;
|
---|
| 488 | }
|
---|
| 489 | }
|
---|
| 490 |
|
---|
| 491 | void GeneFluct3D::WritePPF(string cfname,bool write_real)
|
---|
| 492 | // On ecrit soit le TArray<r_8> ou le TArray<complex <r_8> >
|
---|
| 493 | {
|
---|
[3155] | 494 | cout<<"--- GeneFluct3D::WritePPF: Writing Cube (real="<<write_real<<") to "<<cfname<<endl;
|
---|
[3141] | 495 | try {
|
---|
| 496 | R_.Info()["NX"] = (int_8)Nx_;
|
---|
| 497 | R_.Info()["NY"] = (int_8)Ny_;
|
---|
| 498 | R_.Info()["NZ"] = (int_8)Nz_;
|
---|
| 499 | R_.Info()["DX"] = (r_8)Dx_;
|
---|
| 500 | R_.Info()["DY"] = (r_8)Dy_;
|
---|
| 501 | R_.Info()["DZ"] = (r_8)Dz_;
|
---|
[3768] | 502 | R_.Info()["ZREFPK"] = (r_8)compute_pk_redsh_ref_;
|
---|
[3271] | 503 | R_.Info()["ZREF"] = (r_8)redsh_ref_;
|
---|
| 504 | R_.Info()["KZREF"] = (r_8)kredsh_ref_;
|
---|
[3768] | 505 | R_.Info()["Growth"] = (r_8)Dref();
|
---|
| 506 | R_.Info()["GrowthPk"] = (r_8)DrefPk();
|
---|
[3141] | 507 | POutPersist pos(cfname.c_str());
|
---|
| 508 | if(write_real) pos << PPFNameTag("rgen") << R_;
|
---|
| 509 | else pos << PPFNameTag("pkgen") << T_;
|
---|
| 510 | } catch (PThrowable & exc) {
|
---|
| 511 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
| 512 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
| 513 | return;
|
---|
| 514 | } catch (...) {
|
---|
| 515 | cout<<" some other exception was caught !"<<endl;
|
---|
| 516 | return;
|
---|
| 517 | }
|
---|
| 518 | }
|
---|
| 519 |
|
---|
| 520 | void GeneFluct3D::ReadPPF(string cfname)
|
---|
| 521 | {
|
---|
[3155] | 522 | cout<<"--- GeneFluct3D::ReadPPF: Reading Cube from "<<cfname<<endl;
|
---|
[3141] | 523 | try {
|
---|
| 524 | bool from_real = true;
|
---|
| 525 | PInPersist pis(cfname.c_str());
|
---|
| 526 | string name_tag_k = "pkgen";
|
---|
| 527 | bool found_tag_k = pis.GotoNameTag("pkgen");
|
---|
| 528 | if(found_tag_k) {
|
---|
[3262] | 529 | cout<<" ...reading spectrum into TArray<complex <r_8> >"<<endl;
|
---|
[3141] | 530 | pis >> PPFNameTag("pkgen") >> T_;
|
---|
| 531 | from_real = false;
|
---|
| 532 | } else {
|
---|
| 533 | cout<<" ...reading space into TArray<r_8>"<<endl;
|
---|
| 534 | pis >> PPFNameTag("rgen") >> R_;
|
---|
| 535 | }
|
---|
[3154] | 536 | setpointers(from_real); // a mettre ici pour relire les DVInfo
|
---|
[3141] | 537 | int_8 nx = R_.Info()["NX"];
|
---|
| 538 | int_8 ny = R_.Info()["NY"];
|
---|
| 539 | int_8 nz = R_.Info()["NZ"];
|
---|
| 540 | r_8 dx = R_.Info()["DX"];
|
---|
| 541 | r_8 dy = R_.Info()["DY"];
|
---|
| 542 | r_8 dz = R_.Info()["DZ"];
|
---|
[3768] | 543 | r_8 pkzref = R_.Info()["ZREFPK"];
|
---|
[3154] | 544 | r_8 zref = R_.Info()["ZREF"];
|
---|
| 545 | r_8 kzref = R_.Info()["KZREF"];
|
---|
[3141] | 546 | setsize(nx,ny,nz,dx,dy,dz);
|
---|
[3154] | 547 | init_fftw();
|
---|
| 548 | SetObservator(zref,kzref);
|
---|
[3330] | 549 | array_allocated_ = true;
|
---|
[3768] | 550 | compute_pk_redsh_ref_ = pkzref;
|
---|
[3141] | 551 | } catch (PThrowable & exc) {
|
---|
| 552 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
| 553 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
| 554 | return;
|
---|
| 555 | } catch (...) {
|
---|
| 556 | cout<<" some other exception was caught !"<<endl;
|
---|
| 557 | return;
|
---|
| 558 | }
|
---|
| 559 | }
|
---|
| 560 |
|
---|
[3281] | 561 | void GeneFluct3D::WriteSlicePPF(string cfname)
|
---|
| 562 | // On ecrit 3 tranches du cube selon chaque axe
|
---|
| 563 | {
|
---|
[3283] | 564 | cout<<"--- GeneFluct3D::WriteSlicePPF: Writing Cube Slices "<<cfname<<endl;
|
---|
[3281] | 565 | try {
|
---|
| 566 |
|
---|
| 567 | POutPersist pos(cfname.c_str());
|
---|
| 568 | TMatrix<r_4> S;
|
---|
| 569 | char str[16];
|
---|
| 570 | long i,j,l;
|
---|
| 571 |
|
---|
| 572 | // Tranches en Z
|
---|
| 573 | for(int s=0;s<3;s++) {
|
---|
| 574 | S.ReSize(Nx_,Ny_);
|
---|
| 575 | if(s==0) l=0; else if(s==1) l=(Nz_+1)/2; else l=Nz_-1;
|
---|
[3289] | 576 | sprintf(str,"z%ld",l);
|
---|
[3281] | 577 | for(i=0;i<Nx_;i++) for(j=0;j<Ny_;j++) S(i,j)=data_[IndexR(i,j,l)];
|
---|
| 578 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
| 579 | }
|
---|
| 580 |
|
---|
| 581 | // Tranches en Y
|
---|
| 582 | for(int s=0;s<3;s++) {
|
---|
| 583 | S.ReSize(Nz_,Nx_);
|
---|
| 584 | if(s==0) j=0; else if(s==1) j=(Ny_+1)/2; else j=Ny_-1;
|
---|
[3289] | 585 | sprintf(str,"y%ld",j);
|
---|
[3281] | 586 | for(i=0;i<Nx_;i++) for(l=0;l<Nz_;l++) S(l,i)=data_[IndexR(i,j,l)];
|
---|
| 587 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
| 588 | }
|
---|
| 589 |
|
---|
| 590 | // Tranches en X
|
---|
| 591 | for(int s=0;s<3;s++) {
|
---|
| 592 | S.ReSize(Nz_,Ny_);
|
---|
| 593 | if(s==0) i=0; else if(s==1) i=(Nx_+1)/2; else i=Nx_-1;
|
---|
[3289] | 594 | sprintf(str,"x%ld",i);
|
---|
[3281] | 595 | for(j=0;j<Ny_;j++) for(l=0;l<Nz_;l++) S(l,j)=data_[IndexR(i,j,l)];
|
---|
| 596 | pos<<PPFNameTag(str)<<S; S.RenewObjId();
|
---|
| 597 | }
|
---|
| 598 |
|
---|
| 599 | } catch (PThrowable & exc) {
|
---|
| 600 | cout<<"Exception : "<<(string)typeid(exc).name()
|
---|
| 601 | <<" - Msg= "<<exc.Msg()<<endl;
|
---|
| 602 | return;
|
---|
| 603 | } catch (...) {
|
---|
| 604 | cout<<" some other exception was caught !"<<endl;
|
---|
| 605 | return;
|
---|
| 606 | }
|
---|
| 607 | }
|
---|
| 608 |
|
---|
[3141] | 609 | //-------------------------------------------------------
|
---|
[3524] | 610 | void GeneFluct3D::NTupleCheck(POutPersist &pos,string ntname,unsigned long nent)
|
---|
| 611 | // Remplit le NTuple "ntname" avec "nent" valeurs du cube (reel ou complex) et l'ecrit dans "pos"
|
---|
| 612 | {
|
---|
| 613 | if(ntname.size()<=0 || nent==0) return;
|
---|
| 614 | int nvar = 0;
|
---|
| 615 | if(array_type==1) nvar = 3;
|
---|
| 616 | else if(array_type==2) nvar = 4;
|
---|
| 617 | else return;
|
---|
[3572] | 618 | const char *vname[4] = {"t","z","re","im"};
|
---|
[3524] | 619 | float xnt[4];
|
---|
| 620 | NTuple nt(nvar,vname);
|
---|
| 621 |
|
---|
| 622 | if(array_type==1) {
|
---|
| 623 | unsigned long nmod = Nx_*Ny_*Nz_/nent; if(nmod==0) nmod=1;
|
---|
| 624 | unsigned long n=0;
|
---|
| 625 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 626 | if(n==nmod) {
|
---|
| 627 | int_8 ip = IndexR(i,j,l);
|
---|
| 628 | xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=data_[ip];
|
---|
| 629 | nt.Fill(xnt);
|
---|
| 630 | n=0;
|
---|
| 631 | }
|
---|
| 632 | n++;
|
---|
| 633 | }
|
---|
| 634 | } else {
|
---|
| 635 | unsigned long nmod = Nx_*Ny_*NCz_/nent; if(nmod==0) nmod=1;
|
---|
| 636 | unsigned long n=0;
|
---|
| 637 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<NCz_;l++) {
|
---|
| 638 | if(n==nmod) {
|
---|
| 639 | xnt[0]=sqrt(i*i+j*j); xnt[1]=l; xnt[2]=T_(l,j,i).real(); xnt[3]=T_(l,j,i).imag();
|
---|
| 640 | nt.Fill(xnt);
|
---|
| 641 | n=0;
|
---|
| 642 | }
|
---|
| 643 | n++;
|
---|
| 644 | }
|
---|
| 645 | }
|
---|
| 646 |
|
---|
| 647 | pos.PutObject(nt,ntname);
|
---|
| 648 | }
|
---|
| 649 |
|
---|
| 650 | //-------------------------------------------------------
|
---|
[3115] | 651 | void GeneFluct3D::Print(void)
|
---|
| 652 | {
|
---|
[3141] | 653 | cout<<"GeneFluct3D(T_alloc="<<array_allocated_<<"):"<<endl;
|
---|
[3115] | 654 | cout<<"Space Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<Nz_<<" ("<<NTz_<<") size="
|
---|
| 655 | <<NRtot_<<endl;
|
---|
| 656 | cout<<" Resol: dx="<<Dx_<<" dy="<<Dy_<<" dz="<<Dz_<<" Mpc"
|
---|
| 657 | <<", dVol="<<dVol_<<", Vol="<<Vol_<<" Mpc^3"<<endl;
|
---|
| 658 | cout<<"Fourier Size : nx="<<Nx_<<" ny="<<Ny_<<" nz="<<NCz_<<endl;
|
---|
| 659 | cout<<" Resol: dkx="<<Dkx_<<" dky="<<Dky_<<" dkz="<<Dkz_<<" Mpc^-1"
|
---|
| 660 | <<", Dk3="<<Dk3_<<" Mpc^-3"<<endl;
|
---|
| 661 | cout<<" (2Pi/k: "<<2.*M_PI/Dkx_<<" "<<2.*M_PI/Dky_<<" "<<2.*M_PI/Dkz_<<" Mpc)"<<endl;
|
---|
| 662 | cout<<" Nyquist: kx="<<Knyqx_<<" ky="<<Knyqy_<<" kz="<<Knyqz_<<" Mpc^-1"
|
---|
| 663 | <<", Kmax="<<GetKmax()<<" Mpc^-1"<<endl;
|
---|
| 664 | cout<<" (2Pi/k: "<<2.*M_PI/Knyqx_<<" "<<2.*M_PI/Knyqy_<<" "<<2.*M_PI/Knyqz_<<" Mpc)"<<endl;
|
---|
[3271] | 665 | cout<<"Redshift "<<redsh_ref_<<" for z axe at k="<<kredsh_ref_<<endl;
|
---|
[3115] | 666 | }
|
---|
| 667 |
|
---|
| 668 | //-------------------------------------------------------
|
---|
[3768] | 669 | void GeneFluct3D::ComputeFourier0(PkSpectrumZ& pk_at_z)
|
---|
[3115] | 670 | // cf ComputeFourier() mais avec autre methode de realisation du spectre
|
---|
| 671 | // (attention on fait une fft pour realiser le spectre)
|
---|
| 672 | {
|
---|
[3768] | 673 | compute_pk_redsh_ref_ = pk_at_z.GetZ();
|
---|
[3115] | 674 | // --- realisation d'un tableau de tirage gaussiens
|
---|
[3768] | 675 | if(lp_>0) cout<<"--- ComputeFourier0 at z="<<compute_pk_redsh_ref_
|
---|
| 676 | <<": before gaussian filling ---"<<endl;
|
---|
[3115] | 677 | // On tient compte du pb de normalisation de FFTW3
|
---|
| 678 | double sntot = sqrt((double)NRtot_);
|
---|
[3129] | 679 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 680 | int_8 ip = IndexR(i,j,l);
|
---|
| 681 | data_[ip] = NorRand()/sntot;
|
---|
[3115] | 682 | }
|
---|
| 683 |
|
---|
| 684 | // --- realisation d'un tableau de tirage gaussiens
|
---|
[3155] | 685 | if(lp_>0) cout<<"...before fft real ---"<<endl;
|
---|
[3518] | 686 | GEN3D_FFTW_EXECUTE(pf_);
|
---|
[3115] | 687 |
|
---|
| 688 | // --- On remplit avec une realisation
|
---|
[3157] | 689 | if(lp_>0) cout<<"...before Fourier realization filling"<<endl;
|
---|
[3518] | 690 | T_(0,0,0) = complex<GEN3D_TYPE>(0.); // on coupe le continue et on l'initialise
|
---|
[3768] | 691 | long lmod = Nx_/20; if(lmod<1) lmod=1;
|
---|
[3129] | 692 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 693 | double kx = Kx(i); kx *= kx;
|
---|
| 694 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<"\t nx-i="<<Nx_-i<<"\t kx="<<kx<<"\t pk="<<pk_at_z(kx)<<endl;
|
---|
[3129] | 695 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 696 | double ky = Ky(j); ky *= ky;
|
---|
[3129] | 697 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 698 | double kz = Kz(l); kz *= kz;
|
---|
[3115] | 699 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
| 700 | double k = sqrt(kx+ky+kz);
|
---|
| 701 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
[3141] | 702 | double pk = pk_at_z(k)/Vol_;
|
---|
[3115] | 703 | // ici pas de "/2" a cause de la remarque ci-dessus
|
---|
| 704 | T_(l,j,i) *= sqrt(pk);
|
---|
| 705 | }
|
---|
| 706 | }
|
---|
| 707 | }
|
---|
| 708 |
|
---|
[3524] | 709 | array_type = 2;
|
---|
| 710 |
|
---|
[3155] | 711 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
[3115] | 712 | double p = compute_power_carte();
|
---|
[3155] | 713 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
[3115] | 714 |
|
---|
| 715 | }
|
---|
| 716 |
|
---|
| 717 | //-------------------------------------------------------
|
---|
[3768] | 718 | void GeneFluct3D::ComputeFourier(PkSpectrumZ& pk_at_z)
|
---|
[3141] | 719 | // Calcule une realisation du spectre "pk_at_z"
|
---|
[3115] | 720 | // Attention: dans TArray le premier indice varie le + vite
|
---|
| 721 | // Explication normalisation: see Coles & Lucchin, Cosmology, p264-265
|
---|
| 722 | // FFTW3: on note N=Nx*Ny*Nz
|
---|
| 723 | // f --(FFT)--> F = TF(f) --(FFT^-1)--> fb = TF^-1(F) = TF^-1(TF(f))
|
---|
| 724 | // sum(f(x_i)^2) = S
|
---|
| 725 | // sum(F(nu_i)^2) = S*N
|
---|
| 726 | // sum(fb(x_i)^2) = S*N^2
|
---|
| 727 | {
|
---|
| 728 | // --- RaZ du tableau
|
---|
[3518] | 729 | T_ = complex<GEN3D_TYPE>(0.);
|
---|
[3768] | 730 | compute_pk_redsh_ref_ = pk_at_z.GetZ();
|
---|
[3115] | 731 |
|
---|
| 732 | // --- On remplit avec une realisation
|
---|
[3768] | 733 | if(lp_>0) cout<<"--- ComputeFourier at z="<<compute_pk_redsh_ref_<<" ---"<<endl;
|
---|
| 734 | long lmod = Nx_/20; if(lmod<1) lmod=1;
|
---|
[3129] | 735 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 736 | double kx = Kx(i); kx *= kx;
|
---|
| 737 | if(lp_>0 && i%lmod==0) cout<<"i="<<i<<"\t nx-i="<<Nx_-i<<"\t kx="<<kx<<"\t pk="<<pk_at_z(kx)<<endl;
|
---|
[3129] | 738 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 739 | double ky = Ky(j); ky *= ky;
|
---|
[3129] | 740 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 741 | double kz = Kz(l); kz *= kz;
|
---|
[3115] | 742 | if(i==0 && j==0 && l==0) continue; // Suppression du continu
|
---|
| 743 | double k = sqrt(kx+ky+kz);
|
---|
| 744 | // cf normalisation: Peacock, Cosmology, formule 16.38 p504
|
---|
[3141] | 745 | double pk = pk_at_z(k)/Vol_;
|
---|
[3115] | 746 | // Explication de la division par 2: voir perandom.cc
|
---|
| 747 | // ou egalement Coles & Lucchin, Cosmology formula 13.7.2 p279
|
---|
[3617] | 748 | T_(l,j,i) = ComplexGaussianRand(sqrt(pk/2.));
|
---|
[3115] | 749 | }
|
---|
| 750 | }
|
---|
| 751 | }
|
---|
| 752 |
|
---|
[3524] | 753 | array_type = 2;
|
---|
[3115] | 754 | manage_coefficients(); // gros effet pour les spectres que l'on utilise !
|
---|
| 755 |
|
---|
[3155] | 756 | if(lp_>0) cout<<"...computing power"<<endl;
|
---|
[3115] | 757 | double p = compute_power_carte();
|
---|
[3155] | 758 | if(lp_>0) cout<<"Puissance dans la realisation: "<<p<<endl;
|
---|
[3115] | 759 |
|
---|
| 760 | }
|
---|
| 761 |
|
---|
[3129] | 762 | long GeneFluct3D::manage_coefficients(void)
|
---|
[3115] | 763 | // Take into account the real and complexe conjugate coefficients
|
---|
| 764 | // because we want a realization of a real data in real space
|
---|
[3768] | 765 | // On ecrit que: P(k_x,k_y,k_z) = P(-k_x,-k_y,-k_z)
|
---|
| 766 | // avec k_x = i, -k_x = N_x - i etc...
|
---|
[3115] | 767 | {
|
---|
[3155] | 768 | if(lp_>1) cout<<"...managing coefficients"<<endl;
|
---|
[3141] | 769 | check_array_alloc();
|
---|
[3115] | 770 |
|
---|
| 771 | // 1./ Le Continu et Nyquist sont reels
|
---|
[3129] | 772 | long nreal = 0;
|
---|
| 773 | for(long kk=0;kk<2;kk++) {
|
---|
| 774 | long k=0; // continu
|
---|
[3115] | 775 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
[3129] | 776 | for(long jj=0;jj<2;jj++) {
|
---|
| 777 | long j=0;
|
---|
[3115] | 778 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
[3129] | 779 | for(long ii=0;ii<2;ii++) {
|
---|
| 780 | long i=0;
|
---|
[3115] | 781 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
[3141] | 782 | int_8 ip = IndexC(i,j,k);
|
---|
| 783 | //cout<<"i="<<i<<" j="<<j<<" k="<<k<<" = ("<<fdata_[ip][0]<<","<<fdata_[ip][1]<<")"<<endl;
|
---|
| 784 | fdata_[ip][1] = 0.; fdata_[ip][0] *= M_SQRT2;
|
---|
[3115] | 785 | nreal++;
|
---|
| 786 | }
|
---|
| 787 | }
|
---|
| 788 | }
|
---|
[3155] | 789 | if(lp_>1) cout<<"Number of forced real number ="<<nreal<<endl;
|
---|
[3115] | 790 |
|
---|
| 791 | // 2./ Les elements complexe conjugues (tous dans le plan k=0,Nyquist)
|
---|
| 792 |
|
---|
| 793 | // a./ les lignes et colonnes du continu et de nyquist
|
---|
[3129] | 794 | long nconj1 = 0;
|
---|
| 795 | for(long kk=0;kk<2;kk++) {
|
---|
| 796 | long k=0; // continu
|
---|
[3115] | 797 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
[3129] | 798 | for(long jj=0;jj<2;jj++) { // selon j
|
---|
| 799 | long j=0;
|
---|
[3115] | 800 | if(jj==1) {if( Ny_%2!=0) continue; else j = Ny_/2;}
|
---|
[3129] | 801 | for(long i=1;i<(Nx_+1)/2;i++) {
|
---|
[3141] | 802 | int_8 ip = IndexC(i,j,k);
|
---|
| 803 | int_8 ip1 = IndexC(Nx_-i,j,k);
|
---|
| 804 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
[3115] | 805 | nconj1++;
|
---|
| 806 | }
|
---|
| 807 | }
|
---|
[3129] | 808 | for(long ii=0;ii<2;ii++) {
|
---|
| 809 | long i=0;
|
---|
[3115] | 810 | if(ii==1) {if( Nx_%2!=0) continue; else i = Nx_/2;}
|
---|
[3129] | 811 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
[3141] | 812 | int_8 ip = IndexC(i,j,k);
|
---|
| 813 | int_8 ip1 = IndexC(i,Ny_-j,k);
|
---|
| 814 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
[3115] | 815 | nconj1++;
|
---|
| 816 | }
|
---|
| 817 | }
|
---|
| 818 | }
|
---|
[3155] | 819 | if(lp_>1) cout<<"Number of forced conjugate on cont+nyq ="<<nconj1<<endl;
|
---|
[3115] | 820 |
|
---|
| 821 | // b./ les lignes et colonnes hors continu et de nyquist
|
---|
[3129] | 822 | long nconj2 = 0;
|
---|
| 823 | for(long kk=0;kk<2;kk++) {
|
---|
| 824 | long k=0; // continu
|
---|
[3115] | 825 | if(kk==1) {if(Nz_%2!=0) continue; else k = Nz_/2;} // Nyquist
|
---|
[3129] | 826 | for(long j=1;j<(Ny_+1)/2;j++) {
|
---|
[3115] | 827 | if(Ny_%2==0 && j==Ny_/2) continue; // on ne retraite pas nyquist en j
|
---|
[3129] | 828 | for(long i=1;i<Nx_;i++) {
|
---|
[3115] | 829 | if(Nx_%2==0 && i==Nx_/2) continue; // on ne retraite pas nyquist en i
|
---|
[3141] | 830 | int_8 ip = IndexC(i,j,k);
|
---|
| 831 | int_8 ip1 = IndexC(Nx_-i,Ny_-j,k);
|
---|
| 832 | fdata_[ip1][0] = fdata_[ip][0]; fdata_[ip1][1] = -fdata_[ip][1];
|
---|
[3115] | 833 | nconj2++;
|
---|
| 834 | }
|
---|
| 835 | }
|
---|
| 836 | }
|
---|
[3155] | 837 | if(lp_>1) cout<<"Number of forced conjugate hors cont+nyq ="<<nconj2<<endl;
|
---|
[3115] | 838 |
|
---|
[3768] | 839 | if(lp_>1) cout<<"Check: ddl= "<<NRtot_<<" =?= "<<2*(Nx_*Ny_*NCz_-nconj1-nconj2)-nreal<<endl;
|
---|
[3115] | 840 |
|
---|
| 841 | return nreal+nconj1+nconj2;
|
---|
| 842 | }
|
---|
| 843 |
|
---|
| 844 | double GeneFluct3D::compute_power_carte(void)
|
---|
| 845 | // Calcul la puissance de la realisation du spectre Pk
|
---|
| 846 | {
|
---|
[3141] | 847 | check_array_alloc();
|
---|
| 848 |
|
---|
[3115] | 849 | double s2 = 0.;
|
---|
[3129] | 850 | for(long l=0;l<NCz_;l++)
|
---|
| 851 | for(long j=0;j<Ny_;j++)
|
---|
| 852 | for(long i=0;i<Nx_;i++) s2 += MODULE2(T_(l,j,i));
|
---|
[3115] | 853 |
|
---|
| 854 | double s20 = 0.;
|
---|
[3129] | 855 | for(long j=0;j<Ny_;j++)
|
---|
| 856 | for(long i=0;i<Nx_;i++) s20 += MODULE2(T_(0,j,i));
|
---|
[3115] | 857 |
|
---|
| 858 | double s2n = 0.;
|
---|
| 859 | if(Nz_%2==0)
|
---|
[3129] | 860 | for(long j=0;j<Ny_;j++)
|
---|
| 861 | for(long i=0;i<Nx_;i++) s2n += MODULE2(T_(NCz_-1,j,i));
|
---|
[3115] | 862 |
|
---|
| 863 | return 2.*s2 -s20 -s2n;
|
---|
| 864 | }
|
---|
| 865 |
|
---|
| 866 | //-------------------------------------------------------------------
|
---|
| 867 | void GeneFluct3D::FilterByPixel(void)
|
---|
| 868 | // Filtrage par la fonction fenetre du pixel (parallelepipede)
|
---|
[3120] | 869 | // TF = 1/(dx*dy*dz)*Int[{-dx/2,dx/2},{-dy/2,dy/2},{-dz/2,dz/2}]
|
---|
[3115] | 870 | // e^(ik_x*x) e^(ik_y*y) e^(ik_z*z) dxdydz
|
---|
[3120] | 871 | // = 2/(k_x*dx) * sin(k_x*dx/2) * (idem y) * (idem z)
|
---|
| 872 | // Gestion divergence en 0: sin(y)/y = 1 - y^2/6*(1-y^2/20)
|
---|
| 873 | // avec y = k_x*dx/2
|
---|
[3115] | 874 | {
|
---|
[3155] | 875 | if(lp_>0) cout<<"--- FilterByPixel ---"<<endl;
|
---|
[3141] | 876 | check_array_alloc();
|
---|
| 877 |
|
---|
[3129] | 878 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 879 | double kx = Kx(i) *Dx_/2;
|
---|
[3330] | 880 | double pk_x = pixelfilter(kx);
|
---|
[3129] | 881 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 882 | double ky = Ky(j) *Dy_/2;
|
---|
[3330] | 883 | double pk_y = pixelfilter(ky);
|
---|
[3129] | 884 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 885 | double kz = Kz(l) *Dz_/2;
|
---|
[3141] | 886 | double pk_z = pixelfilter(kz);
|
---|
| 887 | T_(l,j,i) *= pk_x*pk_y*pk_z;
|
---|
[3115] | 888 | }
|
---|
| 889 | }
|
---|
| 890 | }
|
---|
| 891 |
|
---|
| 892 | }
|
---|
| 893 |
|
---|
| 894 | //-------------------------------------------------------------------
|
---|
[3768] | 895 | void GeneFluct3D::ToVelTrans(void)
|
---|
| 896 | {
|
---|
| 897 | cout<<"-------------------- TO BE VERIFIED ToVelTrans"<<endl;
|
---|
| 898 |
|
---|
| 899 | double zpk = compute_pk_redsh_ref_;
|
---|
| 900 | double dpsd = -cosmo_->H(zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
| 901 | if(lp_>0) cout<<"--- ToVelTrans --- at z="<<zpk<<", D'/D="<<dpsd
|
---|
| 902 | <<" (km/s)/Mpc (comp. for dz="<<good_dzinc_<<")"<<endl;
|
---|
| 903 | check_array_alloc();
|
---|
| 904 |
|
---|
| 905 | for(long i=0;i<Nx_;i++) {
|
---|
| 906 | double kx = Kx(i);
|
---|
| 907 | for(long j=0;j<Ny_;j++) {
|
---|
| 908 | double ky = Ky(j);
|
---|
| 909 | double kt2 = kx*kx + ky*ky;
|
---|
| 910 | for(long l=0;l<NCz_;l++) {
|
---|
| 911 | double kz = Kz(l);
|
---|
| 912 | double k2 = kt2 + kz*kz;
|
---|
| 913 | if(k2<=0.) continue;
|
---|
| 914 | T_(l,j,i) *= complex<double>(0.,dpsd*kz/k2);
|
---|
| 915 | }
|
---|
| 916 | }
|
---|
| 917 | }
|
---|
| 918 | }
|
---|
| 919 |
|
---|
| 920 | //-------------------------------------------------------------------
|
---|
[3331] | 921 | void GeneFluct3D::ApplyGrowthFactor(int type_evol)
|
---|
[3157] | 922 | // Apply Growth to real space
|
---|
| 923 | // Using the correspondance between redshift and los comoving distance
|
---|
| 924 | // describe in vector "zred_" "loscom_"
|
---|
[3516] | 925 | // type_evol = 1 : evolution avec la distance a l'observateur
|
---|
| 926 | // 2 : evolution avec la distance du plan Z
|
---|
[3331] | 927 | // (tous les pixels d'un plan Z sont mis au meme redshift z que celui du milieu)
|
---|
[3157] | 928 | {
|
---|
[3331] | 929 | if(lp_>0) cout<<"--- ApplyGrowthFactor: evol="<<type_evol<<endl;
|
---|
[3157] | 930 | check_array_alloc();
|
---|
| 931 |
|
---|
| 932 | if(growth_ == NULL) {
|
---|
[3572] | 933 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
|
---|
[3199] | 934 | cout<<bla<<endl; throw ParmError(bla);
|
---|
[3157] | 935 | }
|
---|
[3331] | 936 | if(type_evol<1 || type_evol>2) {
|
---|
[3572] | 937 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
|
---|
[3331] | 938 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 939 | }
|
---|
[3157] | 940 |
|
---|
[3199] | 941 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
[3157] | 942 |
|
---|
| 943 | for(long i=0;i<Nx_;i++) {
|
---|
[3331] | 944 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
[3157] | 945 | for(long j=0;j<Ny_;j++) {
|
---|
[3331] | 946 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
[3157] | 947 | for(long l=0;l<Nz_;l++) {
|
---|
[3331] | 948 | double dz = DZcom(l);
|
---|
| 949 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
| 950 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
[3768] | 951 | double z = interpinv(dz); // interpolation par morceau
|
---|
| 952 | double dzgr = (*growth_)(z);
|
---|
[3157] | 953 | int_8 ip = IndexR(i,j,l);
|
---|
| 954 | data_[ip] *= dzgr;
|
---|
| 955 | }
|
---|
| 956 | }
|
---|
| 957 | }
|
---|
| 958 |
|
---|
[3768] | 959 | }
|
---|
[3157] | 960 |
|
---|
[3768] | 961 | //-------------------------------------------------------------------
|
---|
| 962 | void GeneFluct3D::ApplyDerGrowthFactor(int type_evol)
|
---|
| 963 | // Apply Conformal derivative of Growth to real space for transverse velocity cube
|
---|
| 964 | {
|
---|
| 965 | cout<<"-------------------- TO BE VERIFIED ApplyDerGrowthFactor"<<endl;
|
---|
| 966 |
|
---|
| 967 | if(lp_>0) cout<<"--- ApplyDerGrowthFactor: evol="<<type_evol<<endl;
|
---|
| 968 | check_array_alloc();
|
---|
| 969 |
|
---|
| 970 | if(growth_ == NULL) {
|
---|
| 971 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: set GrowthFactor first";
|
---|
| 972 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 973 | }
|
---|
| 974 | if(type_evol<1 || type_evol>2) {
|
---|
| 975 | const char *bla = "GeneFluct3D::ApplyGrowthFactor_Error: bad type_evol value";
|
---|
| 976 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 977 | }
|
---|
| 978 |
|
---|
| 979 | double zpk = compute_pk_redsh_ref_;
|
---|
| 980 | double dpsd_orig = - cosmo_->H(zpk) * growth_->DsDz(zpk,good_dzinc_) / (*growth_)(zpk);
|
---|
| 981 |
|
---|
| 982 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
| 983 | InterpFunc interdpd(zredmin_dpsd_,zredmax_dpsd_,dpsdfrzred_);
|
---|
| 984 |
|
---|
| 985 | for(long i=0;i<Nx_;i++) {
|
---|
| 986 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
| 987 | for(long j=0;j<Ny_;j++) {
|
---|
| 988 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
| 989 | for(long l=0;l<Nz_;l++) {
|
---|
| 990 | double dz = DZcom(l);
|
---|
| 991 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
| 992 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
| 993 | double z = interpinv(dz); // interpolation par morceau
|
---|
| 994 | double dpsd = interdpd(z);
|
---|
| 995 | int_8 ip = IndexR(i,j,l);
|
---|
| 996 | data_[ip] *= dpsd / dpsd_orig;
|
---|
| 997 | }
|
---|
| 998 | }
|
---|
| 999 | }
|
---|
| 1000 |
|
---|
[3157] | 1001 | }
|
---|
| 1002 |
|
---|
| 1003 | //-------------------------------------------------------------------
|
---|
[3115] | 1004 | void GeneFluct3D::ComputeReal(void)
|
---|
| 1005 | // Calcule une realisation dans l'espace reel
|
---|
| 1006 | {
|
---|
[3155] | 1007 | if(lp_>0) cout<<"--- ComputeReal ---"<<endl;
|
---|
[3141] | 1008 | check_array_alloc();
|
---|
[3115] | 1009 |
|
---|
| 1010 | // On fait la FFT
|
---|
[3518] | 1011 | GEN3D_FFTW_EXECUTE(pb_);
|
---|
[3524] | 1012 | array_type = 1;
|
---|
[3115] | 1013 | }
|
---|
| 1014 |
|
---|
| 1015 | //-------------------------------------------------------------------
|
---|
| 1016 | void GeneFluct3D::ReComputeFourier(void)
|
---|
| 1017 | {
|
---|
[3155] | 1018 | if(lp_>0) cout<<"--- ReComputeFourier ---"<<endl;
|
---|
[3141] | 1019 | check_array_alloc();
|
---|
[3115] | 1020 |
|
---|
| 1021 | // On fait la FFT
|
---|
[3518] | 1022 | GEN3D_FFTW_EXECUTE(pf_);
|
---|
[3524] | 1023 | array_type = 2;
|
---|
| 1024 |
|
---|
[3115] | 1025 | // On corrige du pb de la normalisation de FFTW3
|
---|
[3594] | 1026 | complex<r_8> v((r_8)NRtot_,0.);
|
---|
[3129] | 1027 | for(long i=0;i<Nx_;i++)
|
---|
| 1028 | for(long j=0;j<Ny_;j++)
|
---|
[3594] | 1029 | for(long l=0;l<NCz_;l++) T_(l,j,i) /= v;
|
---|
[3115] | 1030 | }
|
---|
| 1031 |
|
---|
| 1032 | //-------------------------------------------------------------------
|
---|
[3141] | 1033 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr)
|
---|
| 1034 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
[3115] | 1035 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
| 1036 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
| 1037 | {
|
---|
[3155] | 1038 | if(lp_>0) cout<<"--- ComputeSpectrum ---"<<endl;
|
---|
[3141] | 1039 | check_array_alloc();
|
---|
[3115] | 1040 |
|
---|
[3141] | 1041 | if(herr.NBins()<0) return -1;
|
---|
| 1042 | herr.Zero();
|
---|
[3115] | 1043 |
|
---|
| 1044 | // Attention a l'ordre
|
---|
[3129] | 1045 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 1046 | double kx = Kx(i); kx *= kx;
|
---|
[3129] | 1047 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 1048 | double ky = Ky(j); ky *= ky;
|
---|
[3129] | 1049 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 1050 | double kz = Kz(l);
|
---|
[3330] | 1051 | double k = sqrt(kx+ky+kz*kz);
|
---|
[3115] | 1052 | double pk = MODULE2(T_(l,j,i));
|
---|
[3141] | 1053 | herr.Add(k,pk);
|
---|
[3115] | 1054 | }
|
---|
| 1055 | }
|
---|
| 1056 | }
|
---|
[3150] | 1057 | herr.ToVariance();
|
---|
[3115] | 1058 |
|
---|
| 1059 | // renormalize to directly compare to original spectrum
|
---|
| 1060 | double norm = Vol_;
|
---|
[3141] | 1061 | herr *= norm;
|
---|
[3115] | 1062 |
|
---|
| 1063 | return 0;
|
---|
| 1064 | }
|
---|
| 1065 |
|
---|
[3141] | 1066 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr)
|
---|
| 1067 | {
|
---|
[3155] | 1068 | if(lp_>0) cout<<"--- ComputeSpectrum2D ---"<<endl;
|
---|
[3141] | 1069 | check_array_alloc();
|
---|
| 1070 |
|
---|
| 1071 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
| 1072 | herr.Zero();
|
---|
| 1073 |
|
---|
| 1074 | // Attention a l'ordre
|
---|
| 1075 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 1076 | double kx = Kx(i); kx *= kx;
|
---|
[3141] | 1077 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 1078 | double ky = Ky(j); ky *= ky;
|
---|
[3141] | 1079 | double kt = sqrt(kx+ky);
|
---|
| 1080 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 1081 | double kz = Kz(l);
|
---|
[3141] | 1082 | double pk = MODULE2(T_(l,j,i));
|
---|
| 1083 | herr.Add(kt,kz,pk);
|
---|
| 1084 | }
|
---|
| 1085 | }
|
---|
| 1086 | }
|
---|
[3150] | 1087 | herr.ToVariance();
|
---|
[3141] | 1088 |
|
---|
| 1089 | // renormalize to directly compare to original spectrum
|
---|
| 1090 | double norm = Vol_;
|
---|
| 1091 | herr *= norm;
|
---|
| 1092 |
|
---|
| 1093 | return 0;
|
---|
| 1094 | }
|
---|
| 1095 |
|
---|
[3330] | 1096 | //-------------------------------------------------------------------
|
---|
| 1097 | int GeneFluct3D::ComputeSpectrum(HistoErr& herr,double sigma,bool pixcor)
|
---|
| 1098 | // Compute spectrum from "T" and fill HistoErr "herr"
|
---|
| 1099 | // AVEC la soustraction du niveau de bruit et la correction par filterpixel()
|
---|
| 1100 | // Si on ne fait pas ca, alors on obtient un spectre non-isotrope!
|
---|
| 1101 | //
|
---|
| 1102 | // T : dans le format standard de GeneFuct3D: T(nz,ny,nx)
|
---|
| 1103 | // cad T(kz,ky,kx) avec 0<kz<kz_nyq -ky_nyq<ky<ky_nyq -kx_nyq<kx<kx_nyq
|
---|
| 1104 | {
|
---|
| 1105 | if(lp_>0) cout<<"--- ComputeSpectrum: sigma="<<sigma<<endl;
|
---|
| 1106 | check_array_alloc();
|
---|
| 1107 |
|
---|
| 1108 | if(sigma<=0.) sigma = 0.;
|
---|
| 1109 | double sigma2 = sigma*sigma / (double)NRtot_;
|
---|
| 1110 |
|
---|
| 1111 | if(herr.NBins()<0) return -1;
|
---|
| 1112 | herr.Zero();
|
---|
| 1113 |
|
---|
| 1114 | TVector<r_8> vfz(NCz_);
|
---|
| 1115 | if(pixcor) // kz = l*Dkz_
|
---|
| 1116 | for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
|
---|
| 1117 |
|
---|
| 1118 | // Attention a l'ordre
|
---|
| 1119 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 1120 | double kx = Kx(i);
|
---|
[3330] | 1121 | double fx = (pixcor) ? pixelfilter(kx*Dx_/2): 1.;
|
---|
| 1122 | kx *= kx; fx *= fx;
|
---|
| 1123 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 1124 | double ky = Ky(j);
|
---|
[3330] | 1125 | double fy = (pixcor) ? pixelfilter(ky*Dy_/2): 1.;
|
---|
| 1126 | ky *= ky; fy *= fy;
|
---|
| 1127 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 1128 | double kz = Kz(l);
|
---|
[3330] | 1129 | double k = sqrt(kx+ky+kz*kz);
|
---|
| 1130 | double pk = MODULE2(T_(l,j,i)) - sigma2;
|
---|
| 1131 | double fz = (pixcor) ? vfz(l): 1.;
|
---|
| 1132 | double f = fx*fy*fz;
|
---|
| 1133 | if(f>0.) herr.Add(k,pk/f);
|
---|
| 1134 | }
|
---|
| 1135 | }
|
---|
| 1136 | }
|
---|
| 1137 | herr.ToVariance();
|
---|
[3351] | 1138 | for(int i=0;i<herr.NBins();i++) herr(i) += sigma2;
|
---|
[3330] | 1139 |
|
---|
| 1140 | // renormalize to directly compare to original spectrum
|
---|
| 1141 | double norm = Vol_;
|
---|
| 1142 | herr *= norm;
|
---|
| 1143 |
|
---|
| 1144 | return 0;
|
---|
| 1145 | }
|
---|
| 1146 |
|
---|
| 1147 | int GeneFluct3D::ComputeSpectrum2D(Histo2DErr& herr,double sigma,bool pixcor)
|
---|
| 1148 | // AVEC la soustraction du niveau de bruit et la correction par filterpixel()
|
---|
| 1149 | {
|
---|
| 1150 | if(lp_>0) cout<<"--- ComputeSpectrum2D: sigma="<<sigma<<endl;
|
---|
| 1151 | check_array_alloc();
|
---|
| 1152 |
|
---|
| 1153 | if(sigma<=0.) sigma = 0.;
|
---|
| 1154 | double sigma2 = sigma*sigma / (double)NRtot_;
|
---|
| 1155 |
|
---|
| 1156 | if(herr.NBinX()<0 || herr.NBinY()<0) return -1;
|
---|
| 1157 | herr.Zero();
|
---|
| 1158 |
|
---|
| 1159 | TVector<r_8> vfz(NCz_);
|
---|
| 1160 | if(pixcor) // kz = l*Dkz_
|
---|
| 1161 | for(long l=0;l<NCz_;l++) {vfz(l)=pixelfilter(l*Dkz_ *Dz_/2); vfz(l)*=vfz(l);}
|
---|
| 1162 |
|
---|
| 1163 | // Attention a l'ordre
|
---|
| 1164 | for(long i=0;i<Nx_;i++) {
|
---|
[3768] | 1165 | double kx = Kx(i);
|
---|
[3330] | 1166 | double fx = (pixcor) ? pixelfilter(kx*Dx_/2) : 1.;
|
---|
| 1167 | kx *= kx; fx *= fx;
|
---|
| 1168 | for(long j=0;j<Ny_;j++) {
|
---|
[3768] | 1169 | double ky = Ky(j);
|
---|
[3330] | 1170 | double fy = (pixcor) ? pixelfilter(ky*Dy_/2) : 1.;
|
---|
| 1171 | ky *= ky; fy *= fy;
|
---|
| 1172 | double kt = sqrt(kx+ky);
|
---|
| 1173 | for(long l=0;l<NCz_;l++) {
|
---|
[3768] | 1174 | double kz = Kz(l);
|
---|
[3330] | 1175 | double pk = MODULE2(T_(l,j,i)) - sigma2;
|
---|
| 1176 | double fz = (pixcor) ? vfz(l): 1.;
|
---|
| 1177 | double f = fx*fy*fz;
|
---|
| 1178 | if(f>0.) herr.Add(kt,kz,pk/f);
|
---|
| 1179 | }
|
---|
| 1180 | }
|
---|
| 1181 | }
|
---|
| 1182 | herr.ToVariance();
|
---|
[3351] | 1183 | for(int i=0;i<herr.NBinX();i++)
|
---|
| 1184 | for(int j=0;j<herr.NBinY();j++) herr(i,j) += sigma2;
|
---|
[3330] | 1185 |
|
---|
| 1186 | // renormalize to directly compare to original spectrum
|
---|
| 1187 | double norm = Vol_;
|
---|
| 1188 | herr *= norm;
|
---|
| 1189 |
|
---|
| 1190 | return 0;
|
---|
| 1191 | }
|
---|
| 1192 |
|
---|
[3115] | 1193 | //-------------------------------------------------------
|
---|
[3134] | 1194 | int_8 GeneFluct3D::VarianceFrReal(double R,double& var)
|
---|
[3115] | 1195 | // Recompute MASS variance in spherical top-hat (rayon=R)
|
---|
[3353] | 1196 | // Par definition: SigmaR^2 = <(M-<M>)^2>/<M>^2
|
---|
| 1197 | // ou M = masse dans sphere de rayon R
|
---|
[3354] | 1198 | // --- ATTENTION: la variance calculee a une tres grande dispersion
|
---|
| 1199 | // (surtout si le volume du cube est petit). Pour verifier
|
---|
| 1200 | // que le sigmaR calcule par cette methode est en accord avec
|
---|
| 1201 | // le sigmaR en input, il faut faire plusieurs simulations (~100)
|
---|
| 1202 | // et regarder la moyenne des sigmaR reconstruits
|
---|
[3115] | 1203 | {
|
---|
[3262] | 1204 | if(lp_>0) cout<<"--- VarianceFrReal R="<<R<<endl;
|
---|
[3141] | 1205 | check_array_alloc();
|
---|
| 1206 |
|
---|
[3353] | 1207 | long dnx = long(R/Dx_)+1; if(dnx<=0) dnx = 1;
|
---|
| 1208 | long dny = long(R/Dy_)+1; if(dny<=0) dny = 1;
|
---|
| 1209 | long dnz = long(R/Dz_)+1; if(dnz<=0) dnz = 1;
|
---|
[3155] | 1210 | if(lp_>0) cout<<"dnx="<<dnx<<" dny="<<dny<<" dnz="<<dnz<<endl;
|
---|
[3115] | 1211 |
|
---|
[3353] | 1212 | double sum=0., sum2=0., sn=0., r2 = R*R;
|
---|
| 1213 | int_8 nsum=0;
|
---|
[3115] | 1214 |
|
---|
[3353] | 1215 | for(long i=dnx;i<Nx_-dnx;i+=2*dnx) {
|
---|
| 1216 | for(long j=dny;j<Ny_-dny;j+=2*dny) {
|
---|
| 1217 | for(long l=dnz;l<Nz_-dnz;l+=2*dnz) {
|
---|
| 1218 | double m=0.; int_8 n=0;
|
---|
[3129] | 1219 | for(long ii=i-dnx;ii<=i+dnx;ii++) {
|
---|
[3115] | 1220 | double x = (ii-i)*Dx_; x *= x;
|
---|
[3129] | 1221 | for(long jj=j-dny;jj<=j+dny;jj++) {
|
---|
[3115] | 1222 | double y = (jj-j)*Dy_; y *= y;
|
---|
[3129] | 1223 | for(long ll=l-dnz;ll<=l+dnz;ll++) {
|
---|
[3115] | 1224 | double z = (ll-l)*Dz_; z *= z;
|
---|
| 1225 | if(x+y+z>r2) continue;
|
---|
[3141] | 1226 | int_8 ip = IndexR(ii,jj,ll);
|
---|
[3353] | 1227 | m += 1.+data_[ip]; // 1+drho/rho
|
---|
[3115] | 1228 | n++;
|
---|
| 1229 | }
|
---|
| 1230 | }
|
---|
| 1231 | }
|
---|
[3353] | 1232 | if(n>0) {sum += m; sum2 += m*m; nsum++; sn += n;}
|
---|
| 1233 | //cout<<i<<","<<j<<","<<l<<" n="<<n<<" m="<<m<<" sum="<<sum<<" sum2="<<sum2<<endl;
|
---|
[3115] | 1234 | }
|
---|
| 1235 | }
|
---|
| 1236 | }
|
---|
| 1237 |
|
---|
| 1238 | if(nsum<=1) {var=0.; return nsum;}
|
---|
| 1239 | sum /= nsum;
|
---|
| 1240 | sum2 = sum2/nsum - sum*sum;
|
---|
[3353] | 1241 | sn /= nsum;
|
---|
| 1242 | if(lp_>0) cout<<"...<n>="<<sn<<", nsum="<<nsum<<" <M>="<<sum<<" <(M-<M>)^2>="<<sum2<<endl;
|
---|
| 1243 | var = sum2/(sum*sum); // <dM^2>/<M>^2
|
---|
| 1244 | if(lp_>0) cout<<"...sigmaR^2 = <(M-<M>)^2>/<M>^2 = "<<var
|
---|
| 1245 | <<" -> sigmaR = "<<sqrt(var)<<endl;
|
---|
[3115] | 1246 |
|
---|
| 1247 | return nsum;
|
---|
| 1248 | }
|
---|
| 1249 |
|
---|
| 1250 | //-------------------------------------------------------
|
---|
[3134] | 1251 | int_8 GeneFluct3D::NumberOfBad(double vmin,double vmax)
|
---|
[3115] | 1252 | // number of pixels outside of ]vmin,vmax[ extremites exclues
|
---|
| 1253 | // -> vmin and vmax are considered as bad
|
---|
| 1254 | {
|
---|
[3141] | 1255 | check_array_alloc();
|
---|
[3115] | 1256 |
|
---|
[3134] | 1257 | int_8 nbad = 0;
|
---|
[3129] | 1258 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 1259 | int_8 ip = IndexR(i,j,l);
|
---|
| 1260 | double v = data_[ip];
|
---|
[3115] | 1261 | if(v<=vmin || v>=vmax) nbad++;
|
---|
| 1262 | }
|
---|
| 1263 |
|
---|
[3358] | 1264 | if(lp_>0) cout<<"--- NumberOfBad "<<nbad<<" px out of ]"<<vmin<<","<<vmax
|
---|
| 1265 | <<"[ i.e. frac="<<nbad/(double)NRtot_<<endl;
|
---|
[3115] | 1266 | return nbad;
|
---|
| 1267 | }
|
---|
| 1268 |
|
---|
[3320] | 1269 | int_8 GeneFluct3D::MinMax(double& xmin,double& xmax,double vmin,double vmax)
|
---|
| 1270 | // Calcul des valeurs xmin et xmax dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
| 1271 | {
|
---|
| 1272 | bool tstval = (vmax>vmin)? true: false;
|
---|
| 1273 | if(lp_>0) {
|
---|
| 1274 | cout<<"--- MinMax";
|
---|
| 1275 | if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
|
---|
| 1276 | cout<<endl;
|
---|
| 1277 | }
|
---|
| 1278 | check_array_alloc();
|
---|
| 1279 |
|
---|
| 1280 | int_8 n = 0;
|
---|
| 1281 | xmin = xmax = data_[0];
|
---|
| 1282 |
|
---|
| 1283 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 1284 | int_8 ip = IndexR(i,j,l);
|
---|
| 1285 | double x = data_[ip];
|
---|
| 1286 | if(tstval && (x<=vmin || x>=vmax)) continue;
|
---|
| 1287 | if(x<xmin) xmin = x;
|
---|
| 1288 | if(x>xmax) xmax = x;
|
---|
| 1289 | n++;
|
---|
| 1290 | }
|
---|
| 1291 |
|
---|
| 1292 | if(lp_>0) cout<<" n="<<n<<" min="<<xmin<<" max="<<xmax<<endl;
|
---|
| 1293 |
|
---|
| 1294 | return n;
|
---|
| 1295 | }
|
---|
| 1296 |
|
---|
[3261] | 1297 | int_8 GeneFluct3D::MeanSigma2(double& rm,double& rs2,double vmin,double vmax
|
---|
| 1298 | ,bool useout,double vout)
|
---|
| 1299 | // Calcul de mean,sigma2 dans le cube reel avec valeurs ]vmin,vmax[ extremites exclues
|
---|
| 1300 | // useout = false: ne pas utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
| 1301 | // true : utiliser les pixels hors limites pour calculer mean,sigma2
|
---|
| 1302 | // en remplacant leurs valeurs par "vout"
|
---|
[3115] | 1303 | {
|
---|
[3261] | 1304 | bool tstval = (vmax>vmin)? true: false;
|
---|
| 1305 | if(lp_>0) {
|
---|
[3262] | 1306 | cout<<"--- MeanSigma2";
|
---|
| 1307 | if(tstval) cout<<" range=]"<<vmin<<","<<vmax<<"[";
|
---|
[3261] | 1308 | if(useout) cout<<", useout="<<useout<<" vout="<<vout;
|
---|
| 1309 | cout<<endl;
|
---|
| 1310 | }
|
---|
[3141] | 1311 | check_array_alloc();
|
---|
[3115] | 1312 |
|
---|
[3134] | 1313 | int_8 n = 0;
|
---|
[3115] | 1314 | rm = rs2 = 0.;
|
---|
| 1315 |
|
---|
[3129] | 1316 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 1317 | int_8 ip = IndexR(i,j,l);
|
---|
| 1318 | double v = data_[ip];
|
---|
[3261] | 1319 | if(tstval) {
|
---|
| 1320 | if(v<=vmin || v>=vmax) {if(useout) v=vout; else continue;}
|
---|
| 1321 | }
|
---|
[3115] | 1322 | rm += v;
|
---|
| 1323 | rs2 += v*v;
|
---|
| 1324 | n++;
|
---|
| 1325 | }
|
---|
| 1326 |
|
---|
| 1327 | if(n>1) {
|
---|
| 1328 | rm /= (double)n;
|
---|
| 1329 | rs2 = rs2/(double)n - rm*rm;
|
---|
| 1330 | }
|
---|
| 1331 |
|
---|
[3261] | 1332 | if(lp_>0) cout<<" n="<<n<<" m="<<rm<<" s2="<<rs2<<" s="<<sqrt(fabs(rs2))<<endl;
|
---|
| 1333 |
|
---|
[3115] | 1334 | return n;
|
---|
| 1335 | }
|
---|
| 1336 |
|
---|
[3134] | 1337 | int_8 GeneFluct3D::SetToVal(double vmin, double vmax,double val0)
|
---|
[3115] | 1338 | // set to "val0" if out of range ]vmin,vmax[ extremites exclues
|
---|
[3261] | 1339 | // cad set to "val0" if in [vmin,vmax] -> vmin and vmax are set to val0
|
---|
[3115] | 1340 | {
|
---|
[3141] | 1341 | check_array_alloc();
|
---|
[3115] | 1342 |
|
---|
[3134] | 1343 | int_8 nbad = 0;
|
---|
[3129] | 1344 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 1345 | int_8 ip = IndexR(i,j,l);
|
---|
| 1346 | double v = data_[ip];
|
---|
| 1347 | if(v<=vmin || v>=vmax) {data_[ip] = val0; nbad++;}
|
---|
[3115] | 1348 | }
|
---|
| 1349 |
|
---|
[3262] | 1350 | if(lp_>0) cout<<"--- SetToVal "<<nbad<<" px set to="<<val0
|
---|
| 1351 | <<" because out of range=]"<<vmin<<","<<vmax<<"["<<endl;
|
---|
[3115] | 1352 | return nbad;
|
---|
| 1353 | }
|
---|
| 1354 |
|
---|
[3283] | 1355 | void GeneFluct3D::ScaleOffset(double scalecube,double offsetcube)
|
---|
| 1356 | // Replace "V" by "scalecube * ( V + offsetcube )"
|
---|
| 1357 | {
|
---|
[3284] | 1358 | if(lp_>0) cout<<"--- ScaleCube scale="<<scalecube<<" offset="<<offsetcube<<endl;
|
---|
[3283] | 1359 |
|
---|
| 1360 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 1361 | int_8 ip = IndexR(i,j,l);
|
---|
| 1362 | data_[ip] = scalecube * ( data_[ip] + offsetcube );
|
---|
| 1363 | }
|
---|
| 1364 |
|
---|
| 1365 | return;
|
---|
| 1366 | }
|
---|
| 1367 |
|
---|
[3115] | 1368 | //-------------------------------------------------------
|
---|
| 1369 | void GeneFluct3D::TurnFluct2Mass(void)
|
---|
| 1370 | // d_rho/rho -> Mass (add one!)
|
---|
| 1371 | {
|
---|
[3155] | 1372 | if(lp_>0) cout<<"--- TurnFluct2Mass ---"<<endl;
|
---|
[3141] | 1373 | check_array_alloc();
|
---|
| 1374 |
|
---|
[3115] | 1375 |
|
---|
[3129] | 1376 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 1377 | int_8 ip = IndexR(i,j,l);
|
---|
| 1378 | data_[ip] += 1.;
|
---|
[3115] | 1379 | }
|
---|
| 1380 | }
|
---|
| 1381 |
|
---|
[3358] | 1382 | double GeneFluct3D::TurnFluct2MeanNumber(double val_by_mpc3)
|
---|
[3365] | 1383 | // ATTENTION: la gestion des pixels<0 proposee ici induit une perte de variance
|
---|
| 1384 | // sur la carte, le spectre Pk reconstruit sera plus faible!
|
---|
| 1385 | // L'effet sera d'autant plus grand que le nombre de pixels<0 sera grand.
|
---|
[3329] | 1386 | {
|
---|
[3358] | 1387 | if(lp_>0) cout<<"--- TurnFluct2MeanNumber : "<<val_by_mpc3<<" quantity (gal or mass)/Mpc^3"<<endl;
|
---|
[3329] | 1388 |
|
---|
[3358] | 1389 | // First convert dRho/Rho into 1+dRho/Rho
|
---|
| 1390 | int_8 nball = 0; double sumall = 0., sumall2 = 0.;
|
---|
[3329] | 1391 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 1392 | int_8 ip = IndexR(i,j,l);
|
---|
[3358] | 1393 | data_[ip] += 1.;
|
---|
| 1394 | nball++; sumall += data_[ip]; sumall2 += data_[ip]*data_[ip];
|
---|
[3329] | 1395 | }
|
---|
[3358] | 1396 | if(nball>2) {
|
---|
| 1397 | sumall /= (double)nball;
|
---|
| 1398 | sumall2 = sumall2/(double)nball - sumall*sumall;
|
---|
| 1399 | if(lp_>0) cout<<"1+dRho/Rho: mean="<<sumall<<" variance="<<sumall2
|
---|
| 1400 | <<" -> "<<sqrt(fabs(sumall2))<<endl;
|
---|
[3329] | 1401 | }
|
---|
| 1402 |
|
---|
[3358] | 1403 | // Find contribution for positive pixels
|
---|
| 1404 | int_8 nbpos = 0; double sumpos = 0. , sumpos2 = 0.;
|
---|
| 1405 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 1406 | int_8 ip = IndexR(i,j,l);
|
---|
| 1407 | double v = data_[ip];
|
---|
| 1408 | if(data_[ip]>0.) {nbpos++; sumpos += v; sumpos2 += v*v;}
|
---|
| 1409 | }
|
---|
| 1410 | if(nbpos<1) {
|
---|
| 1411 | cout<<"TurnFluct2MeanNumber_Error: nbpos<1"<<endl;
|
---|
| 1412 | throw RangeCheckError("TurnFluct2MeanNumber_Error: nbpos<1");
|
---|
| 1413 | }
|
---|
| 1414 | sumpos2 = sumpos2/nball - sumpos*sumpos/(nball*nball);
|
---|
| 1415 | if(lp_>0)
|
---|
| 1416 | cout<<"1+dRho/Rho with v<0 set to zero: mean="<<sumpos/nball
|
---|
| 1417 | <<" variance="<<sumpos2<<" -> "<<sqrt(fabs(sumpos2))<<endl;
|
---|
| 1418 | cout<<"Sum of positive values: sumpos="<<sumpos
|
---|
| 1419 | <<" (n(v>0) = "<<nbpos<<" frac(v>0)="<<nbpos/(double)NRtot_<<")"<<endl;
|
---|
[3329] | 1420 |
|
---|
[3358] | 1421 | // - Mettre exactement val_by_mpc3*Vol galaxies (ou Msol) dans notre survey
|
---|
| 1422 | // - Uniquement dans les pixels de masse >0.
|
---|
| 1423 | // - Mise a zero des pixels <0
|
---|
| 1424 | double dn = val_by_mpc3 * Vol_ / sumpos;
|
---|
| 1425 | if(lp_>0) cout<<"...density move from "
|
---|
| 1426 | <<val_by_mpc3*dVol_<<" to "<<dn<<" / pixel"<<endl;
|
---|
| 1427 |
|
---|
[3329] | 1428 | double sum = 0.;
|
---|
| 1429 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 1430 | int_8 ip = IndexR(i,j,l);
|
---|
| 1431 | if(data_[ip]<=0.) data_[ip] = 0.;
|
---|
| 1432 | else {
|
---|
[3349] | 1433 | data_[ip] *= dn;
|
---|
| 1434 | sum += data_[ip];
|
---|
[3329] | 1435 | }
|
---|
| 1436 | }
|
---|
| 1437 |
|
---|
[3358] | 1438 | if(lp_>0) cout<<"...quantity put into survey "<<sum<<" / "<<val_by_mpc3*Vol_<<endl;
|
---|
[3329] | 1439 |
|
---|
| 1440 | return sum;
|
---|
| 1441 | }
|
---|
| 1442 |
|
---|
[3115] | 1443 | double GeneFluct3D::ApplyPoisson(void)
|
---|
| 1444 | // do NOT treate negative or nul mass -> let it as it is
|
---|
| 1445 | {
|
---|
[3155] | 1446 | if(lp_>0) cout<<"--- ApplyPoisson ---"<<endl;
|
---|
[3141] | 1447 | check_array_alloc();
|
---|
| 1448 |
|
---|
[3115] | 1449 | double sum = 0.;
|
---|
[3129] | 1450 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 1451 | int_8 ip = IndexR(i,j,l);
|
---|
| 1452 | double v = data_[ip];
|
---|
[3115] | 1453 | if(v>0.) {
|
---|
[3617] | 1454 | uint_8 dn = PoissonRand(v,10.);
|
---|
[3141] | 1455 | data_[ip] = (double)dn;
|
---|
[3115] | 1456 | sum += (double)dn;
|
---|
| 1457 | }
|
---|
| 1458 | }
|
---|
[3155] | 1459 | if(lp_>0) cout<<sum<<" galaxies put into survey"<<endl;
|
---|
[3115] | 1460 |
|
---|
| 1461 | return sum;
|
---|
| 1462 | }
|
---|
| 1463 |
|
---|
| 1464 | double GeneFluct3D::TurnNGal2Mass(FunRan& massdist,bool axeslog)
|
---|
| 1465 | // do NOT treate negative or nul mass -> let it as it is
|
---|
| 1466 | // INPUT:
|
---|
| 1467 | // massdist : distribution de masse (m*dn/dm)
|
---|
| 1468 | // axeslog = false : retourne la masse
|
---|
| 1469 | // = true : retourne le log10(mass)
|
---|
| 1470 | // RETURN la masse totale
|
---|
| 1471 | {
|
---|
[3155] | 1472 | if(lp_>0) cout<<"--- TurnNGal2Mass ---"<<endl;
|
---|
[3141] | 1473 | check_array_alloc();
|
---|
| 1474 |
|
---|
[3115] | 1475 | double sum = 0.;
|
---|
[3129] | 1476 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
[3141] | 1477 | int_8 ip = IndexR(i,j,l);
|
---|
| 1478 | double v = data_[ip];
|
---|
[3115] | 1479 | if(v>0.) {
|
---|
[3129] | 1480 | long ngal = long(v+0.1);
|
---|
[3141] | 1481 | data_[ip] = 0.;
|
---|
[3129] | 1482 | for(long i=0;i<ngal;i++) {
|
---|
[3115] | 1483 | double m = massdist.RandomInterp(); // massdist.Random();
|
---|
| 1484 | if(axeslog) m = pow(10.,m);
|
---|
[3141] | 1485 | data_[ip] += m;
|
---|
[3115] | 1486 | }
|
---|
[3141] | 1487 | sum += data_[ip];
|
---|
[3115] | 1488 | }
|
---|
| 1489 | }
|
---|
[3155] | 1490 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
[3115] | 1491 |
|
---|
| 1492 | return sum;
|
---|
| 1493 | }
|
---|
| 1494 |
|
---|
[3320] | 1495 | double GeneFluct3D::TurnNGal2MassQuick(SchechterMassDist& schmdist)
|
---|
| 1496 | // idem TurnNGal2Mass mais beaucoup plus rapide
|
---|
| 1497 | {
|
---|
| 1498 | if(lp_>0) cout<<"--- TurnNGal2MassQuick ---"<<endl;
|
---|
| 1499 | check_array_alloc();
|
---|
| 1500 |
|
---|
| 1501 | double sum = 0.;
|
---|
| 1502 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) for(long l=0;l<Nz_;l++) {
|
---|
| 1503 | int_8 ip = IndexR(i,j,l);
|
---|
| 1504 | double v = data_[ip];
|
---|
| 1505 | if(v>0.) {
|
---|
| 1506 | long ngal = long(v+0.1);
|
---|
| 1507 | data_[ip] = schmdist.TirMass(ngal);
|
---|
| 1508 | sum += data_[ip];
|
---|
| 1509 | }
|
---|
| 1510 | }
|
---|
| 1511 | if(lp_>0) cout<<sum<<" MSol HI mass put into survey"<<endl;
|
---|
| 1512 |
|
---|
| 1513 | return sum;
|
---|
| 1514 | }
|
---|
| 1515 |
|
---|
[3349] | 1516 | void GeneFluct3D::AddNoise2Real(double snoise,int type_evol)
|
---|
| 1517 | // add noise to every pixels (meme les <=0 !)
|
---|
| 1518 | // type_evol = 0 : pas d'evolution de la puissance du bruit
|
---|
| 1519 | // 1 : evolution de la puissance du bruit avec la distance a l'observateur
|
---|
| 1520 | // 2 : evolution de la puissance du bruit avec la distance du plan Z
|
---|
| 1521 | // (tous les plans Z sont mis au meme redshift z de leur milieu)
|
---|
| 1522 | {
|
---|
| 1523 | if(lp_>0) cout<<"--- AddNoise2Real: snoise = "<<snoise<<" evol="<<type_evol<<endl;
|
---|
| 1524 | check_array_alloc();
|
---|
| 1525 |
|
---|
| 1526 | if(type_evol<0) type_evol = 0;
|
---|
| 1527 | if(type_evol>2) {
|
---|
[3572] | 1528 | const char *bla = "GeneFluct3D::AddNoise2Real_Error: bad type_evol value";
|
---|
[3349] | 1529 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 1530 | }
|
---|
| 1531 |
|
---|
| 1532 | vector<double> correction;
|
---|
| 1533 | InterpFunc *intercor = NULL;
|
---|
| 1534 |
|
---|
| 1535 | if(type_evol>0) {
|
---|
| 1536 | // Sigma_Noise(en mass) :
|
---|
| 1537 | // Slim ~ 1/sqrt(DNu) * sqrt(nlobe) en W/m^2Hz
|
---|
| 1538 | // Flim ~ sqrt(DNu) * sqrt(nlobe) en W/m^2
|
---|
| 1539 | // Mlim ~ sqrt(DNu) * (Dlum)^2 * sqrt(nlobe) en Msol
|
---|
| 1540 | // nlobe ~ 1/Dtrcom^2
|
---|
| 1541 | // Mlim ~ sqrt(DNu) * (Dlum)^2 / Dtrcom
|
---|
[3662] | 1542 | if(cosmo_ == NULL || redsh_ref_<0. || loscom2zred_.size()<1) {
|
---|
[3572] | 1543 | const char *bla = "GeneFluct3D::AddNoise2Real_Error: set Observator and Cosmology first";
|
---|
[3349] | 1544 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 1545 | }
|
---|
| 1546 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
| 1547 | long nsz = loscom2zred_.size(), nszmod=((nsz>10)? nsz/10: 1);
|
---|
| 1548 | for(long i=0;i<nsz;i++) {
|
---|
| 1549 | double d = interpinv.X(i);
|
---|
| 1550 | double zred = interpinv(d);
|
---|
| 1551 | double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
|
---|
| 1552 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
| 1553 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
| 1554 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
| 1555 | double corr = sqrt(dnu/dnu_ref_) * pow(dlum/dlum_ref_,2.) * dtrc_ref_/dtrc;
|
---|
| 1556 | if(lp_>0 && (i==0 || i==nsz-1 || i%nszmod==0))
|
---|
| 1557 | cout<<"i="<<i<<" d="<<d<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
|
---|
| 1558 | <<" dtrc="<<dtrc<<" dlum="<<dlum<<" -> cor="<<corr<<endl;
|
---|
| 1559 | correction.push_back(corr);
|
---|
| 1560 | }
|
---|
| 1561 | intercor = new InterpFunc(loscom2zred_min_,loscom2zred_max_,correction);
|
---|
| 1562 | }
|
---|
| 1563 |
|
---|
| 1564 | double corrlim[2] = {1.,1.};
|
---|
| 1565 | for(long i=0;i<Nx_;i++) {
|
---|
| 1566 | double dx2 = DXcom(i); dx2 *= dx2;
|
---|
| 1567 | for(long j=0;j<Ny_;j++) {
|
---|
| 1568 | double dy2 = DYcom(j); dy2 *= dy2;
|
---|
| 1569 | for(long l=0;l<Nz_;l++) {
|
---|
| 1570 | double corr = 1.;
|
---|
| 1571 | if(type_evol>0) {
|
---|
| 1572 | double dz = DZcom(l);
|
---|
| 1573 | if(type_evol==1) dz = sqrt(dx2+dy2+dz*dz);
|
---|
| 1574 | else dz = fabs(dz); // tous les plans Z au meme redshift
|
---|
| 1575 | corr = (*intercor)(dz);
|
---|
| 1576 | if(corr<corrlim[0]) corrlim[0]=corr; else if(corr>corrlim[1]) corrlim[1]=corr;
|
---|
| 1577 | }
|
---|
| 1578 | int_8 ip = IndexR(i,j,l);
|
---|
| 1579 | data_[ip] += snoise*corr*NorRand();
|
---|
| 1580 | }
|
---|
| 1581 | }
|
---|
| 1582 | }
|
---|
| 1583 | if(type_evol>0)
|
---|
| 1584 | cout<<"correction factor range: ["<<corrlim[0]<<","<<corrlim[1]<<"]"<<endl;
|
---|
| 1585 |
|
---|
| 1586 | if(intercor!=NULL) delete intercor;
|
---|
| 1587 | }
|
---|
| 1588 |
|
---|
| 1589 | } // Fin namespace SOPHYA
|
---|
| 1590 |
|
---|
| 1591 |
|
---|
| 1592 |
|
---|
| 1593 |
|
---|
| 1594 | /*********************************************************************
|
---|
[3199] | 1595 | void GeneFluct3D::AddAGN(double lfjy,double lsigma,double powlaw)
|
---|
[3196] | 1596 | // Add AGN flux into simulation:
|
---|
| 1597 | // --- Procedure:
|
---|
| 1598 | // 1. lancer "cmvdefsurv" avec les parametres du survey
|
---|
[3199] | 1599 | // (au redshift de reference du survey)
|
---|
[3196] | 1600 | // et recuperer l'angle solide "angsol sr" du pixel elementaire
|
---|
| 1601 | // au centre du cube.
|
---|
| 1602 | // 2. lancer "cmvtstagn" pour cet angle solide -> cmvtstagn.ppf
|
---|
| 1603 | // 3. regarder l'histo "hlfang" et en deduire un equivalent gaussienne
|
---|
| 1604 | // cad une moyenne <log10(S)> et un sigma "sig"
|
---|
[3199] | 1605 | // Attention: la distribution n'est pas gaussienne les "mean,sigma"
|
---|
| 1606 | // de l'histo ne sont pas vraiment ce que l'on veut
|
---|
[3196] | 1607 | // --- Limitations actuelle du code:
|
---|
[3271] | 1608 | // . les AGN sont supposes evoluer avec la meme loi de puissance pour tout theta,phi
|
---|
[3199] | 1609 | // . le flux des AGN est mis dans une colonne Oz (indice k) et pas sur la ligne de visee
|
---|
| 1610 | // . la distribution est approximee a une gaussienne
|
---|
| 1611 | // ... C'est une approximation pour un observateur loin du centre du cube
|
---|
| 1612 | // et pour un cube peu epais / distance observateur
|
---|
[3196] | 1613 | // --- Parametres de la routine:
|
---|
[3271] | 1614 | // llfy : c'est le <log10(S)> du flux depose par les AGN
|
---|
| 1615 | // dans l'angle solide du pixel elementaire de reference du cube
|
---|
| 1616 | // lsigma : c'est le sigma de la distribution des log10(S)
|
---|
| 1617 | // powlaw : c'est la pente de la distribution cad que le flux "lmsol"
|
---|
[3199] | 1618 | // et considere comme le flux a 1.4GHz et qu'on suppose une loi
|
---|
| 1619 | // F(nu) = (1.4GHz/nu)^powlaw * F(1.4GHz)
|
---|
[3196] | 1620 | // - Comme on est en echelle log10():
|
---|
| 1621 | // on tire log10(Msol) + X
|
---|
| 1622 | // ou X est une realisation sur une gaussienne de variance "sig^2"
|
---|
| 1623 | // La masse realisee est donc: Msol*10^X
|
---|
| 1624 | // - Pas de probleme de pixel negatif car on a une multiplication!
|
---|
| 1625 | {
|
---|
[3199] | 1626 | if(lp_>0) cout<<"--- AddAGN: <log10(S Jy)> = "<<lfjy<<" , sigma = "<<lsigma<<endl;
|
---|
[3196] | 1627 | check_array_alloc();
|
---|
| 1628 |
|
---|
[3271] | 1629 | if(cosmo_ == NULL || redsh_ref_<0.| loscom2zred_.size()<1) {
|
---|
[3199] | 1630 | char *bla = "GeneFluct3D::AddAGN_Error: set Observator and Cosmology first";
|
---|
| 1631 | cout<<bla<<endl; throw ParmError(bla);
|
---|
| 1632 | }
|
---|
[3196] | 1633 |
|
---|
[3271] | 1634 | // Le flux des AGN en Jy et en mass solaire
|
---|
| 1635 | double fagnref = pow(10.,lfjy)*(dnu_ref_*1.e9); // Jy.Hz = W/m^2
|
---|
| 1636 | double magnref = FluxHI2Msol(fagnref*Jansky2Watt_cst,dlum_ref_); // Msol
|
---|
| 1637 | if(lp_>0)
|
---|
| 1638 | cout<<"Au pixel de ref: fagnref="<<fagnref
|
---|
| 1639 | <<" Jy.Hz (a 1.4GHz), magnref="<<magnref<<" Msol"<<endl;
|
---|
[3196] | 1640 |
|
---|
[3199] | 1641 | if(powlaw!=0.) {
|
---|
[3271] | 1642 | // 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)
|
---|
| 1643 | magnref *= pow(1/(1.+redsh_ref_),powlaw);
|
---|
[3199] | 1644 | if(lp_>0) cout<<" powlaw="<<powlaw<<" -> change magnref to "<<magnref<<" Msol"<<endl;
|
---|
| 1645 | }
|
---|
| 1646 |
|
---|
| 1647 | // Les infos en fonction de l'indice "l" selon Oz
|
---|
| 1648 | vector<double> correction;
|
---|
| 1649 | InterpFunc interpinv(loscom2zred_min_,loscom2zred_max_,loscom2zred_);
|
---|
[3271] | 1650 | long nzmod = ((Nz_>10)?Nz_/10:1);
|
---|
[3199] | 1651 | for(long l=0;l<Nz_;l++) {
|
---|
[3271] | 1652 | double z = fabs(DZcom(l));
|
---|
[3199] | 1653 | double zred = interpinv(z);
|
---|
[3271] | 1654 | double dtrc = cosmo_->Dtrcom(zred); // pour variation angle solide
|
---|
[3199] | 1655 | double dlum = cosmo_->Dlum(zred); // pour variation conversion mass HI
|
---|
| 1656 | double dred = Dz_/(cosmo_->Dhubble()/cosmo_->E(zred));
|
---|
| 1657 | double dnu = Fr_HyperFin_Par *dred/pow(1.+zred,2.); // pour variation dNu
|
---|
[3271] | 1658 | // on a: Mass ~ DNu * Dlum^2 / Dtrcom^2
|
---|
| 1659 | double corr = dnu/dnu_ref_*pow(dtrc_ref_/dtrc*dlum/dlum_ref_,2.);
|
---|
| 1660 | // 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)
|
---|
| 1661 | if(powlaw!=0.) corr *= pow((1.+redsh_ref_)/(1.+zred),powlaw);
|
---|
[3199] | 1662 | correction.push_back(corr);
|
---|
[3271] | 1663 | if(lp_>0 && (l==0 || l==Nz_-1 || l%nzmod==0)) {
|
---|
| 1664 | cout<<"l="<<l<<" z="<<z<<" red="<<zred<<" dred="<<dred<<" dnu="<<dnu
|
---|
| 1665 | <<" dtrc="<<dtrc<<" dlum="<<dlum
|
---|
| 1666 | <<" -> cor="<<corr<<endl;
|
---|
[3199] | 1667 | }
|
---|
| 1668 | }
|
---|
| 1669 |
|
---|
| 1670 | double sum=0., sum2=0., nsum=0.;
|
---|
| 1671 | for(long i=0;i<Nx_;i++) for(long j=0;j<Ny_;j++) {
|
---|
| 1672 | double a = lsigma*NorRand();
|
---|
| 1673 | a = magnref*pow(10.,a);
|
---|
| 1674 | // On met le meme tirage le long de Oz (indice k)
|
---|
| 1675 | for(long l=0;l<Nz_;l++) {
|
---|
| 1676 | int_8 ip = IndexR(i,j,l);
|
---|
| 1677 | data_[ip] += a*correction[l];
|
---|
| 1678 | }
|
---|
| 1679 | sum += a; sum2 += a*a; nsum += 1.;
|
---|
| 1680 | }
|
---|
| 1681 |
|
---|
| 1682 | if(lp_>0 && nsum>1.) {
|
---|
[3196] | 1683 | sum /= nsum;
|
---|
| 1684 | sum2 = sum2/nsum - sum*sum;
|
---|
| 1685 | cout<<"...Mean mass="<<sum<<" Msol , s^2="<<sum2<<" s="<<sqrt(fabs(sum2))<<endl;
|
---|
| 1686 | }
|
---|
| 1687 |
|
---|
| 1688 | }
|
---|
[3349] | 1689 | *********************************************************************/
|
---|