TFFTReal.cxx

Go to the documentation of this file.
00001 // @(#)root/fft:$Id: TFFTReal.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 // TFFTReal                                                       
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 //
00018 // Computes transforms called r2r in FFTW manual: 
00019 // - transforms of real input and output in "halfcomplex" format i.e. 
00020 //   real and imaginary parts for a transform of size n stored as 
00021 //   (r0, r1, r2, ..., rn/2, i(n+1)/2-1, ..., i2, i1)
00022 // - discrete Hartley transform
00023 // - sine and cosine transforms (DCT-I,II,III,IV and DST-I,II,III,IV)
00024 // For the detailed information on the computed
00025 // transforms please refer to the FFTW manual, chapter "What FFTW really computes".
00026 //
00027 // How to use it:
00028 // 1) Create an instance of TFFTReal - this will allocate input and output
00029 //    arrays (unless an in-place transform is specified)
00030 // 2) Run the Init() function with the desired flags and settings (see function
00031 //    comments for possible kind parameters)
00032 // 3) Set the data (via SetPoints()or SetPoint() functions)
00033 // 4) Run the Transform() function
00034 // 5) Get the output (via GetPoints() or GetPoint() functions)
00035 // 6) Repeat steps 3)-5) as needed
00036 // For a transform of the same size, but of different kind (or with different flags), 
00037 // rerun the Init() function and continue with steps 3)-5)
00038 //
00039 // NOTE: 1) running Init() function will overwrite the input array! Don't set any data
00040 //          before running the Init() function!
00041 //       2) FFTW computes unnormalized transform, so doing a transform followed by 
00042 //          its inverse will lead to the original array scaled BY:
00043 //          - transform size (N) for R2HC, HC2R, DHT transforms
00044 //          - 2*(N-1) for DCT-I (REDFT00)
00045 //          - 2*(N+1) for DST-I (RODFT00)
00046 //          - 2*N for the remaining transforms
00047 // Transform inverses:
00048 // R2HC<-->HC2R
00049 // DHT<-->DHT
00050 // DCT-I<-->DCT-I
00051 // DCT-II<-->DCT-III
00052 // DCT-IV<-->DCT-IV
00053 // DST-I<-->DST-I
00054 // DST-II<-->DST-III
00055 // DST-IV<-->DST-IV
00056 // 
00057 //////////////////////////////////////////////////////////////////////////
00058 
00059 #include "TFFTReal.h"
00060 #include "fftw3.h"
00061 
00062 ClassImp(TFFTReal)
00063 
00064 //_____________________________________________________________________________
00065 TFFTReal::TFFTReal()
00066 {
00067 //default
00068 
00069    fIn    = 0;
00070    fOut   = 0;
00071    fPlan  = 0;
00072    fN     = 0;
00073    fKind  = 0;
00074    fFlags = 0;
00075    fNdim = 0; 
00076    fTotalSize = 0;
00077 }
00078 
00079 //_____________________________________________________________________________
00080 TFFTReal::TFFTReal(Int_t n, Bool_t inPlace)
00081 {
00082 //For 1d transforms
00083 //n here is the physical size of the transform (see FFTW manual for more details)
00084 
00085    fIn = fftw_malloc(sizeof(Double_t)*n);
00086    if (inPlace) fOut = 0;
00087    else fOut = fftw_malloc(sizeof(Double_t)*n);
00088 
00089    fPlan = 0;
00090    fNdim = 1;
00091    fN = new Int_t[1];
00092    fN[0] = n;
00093    fKind = 0;
00094    fTotalSize = n;
00095    fFlags = 0;
00096 }
00097 
00098 //_____________________________________________________________________________
00099 TFFTReal::TFFTReal(Int_t ndim, Int_t *n, Bool_t inPlace)
00100 {
00101 //For multidimensional transforms
00102 //1st parameter is the # of dimensions,
00103 //2nd is the sizes (physical) of the transform in each dimension
00104 
00105    fTotalSize = 1;
00106    fNdim = ndim;
00107    fN = new Int_t[ndim];
00108    fKind = 0;
00109    fPlan = 0;
00110    fFlags = 0; 
00111    for (Int_t i=0; i<ndim; i++){
00112       fTotalSize*=n[i];
00113       fN[i] = n[i];
00114    }
00115    fIn = fftw_malloc(sizeof(Double_t)*fTotalSize);
00116    if (!inPlace)
00117       fOut = fftw_malloc(sizeof(Double_t)*fTotalSize);
00118    else
00119       fOut = 0;
00120 }
00121 
00122 //_____________________________________________________________________________
00123 TFFTReal::~TFFTReal()
00124 {
00125 //clean-up
00126 
00127    fftw_destroy_plan((fftw_plan)fPlan);
00128    fPlan = 0;
00129    fftw_free(fIn);
00130    fIn = 0;
00131    if (fOut){
00132       fftw_free(fOut);
00133    }
00134    fOut = 0;
00135    if (fN)
00136       delete [] fN;
00137    fN = 0;
00138    if (fKind)
00139       fftw_free((fftw_r2r_kind*)fKind);
00140    fKind = 0;
00141 }
00142 
00143 //_____________________________________________________________________________
00144 void TFFTReal::Init( Option_t* flags,Int_t /*sign*/, const Int_t *kind)
00145 {
00146 //Creates the fftw-plan
00147 //
00148 //NOTE:  input and output arrays are overwritten during initialisation,
00149 //       so don't set any points, before running this function!!!!!
00150 //
00151 //1st parameter:
00152 //  Possible flag_options:
00153 //  "ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
00154 //       performance
00155 //  "M" (from "measure") - some time spend in finding the optimal way to do the transform
00156 //  "P" (from "patient") - more time spend in finding the optimal way to do the transform
00157 //  "EX" (from "exhaustive") - the most optimal way is found
00158 //  This option should be chosen depending on how many transforms of the same size and
00159 //  type are going to be done. Planning is only done once, for the first transform of this
00160 //  size and type.
00161 //2nd parameter is dummy and doesn't need to be specified
00162 //3rd parameter- transform kind for each dimension
00163 //     4 different kinds of sine and cosine transforms are available
00164 //     DCT-I   - kind=0
00165 //     DCT-II  - kind=1
00166 //     DCT-III - kind=2
00167 //     DCT-IV  - kind=3
00168 //     DST-I   - kind=4
00169 //     DST-II  - kind=5
00170 //     DSTIII  - kind=6
00171 //     DSTIV   - kind=7
00172 
00173    if (!fKind)
00174       fKind = (fftw_r2r_kind*)fftw_malloc(sizeof(fftw_r2r_kind)*fNdim);
00175 
00176    if (MapOptions(kind)){
00177       if (fOut)
00178          fPlan = (void*)fftw_plan_r2r(fNdim, fN, (Double_t*)fIn, (Double_t*)fOut, (fftw_r2r_kind*)fKind, MapFlag(flags));
00179       else
00180          fPlan = (void*)fftw_plan_r2r(fNdim, fN, (Double_t*)fIn, (Double_t*)fIn, (fftw_r2r_kind*)fKind, MapFlag(flags));
00181       fFlags = flags;
00182    }
00183 }
00184 
00185 //_____________________________________________________________________________
00186 void TFFTReal::Transform()
00187 {
00188 //Computes the transform, specified in Init() function
00189 
00190    if (fPlan)
00191       fftw_execute((fftw_plan)fPlan);
00192    else {
00193       Error("Transform", "transform hasn't been initialised");
00194       return;
00195    }
00196 }
00197 
00198 //_____________________________________________________________________________
00199 Option_t *TFFTReal::GetType() const
00200 {
00201 //Returns the type of the transform
00202 
00203    if (!fKind) {
00204       Error("GetType", "Type not defined yet (kind not set)");
00205       return "";
00206    }
00207    if (((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC) return "R2HC";
00208    if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R) return "HC2R";
00209    if (((fftw_r2r_kind*)fKind)[0]==FFTW_DHT) return "DHT";
00210    else return "R2R";
00211 }
00212 
00213 //_____________________________________________________________________________
00214 void TFFTReal::GetPoints(Double_t *data, Bool_t fromInput) const
00215 {
00216 //Copies the output (or input) points into the provided array, that should
00217 //be big enough
00218 
00219    const Double_t * array = GetPointsReal(fromInput);
00220    if (!array) return; 
00221    std::copy(array, array+fTotalSize, data);
00222 }
00223 
00224 //_____________________________________________________________________________
00225 Double_t TFFTReal::GetPointReal(Int_t ipoint, Bool_t fromInput) const
00226 {
00227 //For 1d tranforms. Returns point #ipoint
00228 
00229    if (ipoint<0 || ipoint>fTotalSize){
00230       Error("GetPointReal", "No such point");
00231       return 0;
00232    }
00233    const Double_t * array = GetPointsReal(fromInput);
00234    return ( array ) ? array[ipoint] : 0; 
00235 }
00236 
00237 //_____________________________________________________________________________
00238 Double_t TFFTReal::GetPointReal(const Int_t *ipoint, Bool_t fromInput) const
00239 {
00240 //For multidim.transforms. Returns point #ipoint
00241 
00242    Int_t ireal = ipoint[0];
00243    for (Int_t i=0; i<fNdim-1; i++)
00244       ireal=fN[i+1]*ireal + ipoint[i+1];
00245 
00246    const Double_t * array = GetPointsReal(fromInput);
00247    return ( array ) ? array[ireal] : 0; 
00248 }
00249 
00250 //_____________________________________________________________________________
00251 void TFFTReal::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00252 {
00253 //Only for input of HC2R and output of R2HC
00254 
00255    const Double_t * array = GetPointsReal(fromInput);
00256    if (!array) return; 
00257    if ( ( ((fftw_r2r_kind*)fKind)[0]==FFTW_R2HC && !fromInput ) ||
00258         ( ((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R &&  fromInput ) )
00259    {
00260       if (ipoint<fN[0]/2+1){
00261          re = array[ipoint];
00262          im = array[fN[0]-ipoint];
00263       } else {
00264          re = array[fN[0]-ipoint];
00265          im = -array[ipoint];
00266       }
00267       if ((fN[0]%2)==0 && ipoint==fN[0]/2) im = 0;
00268    }
00269 }
00270 //_____________________________________________________________________________
00271 void TFFTReal::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00272 {
00273 //Only for input of HC2R and output of R2HC and for 1d
00274 
00275    GetPointComplex(ipoint[0], re, im, fromInput);
00276 }
00277 
00278 //_____________________________________________________________________________
00279 Double_t* TFFTReal::GetPointsReal(Bool_t fromInput) const
00280 {
00281 //Returns the output (or input) array
00282 
00283    // we have 4 different cases
00284    // fromInput = false; fOut = !NULL (transformed is not in place) : return fOut 
00285    // fromInput = false; fOut = NULL (transformed is in place) : return fIn
00286    // fromInput = true; fOut = !NULL :   return fIn
00287    // fromInput = true; fOut = NULL return an error since input array is overwritten
00288 
00289    if (!fromInput && fOut) 
00290       return (Double_t*)fOut;
00291    else if (fromInput && !fOut) { 
00292       Error("GetPointsReal","Input array was destroyed"); 
00293       return 0; 
00294    }
00295    //R__ASSERT(fIn);
00296    return (Double_t*)fIn;
00297 }
00298 
00299 //_____________________________________________________________________________
00300 void TFFTReal::SetPoint(Int_t ipoint, Double_t re, Double_t im)
00301 {
00302    if (ipoint<0 || ipoint>fTotalSize){
00303       Error("SetPoint", "illegal point index");
00304       return;
00305    }
00306    if (((fftw_r2r_kind*)fKind)[0]==FFTW_HC2R){
00307       if ((fN[0]%2)==0 && ipoint==fN[0]/2)
00308          ((Double_t*)fIn)[ipoint] = re;
00309       else {
00310          ((Double_t*)fIn)[ipoint] = re;
00311          ((Double_t*)fIn)[fN[0]-ipoint]=im;
00312       }
00313    }
00314    else
00315       ((Double_t*)fIn)[ipoint]=re;
00316 }
00317 
00318 //_____________________________________________________________________________
00319 void TFFTReal::SetPoint(const Int_t *ipoint, Double_t re, Double_t /*im*/)
00320 {
00321 //Since multidimensional R2HC and HC2R transforms are not supported,
00322 //third parameter is dummy
00323 
00324    Int_t ireal = ipoint[0];
00325    for (Int_t i=0; i<fNdim-1; i++)
00326       ireal=fN[i+1]*ireal + ipoint[i+1];
00327    if (ireal < 0 || ireal >fTotalSize){
00328       Error("SetPoint", "illegal point index");
00329       return;
00330    }
00331    ((Double_t*)fIn)[ireal]=re;
00332 }
00333 
00334 //_____________________________________________________________________________
00335 void TFFTReal::SetPoints(const Double_t *data)
00336 {
00337 //Sets all points
00338 
00339    for (Int_t i=0; i<fTotalSize; i++)
00340       ((Double_t*)fIn)[i] = data[i];
00341 }
00342 
00343 //_____________________________________________________________________________
00344 Int_t TFFTReal::MapOptions(const Int_t *kind)
00345 {
00346 //transfers the r2r_kind parameters to fftw type
00347 
00348    if (kind[0] == 10){
00349       if (fNdim>1){
00350          Error("Init", "Multidimensional R2HC transforms are not supported, use R2C interface instead");
00351          return 0;
00352       }
00353       ((fftw_r2r_kind*)fKind)[0] = FFTW_R2HC;
00354    }
00355    else if (kind[0] == 11) {
00356       if (fNdim>1){
00357          Error("Init", "Multidimensional HC2R transforms are not supported, use C2R interface instead");
00358          return 0;
00359       }
00360       ((fftw_r2r_kind*)fKind)[0] = FFTW_HC2R;
00361    }
00362    else if (kind[0] == 12) {
00363       for (Int_t i=0; i<fNdim; i++)
00364          ((fftw_r2r_kind*)fKind)[i] = FFTW_DHT;
00365    }
00366    else {
00367       for (Int_t i=0; i<fNdim; i++){
00368          switch (kind[i]) {
00369          case 0: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT00;  break;
00370          case 1: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT01;  break;
00371          case 2: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT10;  break;
00372          case 3: ((fftw_r2r_kind*)fKind)[i] = FFTW_REDFT11;  break;
00373          case 4: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT00;  break;
00374          case 5: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT01;  break;
00375          case 6: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT10;  break;
00376          case 7: ((fftw_r2r_kind*)fKind)[i] = FFTW_RODFT11;  break;
00377          default:
00378          ((fftw_r2r_kind*)fKind)[i] = FFTW_R2HC; break;
00379          }
00380       }
00381    }
00382    return 1;
00383 }
00384 
00385 //_____________________________________________________________________________
00386 UInt_t TFFTReal::MapFlag(Option_t *flag)
00387 {
00388 //allowed options:
00389 //"ES" - FFTW_ESTIMATE
00390 //"M" - FFTW_MEASURE
00391 //"P" - FFTW_PATIENT
00392 //"EX" - FFTW_EXHAUSTIVE
00393 
00394    TString opt = flag;
00395    opt.ToUpper();
00396    if (opt.Contains("ES"))
00397       return FFTW_ESTIMATE;
00398    if (opt.Contains("M"))
00399       return FFTW_MEASURE;
00400    if (opt.Contains("P"))
00401       return FFTW_PATIENT;
00402    if (opt.Contains("EX"))
00403       return FFTW_EXHAUSTIVE;
00404    return FFTW_ESTIMATE;
00405 }

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