TFFTComplexReal.cxx

Go to the documentation of this file.
00001 // @(#)root/fft:$Id: TFFTComplexReal.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 // TFFTComplexReal                                                          
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 the inverse of the real-to-complex transforms (class TFFTRealComplex)
00020 // taking complex input (storing the non-redundant half of a logically Hermitian array)
00021 // to real output (see FFTW manual for more details)
00022 // 
00023 // How to use it:
00024 // 1) Create an instance of TFFTComplexReal - this will allocate input and output
00025 //    arrays (unless an in-place transform is specified)
00026 // 2) Run the Init() function with the desired flags and settings
00027 // 3) Set the data (via SetPoints(), SetPoint() or SetPointComplex() functions)
00028 // 4) Run the Transform() function
00029 // 5) Get the output (via GetPoints(), GetPoint() or GetPointReal() functions)
00030 // 6) Repeat steps 3)-5) as needed
00031 //
00032 // For a transform of the same size, but with different flags, rerun the Init()
00033 // function and continue with steps 3)-5)
00034 // NOTE: 1) running Init() function will overwrite the input array! Don't set any data
00035 //          before running the Init() function
00036 //       2) FFTW computes unnormalized transform, so doing a transform followed by 
00037 //          its inverse will lead to the original array scaled by the transform size
00038 //
00039 //       3) In Complex to Real transform the input array is destroyed. It cannot then 
00040 //          be retrieved when using the Get's methods. 
00041 //                                                                      
00042 //////////////////////////////////////////////////////////////////////////
00043 
00044 #include "TFFTComplexReal.h"
00045 #include "fftw3.h"
00046 #include "TComplex.h"
00047 
00048 
00049 ClassImp(TFFTComplexReal)
00050 
00051 //_____________________________________________________________________________
00052 TFFTComplexReal::TFFTComplexReal()
00053 {
00054 //default
00055 
00056    fIn   = 0;
00057    fOut  = 0;
00058    fPlan = 0;
00059    fN    = 0;
00060    fTotalSize = 0;
00061    fFlags = 0;
00062    fNdim = 0; 
00063 }
00064 
00065 //_____________________________________________________________________________
00066 TFFTComplexReal::TFFTComplexReal(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    if (!inPlace){
00072       fIn = fftw_malloc(sizeof(fftw_complex)*(n/2+1));
00073       fOut = fftw_malloc(sizeof(Double_t)*n);
00074 
00075    } else {
00076       fIn = fftw_malloc(sizeof(fftw_complex)*(n/2+1));
00077       fOut = 0;
00078    }
00079    fNdim = 1;
00080    fN = new Int_t[1];
00081    fN[0] = n;
00082    fPlan = 0;
00083    fTotalSize = n;
00084    fFlags = 0;
00085 }
00086 
00087 //_____________________________________________________________________________
00088 TFFTComplexReal::TFFTComplexReal(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    fNdim = ndim;
00094    fTotalSize = 1;
00095    fN = new Int_t[fNdim];
00096    for (Int_t i=0; i<fNdim; i++){
00097       fN[i] = n[i];
00098       fTotalSize*=n[i];
00099    }
00100    Int_t sizein = Int_t(Double_t(fTotalSize)*(n[ndim-1]/2+1)/n[ndim-1]);
00101    if (!inPlace){
00102       fIn = fftw_malloc(sizeof(fftw_complex)*sizein);
00103       fOut = fftw_malloc(sizeof(Double_t)*fTotalSize);
00104    } else {
00105       fIn = fftw_malloc(sizeof(fftw_complex)*sizein);
00106       fOut = 0;
00107    }
00108    fPlan = 0;
00109    fFlags = 0;
00110 }
00111 
00112 
00113 //_____________________________________________________________________________
00114 TFFTComplexReal::~TFFTComplexReal()
00115 {
00116 //Destroys the data arrays and the plan. However, some plan information stays around
00117 //until the root session is over, and is reused if other plans of the same size are
00118 //created
00119 
00120    fftw_destroy_plan((fftw_plan)fPlan);
00121    fPlan = 0;
00122    fftw_free((fftw_complex*)fIn);
00123    if (fOut)
00124       fftw_free(fOut);
00125    fIn = 0;
00126    fOut = 0;
00127    if (fN)
00128       delete [] fN;
00129    fN = 0;
00130 }
00131 
00132 //_____________________________________________________________________________
00133 void TFFTComplexReal::Init( Option_t *flags, Int_t /*sign*/,const Int_t* /*kind*/)
00134 {
00135 //Creates the fftw-plan
00136 //
00137 //NOTE:  input and output arrays are overwritten during initialisation,
00138 //       so don't set any points, before running this function!!!!!
00139 //
00140 //Arguments sign and kind are dummy and not need to be specified
00141 //Possible flag_options:
00142 //"ES" (from "estimate") - no time in preparing the transform, but probably sub-optimal
00143 //   performanc
00144 //"M" (from "measure") - some time spend in finding the optimal way to do the transform
00145 //"P" (from "patient") - more time spend in finding the optimal way to do the transform
00146 //"EX" (from "exhaustive") - the most optimal way is found
00147 //This option should be chosen depending on how many transforms of the same size and
00148 //type are going to be done. Planning is only done once, for the first transform of this
00149 //size and type.
00150 
00151    fFlags = flags;
00152 
00153    if (fOut)
00154       fPlan = (void*)fftw_plan_dft_c2r(fNdim, fN,(fftw_complex*)fIn,(Double_t*)fOut, MapFlag(flags));
00155    else
00156       fPlan = (void*)fftw_plan_dft_c2r(fNdim, fN, (fftw_complex*)fIn, (Double_t*)fIn, MapFlag(flags));
00157 }
00158 
00159 //_____________________________________________________________________________
00160 void TFFTComplexReal::Transform()
00161 {
00162 //Computes the transform, specified in Init() function
00163 
00164    if (fPlan)
00165       fftw_execute((fftw_plan)fPlan);
00166    else {
00167       Error("Transform", "transform was not initialized");
00168       return;
00169    }
00170 }
00171 
00172 //_____________________________________________________________________________
00173 void TFFTComplexReal::GetPoints(Double_t *data, Bool_t fromInput) const
00174 {
00175 //Fills the argument array with the computed transform
00176 // Works only for output (input array is destroyed in a C2R transform)
00177 
00178    if (fromInput){
00179       Error("GetPoints", "Input array has been destroyed");
00180       return;
00181    }
00182    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00183    std::copy(array,array+fTotalSize, data);
00184 }
00185 
00186 //_____________________________________________________________________________
00187 Double_t TFFTComplexReal::GetPointReal(Int_t ipoint, Bool_t fromInput) const
00188 {
00189 //Returns the point #ipoint
00190 // Works only for output (input array is destroyed in a C2R transform)
00191 
00192    if (fromInput){
00193       Error("GetPointReal", "Input array has been destroyed");
00194       return 0;
00195    }
00196    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00197    return array[ipoint];
00198 }
00199 
00200 //_____________________________________________________________________________
00201 Double_t TFFTComplexReal::GetPointReal(const Int_t *ipoint, Bool_t fromInput) const
00202 {
00203 //For multidimensional transforms. Returns the point #ipoint
00204 // Works only for output (input array is destroyed in a C2R transform)
00205 
00206    if (fromInput){
00207       Error("GetPointReal", "Input array has been destroyed");
00208       return 0;
00209    }
00210    Int_t ireal = ipoint[0];
00211    for (Int_t i=0; i<fNdim-1; i++)
00212       ireal=fN[i+1]*ireal + ipoint[i+1];
00213 
00214    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00215    return array[ireal];
00216 }
00217 
00218 
00219 //_____________________________________________________________________________
00220 void TFFTComplexReal::GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00221 {
00222 // Works only for output (input array is destroyed in a C2R transform)
00223 
00224    if (fromInput){
00225       Error("GetPointComplex", "Input array has been destroyed");
00226       return;
00227    }
00228    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00229    re = array[ipoint];
00230    im = 0;         
00231 }
00232 
00233 //_____________________________________________________________________________
00234 void TFFTComplexReal::GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput) const
00235 {
00236 //For multidimensional transforms. Returns the point #ipoint
00237 // Works only for output (input array is destroyed in a C2R transform)
00238 
00239    if (fromInput){
00240       Error("GetPointComplex", "Input array has been destroyed");
00241       return;
00242    }
00243    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00244 
00245    Int_t ireal = ipoint[0];
00246    for (Int_t i=0; i<fNdim-1; i++)
00247       ireal=fN[i+1]*ireal + ipoint[i+1];
00248 
00249    re = array[ireal];
00250    im = 0;
00251 }
00252 //_____________________________________________________________________________
00253 Double_t* TFFTComplexReal::GetPointsReal(Bool_t fromInput) const
00254 {
00255 //Returns the array of computed transform
00256 // Works only for output (input array is destroyed in a C2R transform)
00257    
00258    // we have 2 different cases
00259    // fromInput = false; fOut = !NULL (transformed is not in place) : return fOut 
00260    // fromInput = false; fOut = NULL (transformed is in place) : return fIn
00261 
00262    if (fromInput) { 
00263       Error("GetPointsReal","Input array was destroyed"); 
00264       return 0; 
00265    }
00266    return  (Double_t*) ( (fOut) ?  fOut :  fIn );  
00267 }
00268 
00269 //_____________________________________________________________________________
00270 void TFFTComplexReal::GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput) const
00271 {
00272 //Fills the argument array with the computed transform
00273 // Works only for output (input array is destroyed in a C2R transform)
00274 
00275    if (fromInput){
00276       Error("GetPointsComplex", "Input array has been destroyed");
00277       return;
00278    }
00279    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00280    for (Int_t i=0; i<fTotalSize; i++){
00281       re[i] = array[i];
00282       im[i] = 0;
00283    }
00284 }
00285 
00286 //_____________________________________________________________________________
00287 void TFFTComplexReal::GetPointsComplex(Double_t *data, Bool_t fromInput) const
00288 {
00289 //Fills the argument array with the computed transform. 
00290 // Works only for output (input array is destroyed in a C2R transform)
00291 
00292    if (fromInput){
00293       Error("GetPointsComplex", "Input array has been destroyed");
00294       return;
00295    }
00296    const Double_t * array =  (const Double_t*) ( (fOut) ?  fOut :  fIn );  
00297    for (Int_t i=0; i<fTotalSize; i+=2){
00298       data[i] = array[i/2];
00299       data[i+1]=0;
00300    }
00301 }
00302 
00303 //_____________________________________________________________________________
00304 void TFFTComplexReal::SetPoint(Int_t ipoint, Double_t re, Double_t im)
00305 {
00306 //since the input must be complex-Hermitian, if the ipoint > n/2, the according
00307 //point before n/2 is set to (re, -im)
00308 
00309    if (ipoint <= fN[0]/2){
00310       ((fftw_complex*)fIn)[ipoint][0] = re;
00311       ((fftw_complex*)fIn)[ipoint][1] = im;
00312    } else {
00313       ((fftw_complex*)fOut)[2*(fN[0]/2)-ipoint][0] = re;
00314       ((fftw_complex*)fOut)[2*(fN[0]/2)-ipoint][1] = -im;
00315    }
00316 }
00317 
00318 //_____________________________________________________________________________
00319 void TFFTComplexReal::SetPoint(const Int_t *ipoint, Double_t re, Double_t im)
00320 {
00321 //Set the point #ipoint. Since the input is Hermitian, only the first (roughly)half of
00322 //the points have to be set.
00323 
00324    Int_t ireal = ipoint[0];
00325    for (Int_t i=0; i<fNdim-2; i++){
00326       ireal=fN[i+1]*ireal + ipoint[i+1];
00327    } 
00328    ireal = (fN[fNdim-1]/2+1)*ireal+ipoint[fNdim-1];
00329    Int_t realN = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00330 
00331    if (ireal > realN){
00332       Error("SetPoint", "Illegal index value");
00333       return;
00334    }
00335    ((fftw_complex*)fIn)[ireal][0] = re;
00336    ((fftw_complex*)fIn)[ireal][1] = im;
00337 }
00338 
00339 //_____________________________________________________________________________
00340 void TFFTComplexReal::SetPointComplex(Int_t ipoint, TComplex &c)
00341 {
00342 //since the input must be complex-Hermitian, if the ipoint > n/2, the according
00343 //point before n/2 is set to (re, -im)
00344 
00345    if (ipoint <= fN[0]/2){
00346       ((fftw_complex*)fIn)[ipoint][0] = c.Re();
00347       ((fftw_complex*)fIn)[ipoint][1] = c.Im();
00348    } else {
00349       ((fftw_complex*)fOut)[2*(fN[0]/2)-ipoint][0] = c.Re();
00350       ((fftw_complex*)fOut)[2*(fN[0]/2)-ipoint][1] = -c.Im();
00351    }
00352 }
00353 
00354 //_____________________________________________________________________________
00355 void TFFTComplexReal::SetPoints(const Double_t *data)
00356 {
00357 //set all points. the values are copied. points should be ordered as follows:
00358 //[re_0, im_0, re_1, im_1, ..., re_n, im_n)
00359 
00360    Int_t sizein = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00361 
00362    for (Int_t i=0; i<2*(sizein); i+=2){
00363       ((fftw_complex*)fIn)[i/2][0]=data[i];
00364       ((fftw_complex*)fIn)[i/2][1]=data[i+1];
00365    }
00366 }
00367 
00368 //_____________________________________________________________________________
00369 void TFFTComplexReal::SetPointsComplex(const Double_t *re, const Double_t *im)
00370 {
00371 //Set all points. The values are copied.
00372 
00373    Int_t sizein = Int_t(Double_t(fTotalSize)*(fN[fNdim-1]/2+1)/fN[fNdim-1]);
00374    for (Int_t i=0; i<sizein; i++){
00375       ((fftw_complex*)fIn)[i][0]=re[i];
00376       ((fftw_complex*)fIn)[i][1]=im[i];
00377    }
00378 }
00379 
00380 //_____________________________________________________________________________
00381 UInt_t TFFTComplexReal::MapFlag(Option_t *flag)
00382 {
00383 //allowed options:
00384 //"ES" - FFTW_ESTIMATE
00385 //"M" - FFTW_MEASURE
00386 //"P" - FFTW_PATIENT
00387 //"EX" - FFTW_EXHAUSTIVE
00388 
00389    TString opt = flag;
00390    opt.ToUpper();
00391    if (opt.Contains("ES"))
00392       return FFTW_ESTIMATE;
00393    if (opt.Contains("M"))
00394       return FFTW_MEASURE;
00395    if (opt.Contains("P"))
00396       return FFTW_PATIENT;
00397    if (opt.Contains("EX"))
00398       return FFTW_EXHAUSTIVE;
00399    return FFTW_ESTIMATE;
00400 }

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