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