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