ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TFFTReal.h
Go to the documentation of this file.
1 // @(#)root/fft:$Id$
2 // Author: Anna Kreshuk 07/4/2006
3 
4 /*************************************************************************
5  * Copyright (C) 1995-2006, Rene Brun and Fons Rademakers. *
6  * All rights reserved. *
7  * *
8  * For the licensing terms see $ROOTSYS/LICENSE. *
9  * For the list of contributors see $ROOTSYS/README/CREDITS. *
10  *************************************************************************/
11 
12 #ifndef ROOT_TFFTReal
13 #define ROOT_TFFTReal
14 
15 //////////////////////////////////////////////////////////////////////////
16 //
17 // TFFTReal
18 // One of the interface classes to the FFTW package, can be used directly
19 // or via the TVirtualFFT class. Only the basic interface of FFTW is implemented.
20 //
21 // Computes transforms called r2r in FFTW manual:
22 // - transforms of real input and output in "halfcomplex" format i.e.
23 // real and imaginary parts for a transform of size n stored as
24 // (r0, r1, r2, ..., rn/2, i(n+1)/2-1, ..., i2, i1)
25 // - discrete Hartley transform
26 // - sine and cosine transforms (DCT-I,II,III,IV and DST-I,II,III,IV)
27 // For the detailed information on the computed
28 // transforms please refer to the FFTW manual, chapter "What FFTW really computes".
29 //
30 // How to use it:
31 // 1) Create an instance of TFFTReal - this will allocate input and output
32 // arrays (unless an in-place transform is specified)
33 // 2) Run the Init() function with the desired flags and settings (see function
34 // comments for possible kind parameters)
35 // 3) Set the data (via SetPoints()or SetPoint() functions)
36 // 4) Run the Transform() function
37 // 5) Get the output (via GetPoints() or GetPoint() functions)
38 // 6) Repeat steps 3)-5) as needed
39 // For a transform of the same size, but of different kind (or with different flags),
40 // rerun the Init() function and continue with steps 3)-5)
41 //
42 // NOTE: 1) running Init() function will overwrite the input array! Don't set any data
43 // before running the Init() function!
44 // 2) FFTW computes unnormalized transform, so doing a transform followed by
45 // its inverse will lead to the original array scaled BY:
46 // - transform size (N) for R2HC, HC2R, DHT transforms
47 // - 2*(N-1) for DCT-I (REDFT00)
48 // - 2*(N+1) for DST-I (RODFT00)
49 // - 2*N for the remaining transforms
50 // Transform inverses:
51 // R2HC<-->HC2R
52 // DHT<-->DHT
53 // DCT-I<-->DCT-I
54 // DCT-II<-->DCT-III
55 // DCT-IV<-->DCT-IV
56 // DST-I<-->DST-I
57 // DST-II<-->DST-III
58 // DST-IV<-->DST-IV
59 //
60 //////////////////////////////////////////////////////////////////////////
61 
62 #ifndef ROOT_TVirtualFFT
63 #include "TVirtualFFT.h"
64 #endif
65 
66 class TComplex;
67 
68 class TFFTReal: public TVirtualFFT{
69  protected:
70  void *fIn; //input array
71  void *fOut; //output array
72  void *fPlan; //fftw plan (the plan how to compute the transform)
73  Int_t fNdim; //number of dimensions
74  Int_t fTotalSize; //total size of the transform
75  Int_t *fN; //transform sizes in each dimension
76  void *fKind; //transform kinds in each dimension
77  Option_t *fFlags; //transform flags
78 
79  Int_t MapOptions(const Int_t *kind);
80  UInt_t MapFlag(Option_t *flag);
81 
82  public:
83  TFFTReal();
84  TFFTReal(Int_t n, Bool_t inPlace=kFALSE);
85  TFFTReal(Int_t ndim, Int_t *n, Bool_t inPlace=kFALSE);
86  virtual ~TFFTReal();
87 
88  virtual void Init( Option_t *flags,Int_t sign, const Int_t *kind);
89 
90  virtual Int_t GetSize() const {return fTotalSize;}
91  virtual Int_t *GetN() const {return fN;}
92  virtual Int_t GetNdim() const {return fNdim;}
93  virtual Option_t *GetType() const;
94  virtual Int_t GetSign() const {return 0;}
95  virtual Option_t *GetTransformFlag() const {return fFlags;}
96  virtual Bool_t IsInplace() const {if (fOut) return kTRUE; else return kFALSE;}
97 
98  virtual void GetPoints(Double_t *data, Bool_t fromInput = kFALSE) const;
99  virtual Double_t GetPointReal(Int_t ipoint, Bool_t fromInput = kFALSE) const;
100  virtual Double_t GetPointReal(const Int_t *ipoint, Bool_t fromInput = kFALSE) const;
101  virtual void GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const;
102 
103  virtual void GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const;
104 
105  virtual Double_t *GetPointsReal(Bool_t fromInput=kFALSE) const;
106  virtual void GetPointsComplex(Double_t* /*re*/, Double_t* /*im*/, Bool_t /*fromInput = kFALSE*/) const{};
107  virtual void GetPointsComplex(Double_t* /*data*/, Bool_t /*fromInput = kFALSE*/) const {};
108 
109  virtual void SetPoint(Int_t ipoint, Double_t re, Double_t im = 0);
110  virtual void SetPoint(const Int_t *ipoint, Double_t re, Double_t /*im=0*/);
111  virtual void SetPoints(const Double_t *data);
112  virtual void SetPointComplex(Int_t /*ipoint*/, TComplex &/*c*/){};
113  virtual void SetPointsComplex(const Double_t* /*re*/, const Double_t* /*im*/){};
114  virtual void Transform();
115 
116 
117  ClassDef(TFFTReal,0);
118 };
119 
120 #endif
virtual void SetPoints(const Double_t *data)
Sets all points.
Definition: TFFTReal.cxx:342
virtual void GetPointComplex(const Int_t *ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const
Only for input of HC2R and output of R2HC and for 1d.
Definition: TFFTReal.cxx:277
virtual void GetPoints(Double_t *data, Bool_t fromInput=kFALSE) const
Copies the output (or input) points into the provided array, that should be big enough.
Definition: TFFTReal.cxx:221
const char Option_t
Definition: RtypesCore.h:62
Int_t * fN
Definition: TFFTReal.h:75
void * fOut
Definition: TFFTReal.h:71
virtual Double_t * GetPointsReal(Bool_t fromInput=kFALSE) const
Returns the output (or input) array.
Definition: TFFTReal.cxx:285
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual Option_t * GetType() const
Returns the type of the transform.
Definition: TFFTReal.cxx:205
virtual Int_t GetSize() const
Definition: TFFTReal.h:90
virtual Int_t GetSign() const
Definition: TFFTReal.h:94
void * fIn
Definition: TFFTReal.h:70
virtual Int_t GetNdim() const
Definition: TFFTReal.h:92
Int_t MapOptions(const Int_t *kind)
transfers the r2r_kind parameters to fftw type
Definition: TFFTReal.cxx:351
UInt_t MapFlag(Option_t *flag)
allowed options: "ES" - FFTW_ESTIMATE "M" - FFTW_MEASURE "P" - FFTW_PATIENT "EX" - FFTW_EXHAUSTIVE ...
Definition: TFFTReal.cxx:397
virtual Option_t * GetTransformFlag() const
Definition: TFFTReal.h:95
virtual void Init(Option_t *flags, Int_t sign, const Int_t *kind)
Creates the fftw-plan.
Definition: TFFTReal.cxx:171
virtual void SetPointComplex(Int_t, TComplex &)
Definition: TFFTReal.h:112
virtual void GetPointsComplex(Double_t *, Bool_t) const
Definition: TFFTReal.h:107
virtual Double_t GetPointReal(Int_t ipoint, Bool_t fromInput=kFALSE) const
For 1d tranforms. Returns point #ipoint.
Definition: TFFTReal.cxx:231
void * fKind
Definition: TFFTReal.h:76
unsigned int UInt_t
Definition: RtypesCore.h:42
virtual void SetPoint(Int_t ipoint, Double_t re, Double_t im=0)
Definition: TFFTReal.cxx:305
TVirtualFFT is an interface class for Fast Fourier Transforms.
Definition: TVirtualFFT.h:92
Int_t fNdim
Definition: TFFTReal.h:73
virtual ~TFFTReal()
clean-up
Definition: TFFTReal.cxx:125
virtual Int_t * GetN() const
Definition: TFFTReal.h:91
double Double_t
Definition: RtypesCore.h:55
virtual void GetPointsComplex(Double_t *, Double_t *, Bool_t) const
Definition: TFFTReal.h:106
virtual void SetPointsComplex(const Double_t *, const Double_t *)
Definition: TFFTReal.h:113
virtual Bool_t IsInplace() const
Definition: TFFTReal.h:96
Int_t fTotalSize
Definition: TFFTReal.h:74
ClassDef(TFFTReal, 0)
Option_t * fFlags
Definition: TFFTReal.h:77
Int_t sign(Double_t x)
Definition: CsgOps.cxx:89
const Bool_t kTRUE
Definition: Rtypes.h:91
virtual void Transform()
Computes the transform, specified in Init() function.
Definition: TFFTReal.cxx:192
const Int_t n
Definition: legend1.C:16
void * fPlan
Definition: TFFTReal.h:72