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