TFFTRealComplex.cxx

Go to the documentation of this file.
00001 // @(#)root/fft:$Id: TFFTRealComplex.cxx 36024 2010-10-01 16:03:30Z moneta $
00002 // Author: Anna Kreshuk   07/4/2006
00003 
00004 /*************************************************************************
00005  * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers.               *
00006  * All rights reserved.                                                  *
00007  *                                                                       *
00008  * For the licensing terms see $ROOTSYS/LICENSE.                         *
00009  * For the list of contributors see $ROOTSYS/README/CREDITS.             *
00010  *************************************************************************/
00011 
00012 //////////////////////////////////////////////////////////////////////////
00013 //                                                                      
00014 // TFFTRealComplex                                                       
00015 //                                                                      
00016 // One of the interface classes to the FFTW package, can be used directly
00017 // or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
00018 //
00019 // Computes a real input/complex output discrete Fourier transform in 1 or more
00020 // dimensions. However, only out-of-place transforms are now supported for transforms
00021 // in more than 1 dimension. For detailed information about the computed transforms,
00022 // please refer to the FFTW manual
00023 //
00024 // How to use it:
00025 // 1) Create an instance of TFFTRealComplex - this will allocate input and output
00026 //    arrays (unless an in-place transform is specified)
00027 // 2) Run the Init() function with the desired flags and settings (see function
00028 //    comments for possible kind parameters)
00029 // 3) Set the data (via SetPoints()or SetPoint() functions)
00030 // 4) Run the Transform() function
00031 // 5) Get the output (via GetPoints() or GetPoint() functions)
00032 // 6) Repeat steps 3)-5) as needed
00033 // For a transform of the same size, but with different flags, 
00034 // rerun the Init() function and continue with steps 3)-5)
00035 //
00036 // NOTE: 1) running Init() function will overwrite the input array! Don't set any data
00037 //          before running the Init() function
00038 //       2) FFTW computes unnormalized transform, so doing a transform followed by 
00039 //          its inverse will lead to the original array scaled by the transform size
00040 // 
00041 //
00042 //////////////////////////////////////////////////////////////////////////
00043 
00044 #include "TFFTRealComplex.h"
00045 #include "fftw3.h"
00046 #include "TComplex.h"
00047 
00048 
00049 ClassImp(TFFTRealComplex)
00050 
00051 //_____________________________________________________________________________
00052 TFFTRealComplex::TFFTRealComplex()
00053 {
00054 //default
00055 
00056    fIn   = 0;
00057    fOut  = 0;
00058    fPlan = 0;
00059    fN    = 0;
00060    fFlags = 0;
00061    fNdim = 0;
00062    fTotalSize = 0;
00063 }
00064 
00065 //_____________________________________________________________________________
00066 TFFTRealComplex::TFFTRealComplex(Int_t n, Bool_t inPlace)
00067 {
00068 //For 1d transforms
00069 //Allocates memory for the input array, and, if inPlace = kFALSE, for the output array
00070 
00071 
00072    if (!inPlace){
00073       fIn = fftw_malloc(sizeof(Double_t)*n);
00074       fOut = fftw_malloc(sizeof(fftw_complex)*(n/2+1));
00075    } else {
00076       fIn = fftw_malloc(sizeof(Double_t)*(2*(n/2+1)));
00077       fOut = 0;
00078    }
00079    fN = new Int_t[1];
00080    fN[0] = n;
00081    fTotalSize = n;
00082    fNdim = 1;
00083    fPlan = 0;
00084    fFlags = 0;
00085 }
00086 
00087 //_____________________________________________________________________________
00088 TFFTRealComplex::TFFTRealComplex(Int_t ndim, Int_t *n, Bool_t inPlace)
00089 {
00090 //For ndim-dimensional transforms
00091 //Second argurment contains sizes of the transform in each dimension
00092 
00093    if (ndim>1 && inPlace==kTRUE){
00094       Error("TFFTRealComplex", "multidimensional in-place r2c transforms are not implemented yet");
00095       return;
00096    }
00097    fNdim = ndim;
00098    fTotalSize = 1;
00099    fN = new Int_t[fNdim];
00100    for (Int_t i=0; i<fNdim; i++){
00101       fN[i] = n[i];
00102       fTotalSize*=n[i];
00103    }
00104    Int_t sizeout = Int_t(Double_t(fTotalSize)*(n[ndim-1]/2+1)/n[ndim-1]);
00105    if (!inPlace){
00106       fIn = fftw_malloc(sizeof(Double_t)*fTotalSize);
00107       fOut = fftw_malloc(sizeof(fftw_complex)*sizeout);
00108    } else {
00109       fIn = fftw_malloc(sizeof(Double_t)*(2*sizeout));
00110       fOut = 0;
00111    }
00112    fPlan = 0;
00113    fFlags = 0;
00114 }
00115 
00116 //_____________________________________________________________________________
00117 TFFTRealComplex::~TFFTRealComplex()
00118 {
00119 //Destroys the data arrays and the plan. However, some plan information stays around
00120 //until the root session is over, and is reused if other plans of the same size are
00121 //created
00122 
00123    fftw_destroy_plan((fftw_plan)fPlan);
00124    fPlan = 0;
00125    fftw_free(fIn);
00126    fIn = 0;
00127    if (fOut)
00128       fftw_free((fftw_complex*)fOut);
00129    fOut = 0;
00130    if (fN)
00131       delete [] fN;
00132    fN = 0;
00133 }
00134 
00135 //_____________________________________________________________________________
00136 void TFFTRealComplex::Init(Option_t *flags,Int_t /*sign*/, const Int_t* /*kind*/)
00137 {
00138 //Creates the fftw-plan
00139 //
00140 //NOTE:  input and output arrays are overwritten during initialisation,
00141 //       so don't set any points, before running this function!!!!!
00142 //
00143 //Arguments sign and kind are dummy and not need to be specified
00144 //Possible flag_options:
00145 //"ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
00146 //   performanc
00147 //"M" (from "measure") - some time spend in finding the optimal way to do the transform
00148 //"P" (from "patient") - more time spend in finding the optimal way to do the transform
00149 //"EX" (from "exhaustive") - the most optimal way is found
00150 //This option should be chosen depending on how many transforms of the same size and
00151 //type are going to be done. Planning is only done once, for the first transform of this
00152 //size and type.
00153 
00154    fFlags = flags;
00155 
00156    if (fOut)
00157       fPlan = (void*)fftw_plan_dft_r2c(fNdim, fN, (Double_t*)fIn, (fftw_complex*)fOut,MapFlag(flags));
00158    else
00159       fPlan = (void*)fftw_plan_dft_r2c(fNdim, fN, (Double_t*)fIn, (fftw_complex*)fIn, MapFlag(flags));
00160 }
00161 
00162 //_____________________________________________________________________________
00163 void TFFTRealComplex::Transform()
00164 {
00165 //Computes the transform, specified in Init() function
00166 
00167 
00168    if (fPlan){
00169       fftw_execute((fftw_plan)fPlan);
00170    }
00171    else {
00172       Error("Transform", "transform hasn't been initialised");
00173       return;
00174    }
00175 }
00176 
00177 //_____________________________________________________________________________
00178 void TFFTRealComplex::GetPoints(Double_t *data, Bool_t fromInput) const
00179 {
00180 //Fills the array data with the computed transform.
00181 //Only (roughly) a half of the transform is copied (exactly the output of FFTW),
00182 //the rest being Hermitian symmetric with the first half
00183 
00184    if (fromInput){
00185       for (Int_t i=0; i<fTotalSize; i++)
00186          data[i] = ((Double_t*)fIn)[i];
00187    } else {
00188       Int_t realN = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00189       if (fOut){
00190          for (Int_t i=0; i<realN; i+=2){
00191             data[i] = ((fftw_complex*)fOut)[i/2][0];
00192             data[i+1] = ((fftw_complex*)fOut)[i/2][1];
00193          }
00194       }
00195       else {
00196          for (Int_t i=0; i<realN; i++)
00197             data[i] = ((Double_t*)fIn)[i];
00198       }
00199    }
00200 }
00201 
00202 //_____________________________________________________________________________
00203 Double_t TFFTRealComplex::GetPointReal(Int_t ipoint, Bool_t fromInput) const
00204 {
00205 //Returns the real part of the point #ipoint from the output or the point #ipoint
00206 //from the input
00207 
00208    if (fOut && !fromInput){
00209       Warning("GetPointReal", "Output is complex. Only real part returned");
00210       return ((fftw_complex*)fOut)[ipoint][0];
00211    }
00212    else
00213       return ((Double_t*)fIn)[ipoint];
00214 }
00215 
00216 //_____________________________________________________________________________
00217 Double_t TFFTRealComplex::GetPointReal(const Int_t *ipoint, Bool_t fromInput) const
00218 {
00219 //Returns the real part of the point #ipoint from the output or the point #ipoint
00220 //from the input
00221 
00222    Int_t ireal = ipoint[0];
00223    for (Int_t i=0; i<fNdim-1; i++)
00224       ireal=fN[i+1]*ireal + ipoint[i+1];
00225 
00226     if (fOut && !fromInput){
00227       Warning("GetPointReal", "Output is complex. Only real part returned");
00228       return ((fftw_complex*)fOut)[ireal][0];
00229    }
00230    else
00231       return ((Double_t*)fIn)[ireal];
00232 }
00233 
00234 
00235 //_____________________________________________________________________________
00236 void TFFTRealComplex::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00237 {
00238 //Returns the point #ipoint.
00239 //For 1d, if ipoint > fN/2+1 (the point is in the Hermitian symmetric part), it is still
00240 //returned. For >1d, only the first (roughly)half of points can be returned
00241 //For 2d, see function GetPointComplex(Int_t *ipoint,...)
00242 
00243    if (fromInput){
00244       re = ((Double_t*)fIn)[ipoint];
00245    } else {
00246       if (fNdim==1){
00247          if (fOut){
00248             if (ipoint<fN[0]/2+1){
00249                re = ((fftw_complex*)fOut)[ipoint][0];
00250                im = ((fftw_complex*)fOut)[ipoint][1];
00251             } else {
00252                re = ((fftw_complex*)fOut)[fN[0]-ipoint][0];
00253                im = -((fftw_complex*)fOut)[fN[0]-ipoint][1];
00254             }
00255          } else {
00256             if (ipoint<fN[0]/2+1){
00257                re = ((Double_t*)fIn)[2*ipoint];
00258                im = ((Double_t*)fIn)[2*ipoint+1];
00259             } else {
00260                re = ((Double_t*)fIn)[2*(fN[0]-ipoint)];
00261                im = ((Double_t*)fIn)[2*(fN[0]-ipoint)+1];
00262             }
00263          }
00264       }
00265       else {
00266          Int_t realN = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00267          if (ipoint>realN/2){
00268             Error("GetPointComplex", "Illegal index value");
00269             return;
00270          }
00271          if (fOut){
00272             re = ((fftw_complex*)fOut)[ipoint][0];
00273             im = ((fftw_complex*)fOut)[ipoint][1];
00274          } else {
00275             re = ((Double_t*)fIn)[2*ipoint];
00276             im = ((Double_t*)fIn)[2*ipoint+1];
00277          }
00278       }
00279    }
00280 }
00281 //_____________________________________________________________________________
00282 void TFFTRealComplex::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00283 {
00284 //For multidimensional transforms. Returns the point #ipoint.
00285 //In case of transforms of more than 2 dimensions,
00286 //only points from the first (roughly)half are returned, the rest being Hermitian symmetric
00287 
00288 
00289    Int_t ireal = ipoint[0];
00290    for (Int_t i=0; i<fNdim-2; i++)
00291       ireal=fN[i+1]*ireal + ipoint[i+1];
00292    //special treatment of the last dimension
00293    ireal = (fN[fNdim-1]/2+1)*ireal + ipoint[fNdim-1];
00294 
00295    if (fromInput){
00296       re = ((Double_t*)fIn)[ireal];
00297       return;
00298    } 
00299    if (fNdim==1){
00300       if (fOut){
00301          if (ipoint[0] <fN[0]/2+1){
00302             re = ((fftw_complex*)fOut)[ipoint[0]][0];
00303             im = ((fftw_complex*)fOut)[ipoint[0]][1];
00304          } else {
00305             re = ((fftw_complex*)fOut)[fN[0]-ipoint[0]][0];
00306             im = -((fftw_complex*)fOut)[fN[0]-ipoint[0]][1];
00307          }
00308       } else {
00309          if (ireal <fN[0]/2+1){
00310             re = ((Double_t*)fIn)[2*ipoint[0]];
00311             im = ((Double_t*)fIn)[2*ipoint[0]+1];
00312          } else {
00313             re = ((Double_t*)fIn)[2*(fN[0]-ipoint[0])];
00314             im = ((Double_t*)fIn)[2*(fN[0]-ipoint[0])+1];
00315          }
00316       }
00317    }
00318    else if (fNdim==2){
00319       if (fOut){
00320          if (ipoint[1]<fN[1]/2+1){
00321             re = ((fftw_complex*)fOut)[ipoint[1]+(fN[1]/2+1)*ipoint[0]][0];
00322             im = ((fftw_complex*)fOut)[ipoint[1]+(fN[1]/2+1)*ipoint[0]][1];
00323          } else {
00324             if (ipoint[0]==0){
00325                re = ((fftw_complex*)fOut)[fN[1]-ipoint[1]][0];
00326                im = -((fftw_complex*)fOut)[fN[1]-ipoint[1]][1];
00327             } else {
00328                re = ((fftw_complex*)fOut)[(fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0])][0];
00329                im = -((fftw_complex*)fOut)[(fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0])][1];
00330             }
00331          }
00332       } else {
00333          if (ipoint[1]<fN[1]/2+1){
00334             re = ((Double_t*)fIn)[2*(ipoint[1]+(fN[1]/2+1)*ipoint[0])];
00335             im = ((Double_t*)fIn)[2*(ipoint[1]+(fN[1]/2+1)*ipoint[0])+1];
00336          } else {
00337             if (ipoint[0]==0){
00338                re = ((Double_t*)fIn)[2*(fN[1]-ipoint[1])];
00339                im = -((Double_t*)fIn)[2*(fN[1]-ipoint[1])+1];
00340             } else {
00341                re = ((Double_t*)fIn)[2*((fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0]))];
00342                im = -((Double_t*)fIn)[2*((fN[1]-ipoint[1])+(fN[1]/2+1)*(fN[0]-ipoint[0]))+1];
00343             }
00344          }
00345       }
00346    }
00347    else {
00348       
00349       if (fOut){
00350          re = ((fftw_complex*)fOut)[ireal][0];
00351          im = ((fftw_complex*)fOut)[ireal][1];
00352       } else {
00353          re = ((Double_t*)fIn)[2*ireal];
00354          im = ((Double_t*)fIn)[2*ireal+1];
00355       }
00356    }
00357 }
00358 
00359 
00360 //_____________________________________________________________________________
00361 Double_t* TFFTRealComplex::GetPointsReal(Bool_t fromInput) const
00362 {
00363 //Returns the input array// One of the interface classes to the FFTW package, can be used directly
00364 // or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
00365 
00366    if (!fromInput){
00367       Error("GetPointsReal", "Output array is complex");
00368       return 0;
00369    }
00370    return (Double_t*)fIn;
00371 }
00372 
00373 //_____________________________________________________________________________
00374 void TFFTRealComplex::GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput) const
00375 {
00376 //Fills the argument arrays with the real and imaginary parts of the computed transform.
00377 //Only (roughly) a half of the transform is copied, the rest being Hermitian
00378 //symmetric with the first half
00379 
00380    Int_t nreal;
00381    if (fOut && !fromInput){
00382       nreal = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00383       for (Int_t i=0; i<nreal; i++){
00384          re[i] = ((fftw_complex*)fOut)[i][0];
00385          im[i] = ((fftw_complex*)fOut)[i][1];
00386       }
00387    } else if (fromInput) {
00388       for (Int_t i=0; i<fTotalSize; i++){
00389          re[i] = ((Double_t *)fIn)[i];
00390          im[i] = 0;
00391       }
00392    }
00393    else {
00394       nreal = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00395       for (Int_t i=0; i<nreal; i+=2){
00396          re[i/2] = ((Double_t*)fIn)[i];
00397          im[i/2] = ((Double_t*)fIn)[i+1];
00398       }
00399    }
00400 }
00401 
00402 //_____________________________________________________________________________
00403 void TFFTRealComplex::GetPointsComplex(Double_t *data, Bool_t fromInput) const
00404 {
00405 //Fills the argument arrays with the real and imaginary parts of the computed transform.
00406 //Only (roughly) a half of the transform is copied, the rest being Hermitian
00407 //symmetric with the first half
00408 
00409    Int_t nreal;
00410 
00411    if (fOut && !fromInput){
00412       nreal = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00413 
00414       for (Int_t i=0; i<nreal; i+=2){
00415          data[i] = ((fftw_complex*)fOut)[i/2][0];
00416          data[i+1] = ((fftw_complex*)fOut)[i/2][1];
00417       }
00418    } else if (fromInput){
00419       for (Int_t i=0; i<fTotalSize; i+=2){
00420          data[i] = ((Double_t*)fIn)[i/2];
00421          data[i+1] = 0;
00422       }
00423    }
00424    else {
00425 
00426       nreal = 2*Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00427       for (Int_t i=0; i<nreal; i++)
00428          data[i] = ((Double_t*)fIn)[i];
00429    }
00430 }
00431 
00432 //_____________________________________________________________________________
00433 void TFFTRealComplex::SetPoint(Int_t ipoint, Double_t re, Double_t /*im*/)
00434 {
00435 //Set the point #ipoint
00436 
00437    ((Double_t *)fIn)[ipoint] = re;
00438 }
00439 
00440 //_____________________________________________________________________________
00441 void TFFTRealComplex::SetPoint(const Int_t *ipoint, Double_t re, Double_t /*im*/)
00442 {
00443 //For multidimensional transforms. Set the point #ipoint
00444 
00445    Int_t ireal = ipoint[0];
00446    for (Int_t i=0; i<fNdim-1; i++)
00447       ireal=fN[i+1]*ireal + ipoint[i+1];
00448    ((Double_t*)fIn)[ireal]=re;
00449 }
00450 
00451 //_____________________________________________________________________________
00452 void TFFTRealComplex::SetPoints(const Double_t *data)
00453 {
00454 //Set all input points
00455 
00456    for (Int_t i=0; i<fTotalSize; i++){
00457       ((Double_t*)fIn)[i]=data[i];
00458    }
00459 }
00460 
00461 //_____________________________________________________________________________
00462 void TFFTRealComplex::SetPointComplex(Int_t ipoint, TComplex &c)
00463 {
00464 //Sets the point #ipoint (only the real part of the argument is taken)
00465 
00466    ((Double_t *)fIn)[ipoint]=c.Re();
00467 }
00468 
00469 //_____________________________________________________________________________
00470 void TFFTRealComplex::SetPointsComplex(const Double_t *re, const Double_t* /*im*/)
00471 {
00472 //Set all points. Only the real array is used
00473 
00474    for (Int_t i=0; i<fTotalSize; i++)
00475       ((Double_t *)fIn)[i] = re[i];
00476 }
00477 
00478 //_____________________________________________________________________________
00479 UInt_t TFFTRealComplex::MapFlag(Option_t *flag)
00480 {
00481 //allowed options:
00482 //"ES"
00483 //"M"
00484 //"P"
00485 //"EX"
00486 
00487 
00488    TString opt = flag;
00489    opt.ToUpper();
00490    if (opt.Contains("ES"))
00491       return FFTW_ESTIMATE;
00492    if (opt.Contains("M"))
00493       return FFTW_MEASURE;
00494    if (opt.Contains("P"))
00495       return FFTW_PATIENT;
00496    if (opt.Contains("EX"))
00497       return FFTW_EXHAUSTIVE;
00498    return FFTW_ESTIMATE;
00499 }

Generated on Tue Jul 5 14:30:44 2011 for ROOT_528-00b_version by  doxygen 1.5.1