ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FFT.C
Go to the documentation of this file.
1 #include "TH1D.h"
2 #include "TVirtualFFT.h"
3 #include "TF1.h"
4 #include "TCanvas.h"
5 #include "TMath.h"
6 
7 void FFT()
8 {
9 
10 //This tutorial illustrates the Fast Fourier Transforms interface in ROOT.
11 //FFT transform types provided in ROOT:
12 // - "C2CFORWARD" - a complex input/output discrete Fourier transform (DFT)
13 // in one or more dimensions, -1 in the exponent
14 // - "C2CBACKWARD"- a complex input/output discrete Fourier transform (DFT)
15 // in one or more dimensions, +1 in the exponent
16 // - "R2C" - a real-input/complex-output discrete Fourier transform (DFT)
17 // in one or more dimensions,
18 // - "C2R" - inverse transforms to "R2C", taking complex input
19 // (storing the non-redundant half of a logically Hermitian array)
20 // to real output
21 // - "R2HC" - a real-input DFT with output in ¡Èhalfcomplex¡É format,
22 // i.e. real and imaginary parts for a transform of size n stored as
23 // r0, r1, r2, ..., rn/2, i(n+1)/2-1, ..., i2, i1
24 // - "HC2R" - computes the reverse of FFTW_R2HC, above
25 // - "DHT" - computes a discrete Hartley transform
26 // Sine/cosine transforms:
27 // DCT-I (REDFT00 in FFTW3 notation)
28 // DCT-II (REDFT10 in FFTW3 notation)
29 // DCT-III(REDFT01 in FFTW3 notation)
30 // DCT-IV (REDFT11 in FFTW3 notation)
31 // DST-I (RODFT00 in FFTW3 notation)
32 // DST-II (RODFT10 in FFTW3 notation)
33 // DST-III(RODFT01 in FFTW3 notation)
34 // DST-IV (RODFT11 in FFTW3 notation)
35 //First part of the tutorial shows how to transform the histograms
36 //Second part shows how to transform the data arrays directly
37 //Authors: Anna Kreshuk and Jens Hoffmann
38 
39 
40 //********* Histograms ********//
41 
42 
43  //prepare the canvas for drawing
44  TCanvas *myc = new TCanvas("myc", "Fast Fourier Transform", 800, 600);
45  myc->SetFillColor(45);
46  TPad *c1_1 = new TPad("c1_1", "c1_1",0.01,0.67,0.49,0.99);
47  TPad *c1_2 = new TPad("c1_2", "c1_2",0.51,0.67,0.99,0.99);
48  TPad *c1_3 = new TPad("c1_3", "c1_3",0.01,0.34,0.49,0.65);
49  TPad *c1_4 = new TPad("c1_4", "c1_4",0.51,0.34,0.99,0.65);
50  TPad *c1_5 = new TPad("c1_5", "c1_5",0.01,0.01,0.49,0.32);
51  TPad *c1_6 = new TPad("c1_6", "c1_6",0.51,0.01,0.99,0.32);
52  c1_1->Draw();
53  c1_2->Draw();
54  c1_3->Draw();
55  c1_4->Draw();
56  c1_5->Draw();
57  c1_6->Draw();
58  c1_1->SetFillColor(30);
59  c1_1->SetFrameFillColor(42);
60  c1_2->SetFillColor(30);
61  c1_2->SetFrameFillColor(42);
62  c1_3->SetFillColor(30);
63  c1_3->SetFrameFillColor(42);
64  c1_4->SetFillColor(30);
65  c1_4->SetFrameFillColor(42);
66  c1_5->SetFillColor(30);
67  c1_5->SetFrameFillColor(42);
68  c1_6->SetFillColor(30);
69  c1_6->SetFrameFillColor(42);
70 
71  c1_1->cd();
73 
74  //A function to sample
75  TF1 *fsin = new TF1("fsin", "sin(x)+sin(2*x)+sin(0.5*x)+1", 0, 4*TMath::Pi());
76  fsin->Draw();
77 
78  Int_t n=25;
79  TH1D *hsin = new TH1D("hsin", "hsin", n+1, 0, 4*TMath::Pi());
80  Double_t x;
81 
82  //Fill the histogram with function values
83  for (Int_t i=0; i<=n; i++){
84  x = (Double_t(i)/n)*(4*TMath::Pi());
85  hsin->SetBinContent(i+1, fsin->Eval(x));
86  }
87  hsin->Draw("same");
88  fsin->GetXaxis()->SetLabelSize(0.05);
89  fsin->GetYaxis()->SetLabelSize(0.05);
90 
91  c1_2->cd();
92  //Compute the transform and look at the magnitude of the output
93  TH1 *hm =0;
95  hm = hsin->FFT(hm, "MAG");
96  hm->SetTitle("Magnitude of the 1st transform");
97  hm->Draw();
98  //NOTE: for "real" frequencies you have to divide the x-axes range with the range of your function
99  //(in this case 4*Pi); y-axes has to be rescaled by a factor of 1/SQRT(n) to be right: this is not done automatically!
100 
101  hm->SetStats(kFALSE);
102  hm->GetXaxis()->SetLabelSize(0.05);
103  hm->GetYaxis()->SetLabelSize(0.05);
104  c1_3->cd();
105  //Look at the phase of the output
106  TH1 *hp = 0;
107  hp = hsin->FFT(hp, "PH");
108  hp->SetTitle("Phase of the 1st transform");
109  hp->Draw();
110  hp->SetStats(kFALSE);
111  hp->GetXaxis()->SetLabelSize(0.05);
112  hp->GetYaxis()->SetLabelSize(0.05);
113 
114  //Look at the DC component and the Nyquist harmonic:
115  Double_t re, im;
116  //That's the way to get the current transform object:
118  c1_4->cd();
119  //Use the following method to get just one point of the output
120  fft->GetPointComplex(0, re, im);
121  printf("1st transform: DC component: %f\n", re);
122  fft->GetPointComplex(n/2+1, re, im);
123  printf("1st transform: Nyquist harmonic: %f\n", re);
124 
125  //Use the following method to get the full output:
126  Double_t *re_full = new Double_t[n];
127  Double_t *im_full = new Double_t[n];
128  fft->GetPointsComplex(re_full,im_full);
129 
130  //Now let's make a backward transform:
131  TVirtualFFT *fft_back = TVirtualFFT::FFT(1, &n, "C2R M K");
132  fft_back->SetPointsComplex(re_full,im_full);
133  fft_back->Transform();
134  TH1 *hb = 0;
135  //Let's look at the output
136  hb = TH1::TransformHisto(fft_back,hb,"Re");
137  hb->SetTitle("The backward transform result");
138  hb->Draw();
139  //NOTE: here you get at the x-axes number of bins and not real values
140  //(in this case 25 bins has to be rescaled to a range between 0 and 4*Pi;
141  //also here the y-axes has to be rescaled (factor 1/bins)
142  hb->SetStats(kFALSE);
143  hb->GetXaxis()->SetLabelSize(0.05);
144  hb->GetYaxis()->SetLabelSize(0.05);
145  delete fft_back;
146  fft_back=0;
147 
148 //********* Data array - same transform ********//
149 
150  //Allocate an array big enough to hold the transform output
151  //Transform output in 1d contains, for a transform of size N,
152  //N/2+1 complex numbers, i.e. 2*(N/2+1) real numbers
153  //our transform is of size n+1, because the histogram has n+1 bins
154 
155  Double_t *in = new Double_t[2*((n+1)/2+1)];
156  Double_t re_2,im_2;
157  for (Int_t i=0; i<=n; i++){
158  x = (Double_t(i)/n)*(4*TMath::Pi());
159  in[i] = fsin->Eval(x);
160  }
161 
162  //Make our own TVirtualFFT object (using option "K")
163  //Third parameter (option) consists of 3 parts:
164  //-transform type:
165  // real input/complex output in our case
166  //-transform flag:
167  // the amount of time spent in planning
168  // the transform (see TVirtualFFT class description)
169  //-to create a new TVirtualFFT object (option "K") or use the global (default)
170  Int_t n_size = n+1;
171  TVirtualFFT *fft_own = TVirtualFFT::FFT(1, &n_size, "R2C ES K");
172  if (!fft_own) return;
173  fft_own->SetPoints(in);
174  fft_own->Transform();
175 
176  //Copy all the output points:
177  fft_own->GetPoints(in);
178  //Draw the real part of the output
179  c1_5->cd();
180  TH1 *hr = 0;
181  hr = TH1::TransformHisto(fft_own, hr, "RE");
182  hr->SetTitle("Real part of the 3rd (array) tranfsorm");
183  hr->Draw();
184  hr->SetStats(kFALSE);
185  hr->GetXaxis()->SetLabelSize(0.05);
186  hr->GetYaxis()->SetLabelSize(0.05);
187  c1_6->cd();
188  TH1 *him = 0;
189  him = TH1::TransformHisto(fft_own, him, "IM");
190  him->SetTitle("Im. part of the 3rd (array) transform");
191  him->Draw();
192  him->SetStats(kFALSE);
193  him->GetXaxis()->SetLabelSize(0.05);
194  him->GetYaxis()->SetLabelSize(0.05);
195 
196  myc->cd();
197  //Now let's make another transform of the same size
198  //The same transform object can be used, as the size and the type of the transform
199  //haven't changed
200  TF1 *fcos = new TF1("fcos", "cos(x)+cos(0.5*x)+cos(2*x)+1", 0, 4*TMath::Pi());
201  for (Int_t i=0; i<=n; i++){
202  x = (Double_t(i)/n)*(4*TMath::Pi());
203  in[i] = fcos->Eval(x);
204  }
205  fft_own->SetPoints(in);
206  fft_own->Transform();
207  fft_own->GetPointComplex(0, re_2, im_2);
208  printf("2nd transform: DC component: %f\n", re_2);
209  fft_own->GetPointComplex(n/2+1, re_2, im_2);
210  printf("2nd transform: Nyquist harmonic: %f\n", re_2);
211  delete fft_own;
212  delete [] in;
213  delete [] re_full;
214  delete [] im_full;
215 }
216 
static void SetTransform(TVirtualFFT *fft)
static: set the current transfrom to parameter
virtual void SetPointsComplex(const Double_t *re, const Double_t *im)=0
static TVirtualFFT * GetCurrentTransform()
static: return current fgFFT
TVirtualPad * cd(Int_t subpadnumber=0)
Set current canvas & pad.
Definition: TCanvas.cxx:659
int Int_t
Definition: RtypesCore.h:41
const Bool_t kFALSE
Definition: Rtypes.h:92
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1052
static void AddDirectory(Bool_t add=kTRUE)
Sets the flag controlling the automatic add of histograms in memory.
Definition: TH1.cxx:1231
TVirtualPad * cd(Int_t subpadnumber=0)
Set Current pad.
Definition: TPad.cxx:514
Double_t x[n]
Definition: legend1.C:17
virtual void Draw(Option_t *option="")
Draw Pad in Current pad (re-parent pad if necessary).
Definition: TPad.cxx:1192
virtual TH1 * FFT(TH1 *h_output, Option_t *option)
This function allows to do discrete Fourier transforms of TH1 and TH2.
Definition: TH1.cxx:3098
virtual void Draw(Option_t *option="")
Draw this histogram with options.
Definition: TH1.cxx:2878
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
static TH1 * TransformHisto(TVirtualFFT *fft, TH1 *h_output, Option_t *option)
For a given transform (first parameter), fills the histogram (second parameter) with the transform ou...
Definition: TH1.cxx:8647
static TVirtualFFT * FFT(Int_t ndim, Int_t *n, Option_t *option)
Returns a pointer to the FFT of requested size and type.
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
The most important graphics class in the ROOT system.
Definition: TPad.h:46
virtual void Transform()=0
virtual void GetPointsComplex(Double_t *re, Double_t *im, Bool_t fromInput=kFALSE) const =0
TAxis * GetYaxis()
Definition: TH1.h:320
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:613
virtual void SetLabelSize(Float_t size=0.04)
Set size of axis labels The size is expressed in per cent of the pad width.
Definition: TAttAxis.cxx:187
TAxis * GetYaxis() const
Get y axis of the function.
Definition: TF1.cxx:2083
virtual void SetPoints(const Double_t *data)=0
TVirtualFFT is an interface class for Fast Fourier Transforms.
Definition: TVirtualFFT.h:92
Double_t Pi()
Definition: TMath.h:44
virtual void GetPoints(Double_t *data, Bool_t fromInput=kFALSE) const =0
The Canvas class.
Definition: TCanvas.h:48
virtual void GetPointComplex(Int_t ipoint, Double_t &re, Double_t &im, Bool_t fromInput=kFALSE) const =0
double Double_t
Definition: RtypesCore.h:55
ClassImp(TMCParticle) void TMCParticle printf(": p=(%7.3f,%7.3f,%9.3f) ;", fPx, fPy, fPz)
The TH1 histogram class.
Definition: TH1.h:80
TAxis * GetXaxis() const
Get x axis of the function.
Definition: TF1.cxx:2072
1-Dim function class
Definition: TF1.h:149
virtual Double_t Eval(Double_t x, Double_t y=0, Double_t z=0, Double_t t=0) const
Evaluate this function.
Definition: TF1.cxx:1162
void SetFrameFillColor(Color_t color=1)
Definition: TAttPad.h:83
virtual void SetTitle(const char *title)
Change (i.e.
Definition: TH1.cxx:6268
void FFT()
Definition: FFT.C:7
const Int_t n
Definition: legend1.C:16
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320
TAxis * GetXaxis()
Definition: TH1.h:319