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