49 ClassImp(TFFTComplexReal);
54 TFFTComplexReal::TFFTComplexReal()
68 TFFTComplexReal::TFFTComplexReal(Int_t n, Bool_t inPlace)
71 fIn = fftw_malloc(
sizeof(fftw_complex)*(n/2+1));
72 fOut = fftw_malloc(
sizeof(Double_t)*n);
75 fIn = fftw_malloc(
sizeof(fftw_complex)*(n/2+1));
89 TFFTComplexReal::TFFTComplexReal(Int_t ndim, Int_t *n, Bool_t inPlace)
93 fN =
new Int_t[fNdim];
94 for (Int_t i=0; i<fNdim; i++){
98 Int_t sizein = Int_t(Double_t(fTotalSize)*(n[ndim-1]/2+1)/n[ndim-1]);
100 fIn = fftw_malloc(
sizeof(fftw_complex)*sizein);
101 fOut = fftw_malloc(
sizeof(Double_t)*fTotalSize);
103 fIn = fftw_malloc(
sizeof(fftw_complex)*sizein);
115 TFFTComplexReal::~TFFTComplexReal()
117 fftw_destroy_plan((fftw_plan)fPlan);
119 fftw_free((fftw_complex*)fIn);
148 void TFFTComplexReal::Init( Option_t *flags, Int_t ,
const Int_t* )
153 fftw_destroy_plan((fftw_plan)fPlan);
157 fPlan = (
void*)fftw_plan_dft_c2r(fNdim, fN,(fftw_complex*)fIn,(Double_t*)fOut, MapFlag(flags));
159 fPlan = (
void*)fftw_plan_dft_c2r(fNdim, fN, (fftw_complex*)fIn, (Double_t*)fIn, MapFlag(flags));
165 void TFFTComplexReal::Transform()
168 fftw_execute((fftw_plan)fPlan);
170 Error(
"Transform",
"transform was not initialized");
179 void TFFTComplexReal::GetPoints(Double_t *data, Bool_t fromInput)
const
182 Error(
"GetPoints",
"Input array has been destroyed");
185 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
186 std::copy(array,array+fTotalSize, data);
193 Double_t TFFTComplexReal::GetPointReal(Int_t ipoint, Bool_t fromInput)
const
196 Error(
"GetPointReal",
"Input array has been destroyed");
199 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
200 return array[ipoint];
207 Double_t TFFTComplexReal::GetPointReal(
const Int_t *ipoint, Bool_t fromInput)
const
210 Error(
"GetPointReal",
"Input array has been destroyed");
213 Int_t ireal = ipoint[0];
214 for (Int_t i=0; i<fNdim-1; i++)
215 ireal=fN[i+1]*ireal + ipoint[i+1];
217 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
225 void TFFTComplexReal::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput)
const
228 Error(
"GetPointComplex",
"Input array has been destroyed");
231 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
240 void TFFTComplexReal::GetPointComplex(
const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput)
const
243 Error(
"GetPointComplex",
"Input array has been destroyed");
246 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
248 Int_t ireal = ipoint[0];
249 for (Int_t i=0; i<fNdim-1; i++)
250 ireal=fN[i+1]*ireal + ipoint[i+1];
259 Double_t* TFFTComplexReal::GetPointsReal(Bool_t fromInput)
const
266 Error(
"GetPointsReal",
"Input array was destroyed");
269 return (Double_t*) ( (fOut) ? fOut : fIn );
276 void TFFTComplexReal::GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput)
const
279 Error(
"GetPointsComplex",
"Input array has been destroyed");
282 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
283 for (Int_t i=0; i<fTotalSize; i++){
293 void TFFTComplexReal::GetPointsComplex(Double_t *data, Bool_t fromInput)
const
296 Error(
"GetPointsComplex",
"Input array has been destroyed");
299 const Double_t * array = (
const Double_t*) ( (fOut) ? fOut : fIn );
300 for (Int_t i=0; i<fTotalSize; i+=2){
301 data[i] = array[i/2];
310 void TFFTComplexReal::SetPoint(Int_t ipoint, Double_t re, Double_t im)
312 if (ipoint <= fN[0]/2){
313 ((fftw_complex*)fIn)[ipoint][0] = re;
314 ((fftw_complex*)fIn)[ipoint][1] = im;
316 ((fftw_complex*)fOut)[2*(fN[0]/2)-ipoint][0] = re;
317 ((fftw_complex*)fOut)[2*(fN[0]/2)-ipoint][1] = -im;
325 void TFFTComplexReal::SetPoint(
const Int_t *ipoint, Double_t re, Double_t im)
327 Int_t ireal = ipoint[0];
328 for (Int_t i=0; i<fNdim-2; i++){
329 ireal=fN[i+1]*ireal + ipoint[i+1];
331 ireal = (fN[fNdim-1]/2+1)*ireal+ipoint[fNdim-1];
332 Int_t realN = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
335 Error(
"SetPoint",
"Illegal index value");
338 ((fftw_complex*)fIn)[ireal][0] = re;
339 ((fftw_complex*)fIn)[ireal][1] = im;
346 void TFFTComplexReal::SetPointComplex(Int_t ipoint, TComplex &c)
348 if (ipoint <= fN[0]/2){
349 ((fftw_complex*)fIn)[ipoint][0] = c.Re();
350 ((fftw_complex*)fIn)[ipoint][1] = c.Im();
352 ((fftw_complex*)fIn)[2*(fN[0]/2)-ipoint][0] = c.Re();
353 ((fftw_complex*)fIn)[2*(fN[0]/2)-ipoint][1] = -c.Im();
361 void TFFTComplexReal::SetPoints(
const Double_t *data)
363 Int_t sizein = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
365 for (Int_t i=0; i<2*(sizein); i+=2){
366 ((fftw_complex*)fIn)[i/2][0]=data[i];
367 ((fftw_complex*)fIn)[i/2][1]=data[i+1];
374 void TFFTComplexReal::SetPointsComplex(
const Double_t *re,
const Double_t *im)
376 Int_t sizein = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
377 for (Int_t i=0; i<sizein; i++){
378 ((fftw_complex*)fIn)[i][0]=re[i];
379 ((fftw_complex*)fIn)[i][1]=im[i];
390 UInt_t TFFTComplexReal::MapFlag(Option_t *flag)
394 if (opt.Contains(
"ES"))
395 return FFTW_ESTIMATE;
396 if (opt.Contains(
"M"))
398 if (opt.Contains(
"P"))
400 if (opt.Contains(
"EX"))
401 return FFTW_EXHAUSTIVE;
402 return FFTW_ESTIMATE;