TFFTComplex.cxx

Go to the documentation of this file.
00001 // @(#)root/fft:$Id: TFFTComplex.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 // TFFTComplex                                                           
00015 // One of the interface classes to the FFTW package, can be used directly
00016 // or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
00017 // Computes complex input/output discrete Fourier transforms (DFT) 
00018 // in one or more dimensions. For the detailed information on the computed
00019 // transforms please refer to the FFTW manual, chapter "What FFTW really computes".
00020 // 
00021 // How to use it:
00022 // 1) Create an instance of TFFTComplex - this will allocate input and output
00023 //    arrays (unless an in-place transform is specified)
00024 // 2) Run the Init() function with the desired flags and settings
00025 // 3) Set the data (via SetPoints(), SetPoint() or SetPointComplex() functions)
00026 // 4) Run the Transform() function
00027 // 5) Get the output (via GetPoints(), GetPoint() or GetPointComplex() functions)
00028 // 6) Repeat steps 3)-5) as needed
00029 // 
00030 // For a transform of the same size, but with different flags or sign, rerun the Init()
00031 // function and continue with steps 3)-5)
00032 // NOTE: 1) running Init() function will overwrite the input array! Don't set any data
00033 //          before running the Init() function
00034 //       2) FFTW computes unnormalized transform, so doing a transform followed by 
00035 //          its inverse will lead to the original array scaled by the transform size
00036 //                                                                     
00037 //////////////////////////////////////////////////////////////////////////
00038 
00039 #include "TFFTComplex.h"
00040 #include "fftw3.h"
00041 #include "TComplex.h"
00042 
00043 
00044 ClassImp(TFFTComplex)
00045 
00046 //_____________________________________________________________________________
00047 TFFTComplex::TFFTComplex()
00048 {
00049 //default
00050 
00051    fIn   = 0;
00052    fOut  = 0;
00053    fPlan = 0;
00054    fN    = 0;
00055    fFlags = 0; 
00056    fNdim = 0;
00057    fTotalSize = 0; 
00058    fSign = 1;
00059 }
00060 
00061 //_____________________________________________________________________________
00062 TFFTComplex::TFFTComplex(Int_t n, Bool_t inPlace)
00063 {
00064 //For 1d transforms
00065 //Allocates memory for the input array, and, if inPlace = kFALSE, for the output array
00066 
00067    fIn = fftw_malloc(sizeof(fftw_complex) *n);
00068    if (!inPlace)
00069       fOut = fftw_malloc(sizeof(fftw_complex) * n);
00070    else
00071       fOut = 0;
00072    fN    = new Int_t[1];
00073    fN[0] = n;
00074    fTotalSize = n;
00075    fNdim = 1;
00076    fSign = 1;
00077    fPlan = 0;
00078    fFlags = 0; 
00079 }
00080 
00081 //_____________________________________________________________________________
00082 TFFTComplex::TFFTComplex(Int_t ndim, Int_t *n, Bool_t inPlace)
00083 {
00084 //For multidim. transforms
00085 //Allocates memory for the input array, and, if inPlace = kFALSE, for the output array
00086 
00087    fNdim = ndim;
00088    fTotalSize = 1;
00089    fN = new Int_t[fNdim];
00090    for (Int_t i=0; i<fNdim; i++){
00091       fN[i] = n[i];
00092       fTotalSize*=n[i];
00093    }
00094    fIn = fftw_malloc(sizeof(fftw_complex)*fTotalSize);
00095    if (!inPlace)
00096       fOut = fftw_malloc(sizeof(fftw_complex) * fTotalSize);
00097    else
00098       fOut = 0;
00099    fSign = 1;
00100    fPlan = 0;
00101    fFlags = 0; 
00102 }
00103 
00104 //_____________________________________________________________________________
00105 TFFTComplex::~TFFTComplex()
00106 {
00107 //Destroys the data arrays and the plan. However, some plan information stays around
00108 //until the root session is over, and is reused if other plans of the same size are
00109 //created
00110 
00111    fftw_destroy_plan((fftw_plan)fPlan);
00112    fPlan = 0;
00113    fftw_free((fftw_complex*)fIn);
00114    if (fOut)
00115       fftw_free((fftw_complex*)fOut);
00116    if (fN)
00117       delete [] fN;
00118 }
00119 
00120 //_____________________________________________________________________________
00121 void TFFTComplex::Init( Option_t *flags, Int_t sign,const Int_t* /*kind*/)
00122 {
00123 //Creates the fftw-plan
00124 //
00125 //NOTE:  input and output arrays are overwritten during initialisation,
00126 //       so don't set any points, before running this function!!!!!
00127 //
00128 //2nd parameter: +1
00129 //Argument kind is dummy and doesn't need to be specified
00130 //Possible flag_options:
00131 //"ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
00132 //   performance
00133 //"M" (from "measure") - some time spend in finding the optimal way to do the transform
00134 //"P" (from "patient") - more time spend in finding the optimal way to do the transform
00135 //"EX" (from "exhaustive") - the most optimal way is found
00136 //This option should be chosen depending on how many transforms of the same size and
00137 //type are going to be done. Planning is only done once, for the first transform of this
00138 //size and type.
00139 
00140    fSign = sign;
00141    fFlags = flags;
00142    if (fOut)
00143       fPlan = (void*)fftw_plan_dft(fNdim, fN, (fftw_complex*)fIn, (fftw_complex*)fOut, sign,MapFlag(flags));
00144    else
00145       fPlan = (void*)fftw_plan_dft(fNdim, fN, (fftw_complex*)fIn, (fftw_complex*)fIn, sign, MapFlag(flags));
00146 }
00147 
00148 //_____________________________________________________________________________
00149 void TFFTComplex::Transform()
00150 {
00151 //Computes the transform, specified in Init() function
00152 
00153    if (fPlan)
00154       fftw_execute((fftw_plan)fPlan);
00155    else {
00156       Error("Transform", "transform not initialised");
00157       return;
00158    }
00159 }
00160 
00161 //_____________________________________________________________________________
00162 void TFFTComplex::GetPoints(Double_t *data, Bool_t fromInput) const
00163 {
00164 //Copies the output(or input) into the argument array
00165 
00166    if (!fromInput){
00167       for (Int_t i=0; i<2*fTotalSize; i+=2){
00168          data[i] = ((fftw_complex*)fOut)[i/2][0];
00169          data[i+1] = ((fftw_complex*)fOut)[i/2][1];
00170       }
00171    } else {
00172       for (Int_t i=0; i<2*fTotalSize; i+=2){
00173          data[i] = ((fftw_complex*)fIn)[i/2][0];
00174          data[i+1] = ((fftw_complex*)fIn)[i/2][1];
00175       }
00176    }
00177 }
00178 
00179 //_____________________________________________________________________________
00180 void TFFTComplex::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00181 {
00182 //returns real and imaginary parts of the point #ipoint
00183 
00184    if (fOut && !fromInput){
00185       re = ((fftw_complex*)fOut)[ipoint][0];
00186       im = ((fftw_complex*)fOut)[ipoint][1];
00187    } else {
00188       re = ((fftw_complex*)fIn)[ipoint][0];
00189       im = ((fftw_complex*)fIn)[ipoint][1];
00190    }
00191 }
00192 
00193 //_____________________________________________________________________________
00194 void TFFTComplex::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00195 {
00196 //For multidimensional transforms. Returns real and imaginary parts of the point #ipoint
00197 
00198    Int_t ireal = ipoint[0];
00199    for (Int_t i=0; i<fNdim-1; i++)
00200       ireal=fN[i+1]*ireal + ipoint[i+1];
00201 
00202    if (fOut && !fromInput){
00203       re = ((fftw_complex*)fOut)[ireal][0];
00204       im = ((fftw_complex*)fOut)[ireal][1];
00205    } else {
00206       re = ((fftw_complex*)fIn)[ireal][0];
00207       im = ((fftw_complex*)fIn)[ireal][1];
00208    }
00209 }
00210 
00211 //_____________________________________________________________________________
00212 void TFFTComplex::GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput) const
00213 {
00214 //Copies real and imaginary parts of the output (input) into the argument arrays
00215 
00216    if (fOut && !fromInput){
00217       for (Int_t i=0; i<fTotalSize; i++){
00218          re[i] = ((fftw_complex*)fOut)[i][0];
00219          im[i] = ((fftw_complex*)fOut)[i][1];
00220       }
00221    } else {
00222       for (Int_t i=0; i<fTotalSize; i++){
00223          re[i] = ((fftw_complex*)fIn)[i][0];
00224          im[i] = ((fftw_complex*)fIn)[i][1];
00225       }
00226    }
00227 }
00228 
00229 //_____________________________________________________________________________
00230 void TFFTComplex::GetPointsComplex(Double_t *data, Bool_t fromInput) const
00231 {
00232 //Copies the output(input) into the argument array
00233 
00234    if (fOut && !fromInput){
00235       for (Int_t i=0; i<fTotalSize; i+=2){
00236          data[i] = ((fftw_complex*)fOut)[i/2][0];
00237          data[i+1] = ((fftw_complex*)fOut)[i/2][1];
00238       }
00239    } else {
00240       for (Int_t i=0; i<fTotalSize; i+=2){
00241          data[i] = ((fftw_complex*)fIn)[i/2][0];
00242          data[i+1] = ((fftw_complex*)fIn)[i/2][1];
00243       }
00244    }
00245 }
00246 
00247 //_____________________________________________________________________________
00248 void TFFTComplex::SetPoint(Int_t ipoint, Double_t re, Double_t im)
00249 {
00250 //sets real and imaginary parts of point # ipoint
00251 
00252    ((fftw_complex*)fIn)[ipoint][0]=re;
00253    ((fftw_complex*)fIn)[ipoint][1]=im;
00254 }
00255 
00256 //_____________________________________________________________________________
00257 void TFFTComplex::SetPoint(const Int_t *ipoint, Double_t re, Double_t im)
00258 {
00259 //For multidim. transforms. Sets real and imaginary parts of point # ipoint
00260 
00261    Int_t ireal = ipoint[0];
00262    for (Int_t i=0; i<fNdim-1; i++)
00263       ireal=fN[i+1]*ireal + ipoint[i+1];
00264 
00265    ((fftw_complex*)fIn)[ireal][0]=re;
00266    ((fftw_complex*)fIn)[ireal][1]=im;
00267 }
00268 
00269 //_____________________________________________________________________________
00270 void TFFTComplex::SetPointComplex(Int_t ipoint, TComplex &c)
00271 {
00272    ((fftw_complex*)fIn)[ipoint][0] = c.Re();
00273    ((fftw_complex*)fIn)[ipoint][1] = c.Im();
00274 }
00275 
00276 //_____________________________________________________________________________
00277 void TFFTComplex::SetPoints(const Double_t *data)
00278 {
00279 //set all points. the values are copied. points should be ordered as follows:
00280 //[re_0, im_0, re_1, im_1, ..., re_n, im_n)
00281 
00282    for (Int_t i=0; i<2*fTotalSize-1; i+=2){
00283       ((fftw_complex*)fIn)[i/2][0]=data[i];
00284       ((fftw_complex*)fIn)[i/2][1]=data[i+1];
00285    }
00286 }
00287 
00288 //_____________________________________________________________________________
00289 void TFFTComplex::SetPointsComplex(const Double_t *re_data, const Double_t *im_data)
00290 {
00291 //set all points. the values are copied
00292 
00293    if (!fIn){
00294       Error("SetPointsComplex", "Size is not set yet");
00295       return;
00296    }
00297    for (Int_t i=0; i<fTotalSize; i++){
00298       ((fftw_complex*)fIn)[i][0]=re_data[i];
00299       ((fftw_complex*)fIn)[i][1]=im_data[i];
00300    }
00301 }
00302 
00303 //_____________________________________________________________________________
00304 UInt_t TFFTComplex::MapFlag(Option_t *flag)
00305 {
00306 //allowed options:
00307 //"ES" - FFTW_ESTIMATE
00308 //"M" - FFTW_MEASURE
00309 //"P" - FFTW_PATIENT
00310 //"EX" - FFTW_EXHAUSTIVE
00311 
00312    TString opt = flag;
00313    opt.ToUpper();
00314    if (opt.Contains("ES"))
00315       return FFTW_ESTIMATE;
00316    if (opt.Contains("M"))
00317       return FFTW_MEASURE;
00318    if (opt.Contains("P"))
00319       return FFTW_PATIENT;
00320    if (opt.Contains("EX"))
00321       return FFTW_EXHAUSTIVE;
00322    return FFTW_ESTIMATE;
00323 }

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