FFT.C

Go to the documentation of this file.
00001 #include "TH1D.h"
00002 #include "TVirtualFFT.h"
00003 #include "TF1.h"
00004 #include "TCanvas.h"
00005 #include "TMath.h"
00006 
00007 void FFT()
00008 {
00009 
00010 //This tutorial illustrates the Fast Fourier Transforms interface in ROOT.
00011 //FFT transform types provided in ROOT:
00012 // - "C2CFORWARD" - a complex input/output discrete Fourier transform (DFT) 
00013 //                  in one or more dimensions, -1 in the exponent
00014 // - "C2CBACKWARD"- a complex input/output discrete Fourier transform (DFT) 
00015 //                  in one or more dimensions, +1 in the exponent
00016 // - "R2C"        - a real-input/complex-output discrete Fourier transform (DFT)
00017 //                  in one or more dimensions,
00018 // - "C2R"        - inverse transforms to "R2C", taking complex input 
00019 //                  (storing the non-redundant half of a logically Hermitian array) 
00020 //                  to real output
00021 // - "R2HC"       - a real-input DFT with output in ¡Èhalfcomplex¡É format, 
00022 //                  i.e. real and imaginary parts for a transform of size n stored as
00023 //                  r0, r1, r2, ..., rn/2, i(n+1)/2-1, ..., i2, i1
00024 // - "HC2R"       - computes the reverse of FFTW_R2HC, above
00025 // - "DHT"        - computes a discrete Hartley transform
00026 // Sine/cosine transforms:
00027 //  DCT-I  (REDFT00 in FFTW3 notation)
00028 //  DCT-II (REDFT10 in FFTW3 notation)
00029 //  DCT-III(REDFT01 in FFTW3 notation)
00030 //  DCT-IV (REDFT11 in FFTW3 notation)
00031 //  DST-I  (RODFT00 in FFTW3 notation)
00032 //  DST-II (RODFT10 in FFTW3 notation)
00033 //  DST-III(RODFT01 in FFTW3 notation)
00034 //  DST-IV (RODFT11 in FFTW3 notation)
00035 //First part of the tutorial shows how to transform the histograms
00036 //Second part shows how to transform the data arrays directly
00037 //Authors: Anna Kreshuk and Jens Hoffmann
00038 
00039 
00040 //********* Histograms ********//
00041 
00042 
00043    //prepare the canvas for drawing
00044    TCanvas *myc = new TCanvas("myc", "Fast Fourier Transform", 800, 600);
00045    myc->SetFillColor(45);
00046    TPad *c1_1 = new TPad("c1_1", "c1_1",0.01,0.67,0.49,0.99);
00047    TPad *c1_2 = new TPad("c1_2", "c1_2",0.51,0.67,0.99,0.99);
00048    TPad *c1_3 = new TPad("c1_3", "c1_3",0.01,0.34,0.49,0.65);
00049    TPad *c1_4 = new TPad("c1_4", "c1_4",0.51,0.34,0.99,0.65);
00050    TPad *c1_5 = new TPad("c1_5", "c1_5",0.01,0.01,0.49,0.32);
00051    TPad *c1_6 = new TPad("c1_6", "c1_6",0.51,0.01,0.99,0.32);
00052    c1_1->Draw();
00053    c1_2->Draw();
00054    c1_3->Draw();
00055    c1_4->Draw();
00056    c1_5->Draw();
00057    c1_6->Draw();
00058    c1_1->SetFillColor(30);
00059    c1_1->SetFrameFillColor(42);
00060    c1_2->SetFillColor(30);
00061    c1_2->SetFrameFillColor(42);
00062    c1_3->SetFillColor(30);
00063    c1_3->SetFrameFillColor(42);
00064    c1_4->SetFillColor(30);
00065    c1_4->SetFrameFillColor(42);
00066    c1_5->SetFillColor(30);
00067    c1_5->SetFrameFillColor(42);
00068    c1_6->SetFillColor(30);
00069    c1_6->SetFrameFillColor(42);
00070    
00071    c1_1->cd();
00072    TH1::AddDirectory(kFALSE);
00073      
00074    //A function to sample
00075    TF1 *fsin = new TF1("fsin", "sin(x)+sin(2*x)+sin(0.5*x)+1", 0, 4*TMath::Pi());
00076    fsin->Draw();
00077    
00078    Int_t n=25;
00079    TH1D *hsin = new TH1D("hsin", "hsin", n+1, 0, 4*TMath::Pi());
00080    Double_t x;
00081    
00082    //Fill the histogram with function values
00083    for (Int_t i=0; i<=n; i++){
00084       x = (Double_t(i)/n)*(4*TMath::Pi());
00085       hsin->SetBinContent(i+1, fsin->Eval(x));
00086    }
00087    hsin->Draw("same");
00088    fsin->GetXaxis()->SetLabelSize(0.05);
00089    fsin->GetYaxis()->SetLabelSize(0.05);
00090    
00091    c1_2->cd();
00092    //Compute the transform and look at the magnitude of the output
00093    TH1 *hm =0;
00094    TVirtualFFT::SetTransform(0);
00095    hm = hsin->FFT(hm, "MAG");
00096    hm->SetTitle("Magnitude of the 1st transform");
00097    hm->Draw();
00098    //NOTE: for "real" frequencies you have to divide the x-axes range with the range of your function 
00099    //(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!
00100    
00101    hm->SetStats(kFALSE);
00102    hm->GetXaxis()->SetLabelSize(0.05);
00103    hm->GetYaxis()->SetLabelSize(0.05);
00104    c1_3->cd();   
00105    //Look at the phase of the output   
00106    TH1 *hp = 0;
00107    hp = hsin->FFT(hp, "PH");
00108    hp->SetTitle("Phase of the 1st transform");
00109    hp->Draw();
00110    hp->SetStats(kFALSE);
00111    hp->GetXaxis()->SetLabelSize(0.05);
00112    hp->GetYaxis()->SetLabelSize(0.05);
00113    
00114    //Look at the DC component and the Nyquist harmonic:
00115    Double_t re, im;
00116    //That's the way to get the current transform object:
00117    TVirtualFFT *fft = TVirtualFFT::GetCurrentTransform();
00118    c1_4->cd();
00119    //Use the following method to get just one point of the output
00120    fft->GetPointComplex(0, re, im);
00121    printf("1st transform: DC component: %f\n", re);
00122    fft->GetPointComplex(n/2+1, re, im);
00123    printf("1st transform: Nyquist harmonic: %f\n", re);
00124 
00125    //Use the following method to get the full output:
00126    Double_t *re_full = new Double_t[n];
00127    Double_t *im_full = new Double_t[n];
00128    fft->GetPointsComplex(re_full,im_full);
00129   
00130    //Now let's make a backward transform:
00131    TVirtualFFT *fft_back = TVirtualFFT::FFT(1, &n, "C2R M K");
00132    fft_back->SetPointsComplex(re_full,im_full);
00133    fft_back->Transform();
00134    TH1 *hb = 0;
00135    //Let's look at the output
00136    hb = TH1::TransformHisto(fft_back,hb,"Re");
00137    hb->SetTitle("The backward transform result");
00138    hb->Draw();
00139    //NOTE: here you get at the x-axes number of bins and not real values
00140    //(in this case 25 bins has to be rescaled to a range between 0 and 4*Pi; 
00141    //also here the y-axes has to be rescaled (factor 1/bins)
00142    hb->SetStats(kFALSE);
00143    hb->GetXaxis()->SetLabelSize(0.05);
00144    hb->GetYaxis()->SetLabelSize(0.05);
00145    delete fft_back;
00146    fft_back=0;
00147 
00148 //********* Data array - same transform ********//
00149 
00150    //Allocate an array big enough to hold the transform output
00151    //Transform output in 1d contains, for a transform of size N, 
00152    //N/2+1 complex numbers, i.e. 2*(N/2+1) real numbers
00153    //our transform is of size n+1, because the histogram has n+1 bins
00154 
00155    Double_t *in = new Double_t[2*((n+1)/2+1)];
00156    Double_t re_2,im_2;
00157    for (Int_t i=0; i<=n; i++){
00158       x = (Double_t(i)/n)*(4*TMath::Pi());
00159       in[i] =  fsin->Eval(x);
00160    }
00161 
00162    //Make our own TVirtualFFT object (using option "K")
00163    //Third parameter (option) consists of 3 parts:
00164    //-transform type:
00165    // real input/complex output in our case
00166    //-transform flag: 
00167    // the amount of time spent in planning
00168    // the transform (see TVirtualFFT class description)
00169    //-to create a new TVirtualFFT object (option "K") or use the global (default)
00170    Int_t n_size = n+1;
00171    TVirtualFFT *fft_own = TVirtualFFT::FFT(1, &n_size, "R2C ES K");
00172    if (!fft_own) return;
00173    fft_own->SetPoints(in);
00174    fft_own->Transform();
00175 
00176    //Copy all the output points:
00177    fft_own->GetPoints(in);
00178    //Draw the real part of the output
00179    c1_5->cd();
00180    TH1 *hr = 0;
00181    hr = TH1::TransformHisto(fft_own, hr, "RE");
00182    hr->SetTitle("Real part of the 3rd (array) tranfsorm");
00183    hr->Draw();
00184    hr->SetStats(kFALSE);
00185    hr->GetXaxis()->SetLabelSize(0.05);
00186    hr->GetYaxis()->SetLabelSize(0.05);
00187    c1_6->cd();
00188    TH1 *him = 0;
00189    him = TH1::TransformHisto(fft_own, him, "IM");
00190    him->SetTitle("Im. part of the 3rd (array) transform");
00191    him->Draw();
00192    him->SetStats(kFALSE);
00193    him->GetXaxis()->SetLabelSize(0.05);
00194    him->GetYaxis()->SetLabelSize(0.05);
00195 
00196    myc->cd();
00197    //Now let's make another transform of the same size
00198    //The same transform object can be used, as the size and the type of the transform
00199    //haven't changed
00200    TF1 *fcos = new TF1("fcos", "cos(x)+cos(0.5*x)+cos(2*x)+1", 0, 4*TMath::Pi());
00201    for (Int_t i=0; i<=n; i++){
00202       x = (Double_t(i)/n)*(4*TMath::Pi());
00203       in[i] =  fcos->Eval(x);
00204    }
00205    fft_own->SetPoints(in);
00206    fft_own->Transform();
00207    fft_own->GetPointComplex(0, re_2, im_2);
00208    printf("2nd transform: DC component: %f\n", re_2);
00209    fft_own->GetPointComplex(n/2+1, re_2, im_2);
00210    printf("2nd transform: Nyquist harmonic: %f\n", re_2);
00211    delete fft_own;
00212    delete [] in;
00213    delete [] re_full;
00214    delete [] im_full;
00215 }
00216 

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