[3756] | 1 | // Classes to compute 3D power spectrum
|
---|
| 2 | // R. Ansari - Nov 2008, May 2010
|
---|
| 3 |
|
---|
| 4 | #include "specpk.h"
|
---|
| 5 | #include "randr48.h"
|
---|
| 6 |
|
---|
| 7 | //------------------------------------
|
---|
| 8 | // Class SpectralShape
|
---|
| 9 | // -----------------------------------
|
---|
| 10 |
|
---|
| 11 | double Pnu1(double nu)
|
---|
| 12 | {
|
---|
| 13 | return ( sqrt(sqrt(nu)) / ((nu+1.0)/0.2) *
|
---|
| 14 | (1+0.2*cos(2*M_PI*(nu-2.)*0.15)*exp(-nu/50.)) );
|
---|
| 15 | }
|
---|
| 16 |
|
---|
| 17 | double Pnu2(double nu)
|
---|
| 18 | {
|
---|
| 19 | if (nu < 1.e-9) return 0.;
|
---|
| 20 | return ((1.-exp(-nu/0.5))/nu*(1+0.25*cos(2*M_PI*nu*0.1)*exp(-nu/20.)) );
|
---|
| 21 | }
|
---|
| 22 |
|
---|
| 23 |
|
---|
| 24 | double Pnu3(double nu)
|
---|
| 25 | {
|
---|
| 26 | return ( log(nu/100.+1)*(1+sin(2*M_PI*nu/300))*exp(-nu/4000) );
|
---|
| 27 | }
|
---|
| 28 |
|
---|
| 29 |
|
---|
| 30 | double Pnu4(double nu)
|
---|
| 31 | {
|
---|
| 32 | double x = (nu-0.5)/0.05;
|
---|
| 33 | double rc = 2*exp(-x*x);
|
---|
| 34 | x = (nu-3.1)/0.27;
|
---|
| 35 | rc += exp(-x*x);
|
---|
| 36 | x = (nu-7.6)/1.4;
|
---|
| 37 | rc += 0.5*exp(-x*x);
|
---|
| 38 | return ( rc+2.*exp(-x*x) );
|
---|
| 39 | }
|
---|
| 40 |
|
---|
| 41 | //--------------------------------------------------
|
---|
| 42 | // -- SpectralShape class : test P(k) class
|
---|
| 43 | //--------------------------------------------------
|
---|
| 44 | // Constructor
|
---|
| 45 | SpectralShape::SpectralShape(int typ)
|
---|
| 46 | {
|
---|
| 47 | typ_=typ;
|
---|
[3825] | 48 | SetRenormFac();
|
---|
[3756] | 49 | }
|
---|
| 50 |
|
---|
| 51 | // Return the spectral power for a given wave number wk
|
---|
| 52 | double SpectralShape::operator() (double wk)
|
---|
| 53 | {
|
---|
| 54 | wk/=DeuxPI;
|
---|
[3825] | 55 | double retv=1.;
|
---|
[3756] | 56 | switch (typ_)
|
---|
| 57 | {
|
---|
| 58 | case 1:
|
---|
[3825] | 59 | retv=Pnu1(wk);
|
---|
[3756] | 60 | break;
|
---|
| 61 | case 2:
|
---|
[3825] | 62 | retv=Pnu2(wk);
|
---|
[3756] | 63 | break;
|
---|
| 64 | case 3:
|
---|
[3825] | 65 | retv=Pnu3(wk);
|
---|
[3756] | 66 | break;
|
---|
| 67 | case 4:
|
---|
[3825] | 68 | retv=Pnu4(wk);
|
---|
[3756] | 69 | break;
|
---|
| 70 | default :
|
---|
| 71 | {
|
---|
| 72 | // global shape
|
---|
| 73 | double csp = pow( (2*sin(sqrt(sqrt(wk/7.)))),2.);
|
---|
| 74 | if (csp < 0.) return 0.;
|
---|
| 75 |
|
---|
| 76 | // Adding some pics
|
---|
| 77 | double picpos[5] = {75.,150.,225.,300.,375.,};
|
---|
| 78 |
|
---|
| 79 | for(int k=0; k<5; k++) {
|
---|
| 80 | double x0 = picpos[k];
|
---|
| 81 | if ( (wk > x0-25.) && (wk < x0+25.) ) {
|
---|
| 82 | double x = (wk-x0);
|
---|
| 83 | csp *= (1.+0.5*exp(-(x*x)/(2.*5*5)));
|
---|
| 84 | break;
|
---|
| 85 | }
|
---|
| 86 | }
|
---|
[3825] | 87 | retv=csp;
|
---|
[3756] | 88 | }
|
---|
| 89 | break;
|
---|
| 90 | }
|
---|
[3825] | 91 | return retv*renorm_fac;
|
---|
[3756] | 92 | }
|
---|
| 93 | // Return a vector representing the power spectrum (for checking)
|
---|
| 94 | Histo SpectralShape::GetPk(int n)
|
---|
| 95 | {
|
---|
| 96 | if (n<16) n = 256;
|
---|
| 97 | Histo h(0.,1024.*DeuxPI,n);
|
---|
| 98 | for(int k=0; k<h.NBins(); k++) h(k) = Value((k+0.5)*h.BinWidth());
|
---|
| 99 | return h;
|
---|
| 100 | }
|
---|
| 101 |
|
---|
[3825] | 102 | double SpectralShape::Sommek2Pk(double kmax, int n)
|
---|
| 103 | {
|
---|
| 104 | double dk=kmax/(double)n;
|
---|
| 105 | double s=0.;
|
---|
| 106 | for(int i=1; i<=n; i++) {
|
---|
| 107 | double ck=(double)i*dk;
|
---|
| 108 | s += Value(ck)*ck*ck;
|
---|
| 109 | }
|
---|
| 110 | return s*dk*4.*M_PI;
|
---|
| 111 | }
|
---|
[3756] | 112 | //--------------------------------------------------
|
---|
| 113 | // -- Four2DResponse class : test P(k) class
|
---|
| 114 |
|
---|
| 115 | //---------------------------------------------------------------
|
---|
| 116 | // -- Four3DPk class : 3D fourier amplitudes and power spectrum
|
---|
| 117 | //---------------------------------------------------------------
|
---|
| 118 | // Constructeur avec Tableau des coeff. de Fourier en argument
|
---|
| 119 | Four3DPk::Four3DPk(TArray< complex<TF> > & fourcoedd, RandomGeneratorInterface& rg)
|
---|
| 120 | : rg_(rg), fourAmp(fourcoedd)
|
---|
| 121 | {
|
---|
| 122 | SetPrtLevel();
|
---|
| 123 | SetCellSize();
|
---|
| 124 | }
|
---|
| 125 | // Constructor
|
---|
| 126 | Four3DPk::Four3DPk(RandomGeneratorInterface& rg, sa_size_t szx, sa_size_t szy, sa_size_t szz)
|
---|
| 127 | : rg_(rg), fourAmp(szx, szy, szz)
|
---|
| 128 | {
|
---|
| 129 | SetPrtLevel();
|
---|
| 130 | SetCellSize();
|
---|
| 131 | }
|
---|
| 132 |
|
---|
| 133 |
|
---|
| 134 | // Generate mass field Fourier Coefficient
|
---|
| 135 | void Four3DPk::ComputeFourierAmp(SpectralShape& pk)
|
---|
| 136 | {
|
---|
| 137 | // We generate a random gaussian real field
|
---|
| 138 | // fourAmp represent 3-D fourier transform of a real input array.
|
---|
| 139 | // The second half of the array along Y and Z contain negative frequencies
|
---|
| 140 | // double fnorm = 1./sqrt(2.*fourAmp.Size());
|
---|
| 141 | double fnorm = 1.;
|
---|
| 142 | double kxx, kyy, kzz;
|
---|
| 143 | // sa_size_t is large integer type
|
---|
| 144 | for(sa_size_t kz=0; kz<fourAmp.SizeZ(); kz++) {
|
---|
| 145 | kzz = (kz>fourAmp.SizeZ()/2) ? (double)(fourAmp.SizeZ()-kz)*dkz_ : (double)kz*dkz_;
|
---|
| 146 | for(sa_size_t ky=0; ky<fourAmp.SizeY(); ky++) {
|
---|
| 147 | kyy = (ky>fourAmp.SizeY()/2) ? (double)(fourAmp.SizeY()-ky)*dky_ : (double)ky*dky_;
|
---|
| 148 | for(sa_size_t kx=0; kx<fourAmp.SizeX(); kx++) {
|
---|
| 149 | double kxx=(double)kx*dkx_;
|
---|
| 150 | double wk = sqrt(kxx*kxx+kyy*kyy+kzz*kzz);
|
---|
| 151 | double amp = sqrt(pk(wk)*fnorm/2.);
|
---|
| 152 | fourAmp(kx, ky, kz) = complex<TF>(rg_.Gaussian(amp), rg_.Gaussian(amp)); // renormalize fourier coeff usin
|
---|
| 153 | }
|
---|
| 154 | }
|
---|
| 155 | }
|
---|
| 156 | if (prtlev_>0)
|
---|
| 157 | cout << " Four3DPk::ComputeFourierAmp() done ..." << endl;
|
---|
| 158 | }
|
---|
| 159 |
|
---|
| 160 | // Generate mass field Fourier Coefficient
|
---|
| 161 | void Four3DPk::ComputeNoiseFourierAmp(Four2DResponse& resp, bool crmask)
|
---|
| 162 | {
|
---|
| 163 | TMatrix<r_4> mask(fourAmp.SizeY(), fourAmp.SizeX());
|
---|
| 164 | // fourAmp represent 3-D fourier transform of a real input array.
|
---|
| 165 | // The second half of the array along Y and Z contain negative frequencies
|
---|
[3787] | 166 | double kxx, kyy, kzz, rep, amp;
|
---|
[3756] | 167 | // sa_size_t is large integer type
|
---|
| 168 | for(sa_size_t kz=0; kz<fourAmp.SizeZ(); kz++) {
|
---|
[3769] | 169 | kzz = (kz>fourAmp.SizeZ()/2) ? -(double)(fourAmp.SizeZ()-kz)*dkz_ : (double)kz*dkz_;
|
---|
[3756] | 170 | for(sa_size_t ky=0; ky<fourAmp.SizeY(); ky++) {
|
---|
[3769] | 171 | kyy = (ky>fourAmp.SizeY()/2) ? -(double)(fourAmp.SizeY()-ky)*dky_ : (double)ky*dky_;
|
---|
[3756] | 172 | for(sa_size_t kx=0; kx<fourAmp.SizeX(); kx++) {
|
---|
[3787] | 173 | kxx=(double)kx*dkx_;
|
---|
| 174 | rep = resp(kxx, kyy);
|
---|
[3756] | 175 | if (crmask&&(kz==0)) mask(ky,kx)=((rep<1.e-8)?9.e9:(1./rep));
|
---|
| 176 | if (rep<1.e-8) fourAmp(kx, ky, kz) = complex<TF>(9.e9,0.);
|
---|
| 177 | else {
|
---|
[3787] | 178 | amp = 1./sqrt(rep)/sqrt(2.);
|
---|
[3756] | 179 | fourAmp(kx, ky, kz) = complex<TF>(rg_.Gaussian(amp), rg_.Gaussian(amp));
|
---|
| 180 | }
|
---|
| 181 | }
|
---|
| 182 | }
|
---|
| 183 | }
|
---|
| 184 | if (prtlev_>1) fourAmp.Show();
|
---|
| 185 | if (crmask) {
|
---|
| 186 | POutPersist po("mask.ppf");
|
---|
| 187 | po << mask;
|
---|
| 188 | }
|
---|
| 189 | if (prtlev_>0)
|
---|
| 190 | cout << " Four3DPk::ComputeNoiseFourierAmp() done ..." << endl;
|
---|
| 191 | }
|
---|
| 192 |
|
---|
| 193 | // Compute mass field from its Fourier Coefficient
|
---|
| 194 | TArray<TF> Four3DPk::ComputeMassDens()
|
---|
| 195 | {
|
---|
| 196 | TArray<TF> massdens;
|
---|
| 197 | // Backward fourier transform of the fourierAmp array
|
---|
| 198 | FFTWServer ffts(true);
|
---|
| 199 | ffts.setNormalize(true);
|
---|
| 200 | ffts.FFTBackward(fourAmp, massdens, true);
|
---|
| 201 | // cout << " Four3DPk::ComputeMassDens() done NbNeg=" << npbz << " / NPix=" << massDens.Size() << endl;
|
---|
| 202 | cout << " Four3DPk::ComputeMassDens() done NPix=" << massdens.Size() << endl;
|
---|
| 203 | return massdens;
|
---|
| 204 | }
|
---|
| 205 |
|
---|
| 206 | // Compute power spectrum as a function of wave number k
|
---|
| 207 | // cells with amp^2=re^2+im^2>s2cut are ignored
|
---|
| 208 | // Output : power spectrum (profile histogram)
|
---|
[3769] | 209 | HProf Four3DPk::ComputePk(double s2cut, int nbin, double kmin, double kmax)
|
---|
[3756] | 210 | {
|
---|
| 211 | // The second half of the array along Y (matrix rows) contain
|
---|
| 212 | // negative frequencies
|
---|
| 213 | // int nbh = sqrt(fourAmp.SizeX()*fourAmp.SizeX()+fourAmp.SizeY()*fourAmp.SizeY()/4.+fourAmp.SizeZ()*fourAmp.SizeY()/4.);
|
---|
| 214 | // The profile histogram will contain the mean value of FFT amplitude
|
---|
| 215 | // as a function of wave-number k = sqrt((double)(kx*kx+ky*ky))
|
---|
| 216 | // if (nbin < 1) nbin = nbh/2;
|
---|
[3769] | 217 | if ((kmax<0.)||(kmax<kmin)) {
|
---|
| 218 | kmin=0.;
|
---|
| 219 | double maxx=fourAmp.SizeX()*dkx_;
|
---|
| 220 | double maxy=fourAmp.SizeY()*dky_/2;
|
---|
| 221 | double maxz=fourAmp.SizeZ()*dkz_/2;
|
---|
| 222 | kmax=sqrt(maxx*maxx+maxy*maxy+maxz*maxz);
|
---|
| 223 | }
|
---|
| 224 | if (nbin<2) nbin=128;
|
---|
[3756] | 225 | HProf hp(kmin, kmax, nbin);
|
---|
| 226 | hp.SetErrOpt(false);
|
---|
| 227 | ComputePkCumul(hp, s2cut);
|
---|
| 228 | return hp;
|
---|
| 229 | }
|
---|
| 230 |
|
---|
| 231 | // Compute power spectrum as a function of wave number k
|
---|
| 232 | // Cumul dans hp - cells with amp^2=re^2+im^2>s2cut are ignored
|
---|
| 233 | void Four3DPk::ComputePkCumul(HProf& hp, double s2cut)
|
---|
| 234 | {
|
---|
| 235 |
|
---|
| 236 | // fourAmp represent 3-D fourier transform of a real input array.
|
---|
| 237 | // The second half of the array along Y and Z contain negative frequencies
|
---|
| 238 | double kxx, kyy, kzz;
|
---|
| 239 | // sa_size_t is large integer type
|
---|
[3783] | 240 | // We ignore 0th term in all frequency directions ...
|
---|
| 241 | for(sa_size_t kz=1; kz<fourAmp.SizeZ(); kz++) {
|
---|
[3756] | 242 | kzz = (kz > fourAmp.SizeZ()/2) ? (double)(fourAmp.SizeZ()-kz)*dkz_ : (double)kz*dkz_;
|
---|
[3783] | 243 | for(sa_size_t ky=1; ky<fourAmp.SizeY(); ky++) {
|
---|
[3756] | 244 | kyy = (ky > fourAmp.SizeY()/2) ? (double)(fourAmp.SizeY()-ky)*dky_ : (double)ky*dky_;
|
---|
[3783] | 245 | for(sa_size_t kx=1; kx<fourAmp.SizeX(); kx++) { // ignore the 0th coefficient (constant term)
|
---|
[3756] | 246 | double kxx=(double)kx*dkx_;
|
---|
| 247 | complex<TF> za = fourAmp(kx, ky, kz);
|
---|
| 248 | if (za.real()>8.e9) continue;
|
---|
| 249 | double wk = sqrt(kxx*kxx+kyy*kyy+kzz*kzz);
|
---|
| 250 | double amp2 = za.real()*za.real()+za.imag()*za.imag();
|
---|
| 251 | if ((s2cut>1.e-9)&&(amp2>s2cut)) continue;
|
---|
| 252 | hp.Add(wk, amp2);
|
---|
| 253 | }
|
---|
| 254 | }
|
---|
| 255 | }
|
---|
| 256 | return;
|
---|
| 257 | }
|
---|
| 258 |
|
---|
| 259 |
|
---|
| 260 |
|
---|
| 261 | //-----------------------------------------------------
|
---|
| 262 | // -- MassDist2D class : 2D mass distribution
|
---|
| 263 | //-----------------------------------------------------
|
---|
| 264 | // Constructor
|
---|
| 265 | MassDist2D::MassDist2D(GenericFunc& pk, int size, double meandens)
|
---|
| 266 | : pkSpec(pk) , sizeA((size>16)?size:16) , massDens(sizeA, sizeA),
|
---|
| 267 | meanRho(meandens) , fg_fourAmp(false) , fg_massDens(false)
|
---|
| 268 | {
|
---|
| 269 | }
|
---|
| 270 |
|
---|
| 271 | // To the computation job
|
---|
| 272 | void MassDist2D::Compute()
|
---|
| 273 | {
|
---|
| 274 | ComputeFourierAmp();
|
---|
| 275 | ComputeMassDens();
|
---|
| 276 | }
|
---|
| 277 |
|
---|
| 278 | // Generate mass field Fourier Coefficient
|
---|
| 279 | void MassDist2D::ComputeFourierAmp()
|
---|
| 280 | {
|
---|
| 281 | if (fg_fourAmp) return; // job already done
|
---|
| 282 | // We generate a random gaussian real field
|
---|
| 283 | double sigma = 1.;
|
---|
| 284 | // The following line fills the array by gaussian random numbers
|
---|
| 285 | //--Replaced-- massDens = RandomSequence(RandomSequence::Gaussian, 0., sigma);
|
---|
| 286 | // Can be replaced by
|
---|
| 287 | DR48RandGen rg;
|
---|
| 288 | for(sa_size_t ir=0; ir<massDens.NRows(); ir++) {
|
---|
| 289 | for(sa_size_t jc=0; jc<massDens.NCols(); jc++) {
|
---|
| 290 | massDens(ir, jc) = rg.Gaussian(sigma);
|
---|
| 291 | }
|
---|
| 292 | }
|
---|
| 293 | // --- End of random filling
|
---|
| 294 |
|
---|
| 295 | // Compute fourier transform of the random gaussian field -> white noise
|
---|
| 296 | FFTWServer ffts(true);
|
---|
| 297 | ffts.setNormalize(true);
|
---|
| 298 | ffts.FFTForward(massDens, fourAmp);
|
---|
| 299 |
|
---|
| 300 | // fourAmp represent 2-D fourier transform of a real input array.
|
---|
| 301 | // The second half of the array along Y (matrix rows) contain
|
---|
| 302 | // negative frequencies
|
---|
| 303 | // double fnorm = 1./sqrt(2.*fourAmp.Size());
|
---|
| 304 | // PUT smaller value for fnorm and check number of zeros
|
---|
| 305 | double fnorm = 1.;
|
---|
| 306 | // sa_size_t is large integer type
|
---|
| 307 | for(sa_size_t ky=0; ky<fourAmp.NRows(); ky++) {
|
---|
| 308 | double kyy = ky;
|
---|
| 309 | if (ky > fourAmp.NRows()/2) kyy = fourAmp.NRows()-ky; // negative frequencies
|
---|
| 310 | for(sa_size_t kx=0; kx<fourAmp.NCols(); kx++) {
|
---|
| 311 | double wk = sqrt((double)(kx*kx+kyy*kyy));
|
---|
| 312 | double amp = pkSpec(wk)*fnorm;
|
---|
| 313 | fourAmp(ky, kx) *= amp; // renormalize fourier coeff using
|
---|
| 314 | }
|
---|
| 315 | }
|
---|
| 316 | fg_fourAmp = true;
|
---|
| 317 | cout << " MassDist2D::ComputeFourierAmp() done ..." << endl;
|
---|
| 318 | }
|
---|
| 319 |
|
---|
| 320 | // Compute mass field from its Fourier Coefficient
|
---|
| 321 | void MassDist2D::ComputeMassDens()
|
---|
| 322 | {
|
---|
| 323 | if (fg_massDens) return; // job already done
|
---|
| 324 | if (!fg_fourAmp) ComputeFourierAmp(); // Check fourier amp generation
|
---|
| 325 |
|
---|
| 326 | // Backward fourier transform of the fourierAmp array
|
---|
| 327 | FFTWServer ffts(true);
|
---|
| 328 | ffts.setNormalize(true);
|
---|
| 329 | ffts.FFTBackward(fourAmp, massDens, true);
|
---|
| 330 | // We consider that massDens represents delta rho/rho
|
---|
| 331 | // rho = (delta rho/rho + 1) * MeanDensity
|
---|
| 332 | massDens += 1.;
|
---|
| 333 | // We remove negative values
|
---|
| 334 | sa_size_t npbz = 0;
|
---|
| 335 | for (sa_size_t i=0; i<massDens.NRows(); i++)
|
---|
| 336 | for (sa_size_t j=0; j<massDens.NCols(); j++)
|
---|
| 337 | if (massDens(i,j) < 0.) { npbz++; massDens(i,j) = 0.; }
|
---|
| 338 | massDens *= meanRho;
|
---|
| 339 | cout << " MassDist2D::ComputeMassDens() done NbNeg=" << npbz << " / NPix=" << massDens.Size() << endl;
|
---|
| 340 | }
|
---|
| 341 |
|
---|
| 342 | // Compute power spectrum as a function of wave number k
|
---|
| 343 | // Output : power spectrum (profile histogram)
|
---|
| 344 | HProf MassDist2D::ReconstructPk(int nbin)
|
---|
| 345 | {
|
---|
| 346 | // The second half of the array along Y (matrix rows) contain
|
---|
| 347 | // negative frequencies
|
---|
| 348 | int nbh = sqrt(2.0)*fourAmp.NCols();
|
---|
| 349 | // The profile histogram will contain the mean value of FFT amplitude
|
---|
| 350 | // as a function of wave-number k = sqrt((double)(kx*kx+ky*ky))
|
---|
| 351 | if (nbin < 1) nbin = nbh/2;
|
---|
| 352 | HProf hp(-0.5, nbh-0.5, nbin);
|
---|
| 353 | hp.SetErrOpt(false);
|
---|
| 354 |
|
---|
| 355 | for(int ky=0; ky<fourAmp.NRows(); ky++) {
|
---|
| 356 | double kyy = ky;
|
---|
| 357 | if (ky > fourAmp.NRows()/2) kyy = fourAmp.NRows()-ky; // negative frequencies
|
---|
| 358 | for(int kx=0; kx<fourAmp.NCols(); kx++) {
|
---|
| 359 | double wk = sqrt((double)(kx*kx+kyy*kyy));
|
---|
| 360 | complex<r_8> za = fourAmp(ky, kx);
|
---|
| 361 | double amp = sqrt(za.real()*za.real()+za.imag()*za.imag());
|
---|
| 362 | hp.Add(wk, amp);
|
---|
| 363 | }
|
---|
| 364 | }
|
---|
| 365 | return hp;
|
---|
| 366 | }
|
---|
| 367 |
|
---|