Changeset 1391 in Sophya
- Timestamp:
- Feb 9, 2001, 6:09:13 PM (25 years ago)
- Location:
- trunk/SophyaExt/IFFTW
- Files:
-
- 2 edited
Legend:
- Unmodified
- Added
- Removed
-
trunk/SophyaExt/IFFTW/fftwserver.cc
r834 r1391 3 3 #include "FFTW/rfftw.h" 4 4 5 class FFTWServerPlan{ 5 6 #define MAXND_FFTW 5 7 8 class FFTWServerPlan { 6 9 public: 7 10 FFTWServerPlan(int n, fftw_direction dir, bool fgreal=false); 8 FFTWServerPlan(int n x, int ny, fftw_direction dir, bool fgreal=false);11 FFTWServerPlan(int nd, int * nxyz, fftw_direction dir, bool fgreal=false); 9 12 ~FFTWServerPlan(); 10 13 void Recreate(int n); 11 void Recreate(int nx, int ny); 12 13 int _n; 14 int _nx, _ny; 14 void Recreate(int nd, int * nxyz); 15 16 int _n; // Array dimension for 1-d arrays 17 int _nd; // Nb of dimensions for n-d arrays 18 int _nxyz[MAXND_FFTW]; // Array dimensions for n-d arrays 15 19 fftw_direction _dir; 16 20 … … 30 34 pnd = NULL; 31 35 rpnd = NULL; 32 _nx = _ny = -10; 36 _nd = -1; 37 for(int k=0; k<MAXND_FFTW; k++) _nxyz[k] = -10; 33 38 _n = n; 34 39 _dir = dir; … … 36 41 else p = fftw_create_plan(n, dir, FFTW_ESTIMATE); 37 42 } 38 FFTWServerPlan::FFTWServerPlan(int nx, int ny, fftw_direction dir, bool fgreal) 39 { 40 if ( (nx < 1) || (ny <1) ) 41 throw ParmError("FFTWServerPlan: Array size Nx or Ny <= 0 !"); 43 44 FFTWServerPlan::FFTWServerPlan(int nd, int * nxyz, fftw_direction dir, bool fgreal) 45 { 46 int k; 47 if (nd > MAXND_FFTW) 48 throw ParmError("FFTWServerPlan: Array rank (nd) > MAXND_FFTW !"); 42 49 p = NULL; 43 50 rp = NULL; … … 45 52 rpnd = NULL; 46 53 _n = -10; 47 _nx = nx; 48 _ny = ny; 54 _nd = nd; 55 for(k=0; k<nd; k++) { 56 if (nxyz[k] < 1) 57 throw ParmError("FFTWServerPlan: One of the Array size <= 0 !"); 58 _nxyz[k] = nxyz[k]; 59 } 60 for(k=nd; k<MAXND_FFTW; k++) _nxyz[k] = -10; 49 61 _dir = dir; 50 int sz[2]; 51 sz[0]= nx; sz[1] = ny; 52 if (fgreal) rpnd = rfftwnd_create_plan(2, sz, dir, FFTW_ESTIMATE); 53 else pnd = fftwnd_create_plan(2, sz, dir, FFTW_ESTIMATE); 62 if (fgreal) rpnd = rfftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE); 63 else pnd = fftwnd_create_plan(_nd, _nxyz, dir, FFTW_ESTIMATE); 54 64 } 55 65 … … 67 77 if (n < 1) 68 78 throw ParmError("FFTWServerPlan::Recreate(n) n < 0 !"); 69 if ( (_nx > 0) || (_ny > 0))70 throw ParmError("FFTWServerPlan::Recreate(n) Nx or Ny> 0 !");79 if (_nd > 0) 80 throw ParmError("FFTWServerPlan::Recreate(n) Multi-dimensional plan ! > 0 !"); 71 81 if (n == _n) return; 72 82 _n = n; … … 82 92 83 93 void 84 FFTWServerPlan::Recreate(int nx, int ny) 85 { 86 if ( (nx < 1) || (ny <1) ) 87 throw ParmError("FFTWServerPlan:Recreate(nx, ny) size Nx or Ny <= 0 !"); 94 FFTWServerPlan::Recreate(int nd, int * nxyz) 95 { 88 96 if (_n > 0) 89 throw ParmError("FFTWServerPlan::Recreate(nx, ny) N > 0 !"); 90 if ((nx == _nx) && (ny == _ny)) return; 91 _nx = nx; 92 _ny = ny; 93 int sz[2]; 94 sz[0]= nx; sz[1] = ny; 97 throw ParmError("FFTWServerPlan::Recreate(nd, nxyz) 1-dimensional plan !"); 98 int k; 99 if (nd == _nd) { 100 bool samepl = true; 101 for (int k=0; k<nd; k++) 102 if (nxyz[k] != _nxyz[k]) samepl = false; 103 if (samepl) return; 104 } 105 if (nd > MAXND_FFTW) 106 throw ParmError("FFTWServerPlan::Recreate(nd, nxyz) Array rank (nd) > MAXND_FFTW !"); 107 _nd = nd; 108 for(k=0; k<nd; k++) { 109 if (nxyz[k] < 1) 110 throw ParmError("FFTWServerPlan::Recreate(nd, nxyz) One of the Array size <= 0 !"); 111 _nxyz[k] = nxyz[k]; 112 } 113 for(k=nd; k<MAXND_FFTW; k++) _nxyz[k] = -10; 95 114 if (pnd) { 96 115 fftwnd_destroy_plan(pnd); 97 pnd = fftwnd_create_plan( 2, sz,_dir, FFTW_ESTIMATE);116 pnd = fftwnd_create_plan(_nd, _nxyz, _dir, FFTW_ESTIMATE); 98 117 } 99 118 else { 100 119 rfftwnd_destroy_plan(rpnd); 101 rpnd = rfftwnd_create_plan( 2, sz, _dir, FFTW_ESTIMATE);120 rpnd = rfftwnd_create_plan(_nd, _nxyz, _dir, FFTW_ESTIMATE); 102 121 } 103 122 … … 108 127 FFTWServer::FFTWServer() 109 128 : FFTServerInterface("FFTServer using FFTW package") 129 , ckR4(true, false) , ckR8(true, false) 130 110 131 { 111 132 _p1df = NULL; 112 133 _p1db = NULL; 113 _p 2df = NULL;114 _p 2db = NULL;134 _pndf = NULL; 135 _pndb = NULL; 115 136 116 137 _p1drf = NULL; 117 138 _p1drb = NULL; 118 _p 2drf = NULL;119 _p 2drb = NULL;139 _pndrf = NULL; 140 _pndrb = NULL; 120 141 } 121 142 … … 126 147 if (_p1df) delete _p1df ; 127 148 if (_p1db) delete _p1db ; 128 if (_p 2df) delete _p2df ;129 if (_p 2db) delete _p2db ;149 if (_pndf) delete _pndf ; 150 if (_pndb) delete _pndb ; 130 151 131 152 if (_p1drf) delete _p1drf ; 132 153 if (_p1drb) delete _p1drb ; 133 if (_p 2drf) delete _p2drf ;134 if (_p 2drb) delete _p2drb ;154 if (_pndrf) delete _pndrf ; 155 if (_pndrb) delete _pndrb ; 135 156 } 136 157 … … 143 164 /* --Methode-- */ 144 165 void 145 FFTWServer::FFTForward(TVector< complex<double> > const & in, TVector< complex<double> > & out) 146 { 147 if (_p1df) _p1df->Recreate(in.NElts()); 148 else _p1df = new FFTWServerPlan(in.NElts(), FFTW_FORWARD, false); 149 out.ReSize(in.NElts()); 150 fftw_one(_p1df->p, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 151 if(this->getNormalize()) out=out/complex<double>(pow(in.NElts(),0.5),0.); 152 } 153 /* --Methode-- */ 154 void FFTWServer::FFTBackward(TVector< complex<double> > const & in, TVector< complex<double> > & out) 155 { 156 if (_p1db) _p1db->Recreate(in.NElts()); 157 else _p1db = new FFTWServerPlan(in.NElts(), FFTW_BACKWARD, false); 158 out.ReSize(in.NElts()); 159 fftw_one(_p1db->p, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 160 if(this->getNormalize()) out=out/complex<double>(pow(in.NElts(),0.5),0.); 161 162 } 163 164 165 void FFTWServer::FFTForward(TVector< double > const & in, TVector< complex<double> > & out) 166 FFTWServer::FFTForward(TArray< complex<r_8> > const & in, TArray< complex<r_8> > & out) 167 { 168 int rank = ckR8.CheckResize(in, out); 169 if (rank == 1) { // One dimensional transform 170 if (_p1df) _p1df->Recreate(in.Size()); 171 else _p1df = new FFTWServerPlan(in.Size(), FFTW_FORWARD, false); 172 fftw_one(_p1df->p, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 173 } 174 else { // Multi dimensional 175 if (in.NbDimensions() > MAXND_FFTW) 176 throw ParmError("FFTWServer::FFTForward( complex<r_8>, complex<r_8> ) rank > MAXND_FFTW !"); 177 int sz[MAXND_FFTW]; 178 int k1 = 0; 179 int k2 = 0; 180 for(k1=in.NbDimensions()-1; k1>=0; k1--) { 181 sz[k2] = in.Size(k1); k2++; 182 } 183 if (_pndf) _pndf->Recreate(in.NbDimensions(), sz); 184 else _pndf = new FFTWServerPlan(in.NbDimensions(), sz, FFTW_FORWARD, false); 185 fftwnd_one(_pndf->pnd, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 186 } 187 // $CHECK$ Reza 9/2/2001 , Verifier normalisation 188 if(this->getNormalize()) out=out/complex<r_8>(sqrt((double)in.Size()),0.); 189 return; 190 } 191 192 /* --Methode-- */ 193 void FFTWServer::FFTBackward(TArray< complex<r_8> > const & in, TArray< complex<r_8> > & out) 194 { 195 int rank = ckR8.CheckResize(in, out); 196 if (rank == 1) { // One dimensional transform 197 if (_p1db) _p1db->Recreate(in.Size()); 198 else _p1db = new FFTWServerPlan(in.Size(), FFTW_BACKWARD, false); 199 fftw_one(_p1db->p, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 200 } 201 else { // Multi dimensional 202 if (in.NbDimensions() > MAXND_FFTW) 203 throw ParmError("FFTWServer::FFTForward( complex<r_8>, complex<r_8> ) rank > MAXND_FFTW !"); 204 int sz[MAXND_FFTW]; 205 int k1 = 0; 206 int k2 = 0; 207 for(k1=in.NbDimensions()-1; k1>=0; k1--) { 208 sz[k2] = in.Size(k1); k2++; 209 } 210 if (_pndb) _pndb->Recreate(in.NbDimensions(), sz); 211 else _pndb = new FFTWServerPlan(in.NbDimensions(), sz, FFTW_BACKWARD, false); 212 } 213 // $CHECK$ Reza 9/2/2001 , Verifier normalisation 214 if(this->getNormalize()) out=out/complex<r_8>(sqrt((double)in.Size()),0.); 215 return; 216 } 217 218 219 void FFTWServer::FFTForward(TArray< r_8 > const & in, TArray< complex<r_8> > & out) 166 220 { 167 int size = in.NElts()/2; 168 169 if(in.NElts()%2 != 0) { size = in.NElts()/2 +1;} 170 else { size = in.NElts()/2 +1 ;} 171 172 TVector< double > const outTemp(in.NElts()); 173 out.ReSize(size); 174 if (_p1drf) _p1drf->Recreate(in.NElts()); 175 else _p1drf = new FFTWServerPlan(in.NElts(), FFTW_REAL_TO_COMPLEX, true); 176 rfftw_one(_p1drf->rp, (fftw_real *)(in.Data()) , (fftw_real *)(outTemp.Data())); 177 ReShapetoCompl(outTemp, out); 178 if(this->getNormalize()) out=out/complex<double>(pow(in.NElts(),0.5),0.); 179 } 180 181 182 183 void FFTWServer::FFTBackward(TVector< complex<double> > const & in, TVector< double > & out) 184 { 221 int rank = ckR8.CheckResize(in, out); 222 TArray<r_8> outtemp(in, false); 223 224 if (rank == 1) { // One dimensional transform 225 if (_p1drf) _p1drf->Recreate(in.Size()); 226 else _p1drf = new FFTWServerPlan(in.Size(), FFTW_REAL_TO_COMPLEX, true); 227 rfftw_one(_p1drf->rp, (fftw_real *)(in.Data()) , (fftw_real *)(outtemp.Data())); 228 ReShapetoCompl(outtemp, out); 229 } 230 else { // Multi dimensional 231 if (in.NbDimensions() > MAXND_FFTW) 232 throw ParmError("FFTWServer::FFTForward( complex<r_8>, complex<r_8> ) rank > MAXND_FFTW !"); 233 int sz[MAXND_FFTW]; 234 int k1 = 0; 235 int k2 = 0; 236 for(k1=in.NbDimensions()-1; k1>=0; k1--) { 237 sz[k2] = in.Size(k1); k2++; 238 } 239 if (_pndrf) _pndrf->Recreate(in.NbDimensions(), sz); 240 else _pndrf = new FFTWServerPlan(in.NbDimensions(), sz, FFTW_REAL_TO_COMPLEX, true); 241 rfftwnd_one_real_to_complex(_pndrf->rpnd, (fftw_real *)(in.Data()) , 242 (fftw_complex *)(out.Data()) ); 243 } 244 // $CHECK$ Reza 9/2/2001 , Verifier normalisation 245 if(this->getNormalize()) out=out/complex<r_8>(sqrt((double)in.Size()),0.); 246 return; 247 248 } 249 250 251 252 void FFTWServer::FFTBackward(TArray< complex<r_8> > const & in, TArray< r_8 > & out) 253 { 254 throw ParmError("FFTWServer::FFTBackward(TArray< complex<r_8> > ... Not implemented ... !"); 255 /* 185 256 int size; 186 257 if(in(in.NElts()).imag() == 0) { size = 2*in.NElts()-2;} 187 258 else { size = 2*in.NElts()-1;} 188 259 189 T Vector< double> inTemp(size);260 TArray< r_8 > inTemp(size); 190 261 out.ReSize(size); 191 262 … … 196 267 rfftw_one(_p1drb->rp, (fftw_real *)(inTemp.Data()) , (fftw_real *)(out.Data())); 197 268 if(this->getNormalize()) out=out/pow(size,0.5); 198 } 199 200 /* --Methode-- */ 201 void FFTWServer::FFTForward(TMatrix< complex<double> > const & in, TMatrix< complex<double> > & out) 269 */ 270 } 271 272 /* 273 void FFTWServer::FFTForward(TArray< complex<r_8> > const & in, TArray< complex<r_8> > & out) 202 274 { 203 275 out.ReSize(in.NRows(),in.NCols()); 204 276 205 if (_p2df) _p2df->Recreate( in.NRows(),in.NCols()); 206 else _p2df = new FFTWServerPlan( in.NCols(),in.NRows(), FFTW_FORWARD, false); 207 208 fftwnd_one(_p2df->pnd, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 209 if(this->getNormalize()) out=out/complex<double>(pow(in.NRows()*in.NCols(),0.5),0.); 210 } 211 212 /* --Methode-- */ 213 void FFTWServer::FFTBackward(TMatrix< complex<double> > const & in, TMatrix< complex<double> > & out) 214 { 215 if (_p2db) _p2db->Recreate(in.NCols(), in.NRows()); 216 else _p2db = new FFTWServerPlan(in.NCols(), in.NRows(), FFTW_BACKWARD, false); 277 if (_pndf) _pndf->Recreate( in.NRows(),in.NCols()); 278 else _pndf = new FFTWServerPlan( in.NCols(),in.NRows(), FFTW_FORWARD, false); 279 280 fftwnd_one(_pndf->pnd, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 281 if(this->getNormalize()) out=out/complex<r_8>(pow(in.NRows()*in.NCols(),0.5),0.); 282 } 283 284 void FFTWServer::FFTBackward(TArray< complex<r_8> > const & in, TArray< complex<r_8> > & out) 285 { 286 if (_pndb) _pndb->Recreate(in.NCols(), in.NRows()); 287 else _pndb = new FFTWServerPlan(in.NCols(), in.NRows(), FFTW_BACKWARD, false); 217 288 out.ReSize(in.NRows(), in.NCols()); 218 fftwnd_one(_p2db->pnd, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 219 if(this->getNormalize()) out=out/complex<double>(pow(in.NRows()*in.NCols(),0.5),0.); 220 221 } 222 223 224 /* --Methode-- */ 225 void FFTWServer::FFTForward(TMatrix< double > const & in, TMatrix< complex<double> > & out) 226 { 227 228 TMatrix< double > inNew(in.NCols(),in.NRows()); 289 fftwnd_one(_pndb->pnd, (fftw_complex *)(in.Data()) , (fftw_complex *)(out.Data()) ); 290 if(this->getNormalize()) out=out/complex<r_8>(pow(in.NRows()*in.NCols(),0.5),0.); 291 292 } 293 294 295 void FFTWServer::FFTForward(TArray< r_8 > const & in, TArray< complex<r_8> > & out) 296 { 297 298 TArray< r_8 > inNew(in.NCols(),in.NRows()); 229 299 for(int i=0; i<in.NRows(); i++) 230 300 for(int j=0;j<in.NCols(); j++) 231 301 inNew(j,i) = in(i,j); 232 302 233 if (_p 2drf) _p2drf->Recreate(inNew.NRows(),inNew.NCols());234 else _p 2drf = new FFTWServerPlan(inNew.NRows(), inNew.NCols(),FFTW_REAL_TO_COMPLEX, true);303 if (_pndrf) _pndrf->Recreate(inNew.NRows(),inNew.NCols()); 304 else _pndrf = new FFTWServerPlan(inNew.NRows(), inNew.NCols(),FFTW_REAL_TO_COMPLEX, true); 235 305 // rfftwnd_plan p; 236 T Matrix< complex<double> > outTemp;306 TArray< complex<r_8> > outTemp; 237 307 outTemp.ReSize(in.NRows(),in.NCols()); 238 308 239 rfftwnd_one_real_to_complex(_p2drf->rpnd, (fftw_real *)(in.Data()) , (fftw_complex *)(out.Data()) ); 240 } 241 242 /* --Methode-- */ 243 void FFTWServer::FFTBackward(TMatrix< complex<double> > const & in, TMatrix< double > & out) 244 { 245 246 TMatrix< complex<double> > inNew(in.NCols(),in.NRows()); 309 rfftwnd_one_real_to_complex(_pndrf->rpnd, (fftw_real *)(in.Data()) , (fftw_complex *)(out.Data()) ); 310 } 311 312 void FFTWServer::FFTBackward(TArray< complex<r_8> > const & in, TArray< r_8 > & out) 313 { 314 315 TArray< complex<r_8> > inNew(in.NCols(),in.NRows()); 247 316 for(int i=0; i<in.NRows(); i++) 248 317 for(int j=0;j<in.NCols(); j++) 249 318 inNew(j,i) = in(i,j); 250 319 251 if (_p 2drb) _p2drb->Recreate(inNew.NRows(),inNew.NCols());252 else _p 2drb = new FFTWServerPlan(inNew.NRows(), inNew.NCols(),FFTW_COMPLEX_TO_REAL, true);320 if (_pndrb) _pndrb->Recreate(inNew.NRows(),inNew.NCols()); 321 else _pndrb = new FFTWServerPlan(inNew.NRows(), inNew.NCols(),FFTW_COMPLEX_TO_REAL, true); 253 322 // rfftwnd_plan p; 254 323 out.ReSize(in.NRows(),in.NCols()); 255 324 256 rfftwnd_one_complex_to_real(_p 2drb->rpnd, (fftw_complex *)(in.Data()) , (fftw_real *)(out.Data()) );325 rfftwnd_one_complex_to_real(_pndrb->rpnd, (fftw_complex *)(in.Data()) , (fftw_real *)(out.Data()) ); 257 326 cout << " in the function !!!" << endl; 258 327 if(this->getNormalize()) 259 328 { 260 double norm = (double)(in.NRows()*in.NCols());329 r_8 norm = (r_8)(in.NRows()*in.NCols()); 261 330 out=out/norm; 262 331 } 263 332 } 264 333 265 266 /* --Methode-- */ 267 void FFTWServer::ReShapetoReal( TVector< complex<double> > const & in, TVector< double > & out) 268 { 269 int N = in.NElts(); 334 */ 335 336 337 /* --Methode-- */ 338 void FFTWServer::ReShapetoReal( TArray< complex<r_8> > const & in, TArray< r_8 > & out) 339 { 340 int N = in.Size(); 341 /* 270 342 int i; 271 if (in(in. NElts()).imag() == 0)343 if (in(in.Size()).imag() == 0) 272 344 { 273 345 out(0) = in(0).real(); … … 293 365 } 294 366 } 367 */ 368 out[0] = in[0].real(); 369 int k=0; 370 for(int i=1; i<N; i++) { 371 out[i] = in[i].real(); 372 out[N-i] = in[i].imag(); 373 } 374 295 375 // for(int k=0; k<out.NElts(); k++) cout << "ReShapetoReal out " << k << " " << out(k) << endl; 296 376 } … … 298 378 299 379 /* --Methode-- */ 300 void FFTWServer::ReShapetoCompl(T Vector< double > const & in, TVector< complex<double> > & out)301 { 302 int N = in. NElts();380 void FFTWServer::ReShapetoCompl(TArray< r_8 > const & in, TArray< complex<r_8> > & out) 381 { 382 int N = in.Size(); 303 383 // for(int k=0; k<in.NElts(); k++) cout << "ReShapetoCompl in " << k << " " << in(k) << endl; 304 out(0) = complex<double> (in(0),0.); 305 if(in.NElts()%2 !=0) 306 { 307 for(int k=1;k<=N/2+1;k++) 308 { 309 out(k) = complex<double> (in(k),in(N-k)); 310 } 311 } 312 else 313 { 314 for(int k=1;k<N/2;k++) 315 { 316 out(k) = complex<double> (in(k),in(N-k)); 317 } 318 out(N/2) = complex<double> (in(N/2),0.); 319 } 384 out[0] = complex<r_8> (in[0],0.); 385 for(int k=1; k<N+1/2; k++) 386 out[k] = complex<r_8>(in[k], in[N-k]); 387 if (N%2 == 0) out[N/2] = complex<r_8>(in[N/2], 0.); 320 388 // for(int k=0; k<out.NElts(); k++) cout << "ReShapetoCompl out " << k << " " << out(k) << endl; 321 389 } -
trunk/SophyaExt/IFFTW/fftwserver.h
r765 r1391 19 19 virtual FFTServerInterface * Clone(); 20 20 21 // Transformee unidimensionnelle 22 virtual void FFTForward(TVector< complex<double> > const & in, TVector< complex<double> > & out); 23 virtual void FFTBackward(TVector< complex<double> > const & in, TVector< complex<double> > & out); 24 virtual void FFTForward(TVector< double > const & in, TVector< complex<double> > & out); 25 virtual void FFTBackward(TVector< complex<double> > const & in, TVector< double > & out); 26 27 // Transforme 2D 28 virtual void FFTForward(TMatrix< complex<double> > const & in, TMatrix< complex<double> > & out); 29 virtual void FFTBackward(TMatrix< complex<double> > const & in, TMatrix< complex<double> > & out); 30 virtual void FFTForward(TMatrix< double > const & in, TMatrix< complex<double> > & out); 31 virtual void FFTBackward(TMatrix< complex<double> > const & in, TMatrix< double > & out); 21 // Transforme unidimensionnelle , N-dimensionnel 22 virtual void FFTForward(TArray< complex<r_8> > const & in, TArray< complex<r_8> > & out); 23 virtual void FFTBackward(TArray< complex<r_8> > const & in, TArray< complex<r_8> > & out); 24 virtual void FFTForward(TArray< r_8 > const & in, TArray< complex<r_8> > & out); 25 virtual void FFTBackward(TArray< complex<r_8> > const & in, TArray< r_8 > & out); 32 26 33 27 // Methodes statiques pour reordonner les donnees 34 virtual void ReShapetoReal( T Vector< complex<r_8> > const & in, TVector< r_8 > & out);35 virtual void ReShapetoCompl(T Vector< r_8 > const & in, TVector< complex<r_8> > & out);28 virtual void ReShapetoReal( TArray< complex<r_8> > const & in, TArray< r_8 > & out); 29 virtual void ReShapetoCompl(TArray< r_8 > const & in, TArray< complex<r_8> > & out); 36 30 37 31 protected: 38 32 FFTWServerPlan * _p1df; 39 33 FFTWServerPlan * _p1db; 40 FFTWServerPlan * _p 2df;41 FFTWServerPlan * _p 2db;34 FFTWServerPlan * _pndf; 35 FFTWServerPlan * _pndb; 42 36 43 37 FFTWServerPlan * _p1drf; 44 38 FFTWServerPlan * _p1drb; 45 FFTWServerPlan * _p2drf; 46 FFTWServerPlan * _p2drb; 39 FFTWServerPlan * _pndrf; 40 FFTWServerPlan * _pndrb; 41 42 FFTArrayChecker<r_4> ckR4; 43 FFTArrayChecker<r_8> ckR8; 47 44 }; 48 45
Note:
See TracChangeset
for help on using the changeset viewer.