testInterpolation.cxx

Go to the documentation of this file.
00001 //Example of the various methods of interpolation provided by the 
00002 // ROOT::Math::Interpolator class
00003 // 
00004 //Example can also be run in ROOT by doing : 
00005 // 
00006 // root> .x testInterpolation.cxx
00007 //
00008 #include "TGraph.h"
00009 #include "TAxis.h"
00010 #include "TCanvas.h"
00011 #include "TLegend.h"
00012 #include "TApplication.h"
00013 #include "Math/Interpolator.h"
00014 #include <iostream>
00015 #include <string>
00016 #include <cstdlib>
00017 
00018 #include <cmath>
00019 
00020 bool showGraphics = true;
00021 
00022 TGraph *grorig = 0;
00023 
00024 void interpolate( const  ROOT::Math::Interpolator & itp ) { 
00025 
00026    std::string type = itp.Type(); 
00027    std::cout << "\n" << type << "  interpolation:" << std::endl;
00028 
00029    std::cout << "x[i]     y[i]     deriv[i]     deriv2[i]    integral[i] \n" << std::endl; 
00030    // print result of interpolation
00031    const Int_t n = 51;
00032    Int_t i = 0;
00033    Float_t xcoord[n], ycoord[n];
00034    for ( double xi = 0; xi < 10; xi += 0.2) { 
00035       xcoord[i] = xi;
00036       ycoord[i] = itp.Eval(xi);
00037       double dyi = itp.Deriv(xi); 
00038       double dy2i = itp.Deriv2(xi);
00039       double igyi = itp.Integ(0, xi); 
00040       std::cout << xcoord[i]  << "  " << ycoord[i] << "  " << dyi << "  " << dy2i << "  " <<  igyi << std::endl;
00041       i++; 
00042    }
00043 
00044    TGraph *gr = new TGraph(n,xcoord,ycoord);
00045    gr->SetMarkerColor(kBlue);
00046    gr->SetMarkerStyle(7);
00047    gr->Draw("CP");
00048    TLegend * l = new TLegend(0.1,0.7,0.4,0.9); 
00049    l->AddEntry(grorig,"original data"); 
00050    l->AddEntry(gr,type.c_str());
00051    l->Draw();
00052 
00053    return;
00054 }
00055 
00056 
00057 void testInterpolation() {
00058 
00059    // Create data
00060    int n = 10; 
00061    Float_t xorig[11], yorig[11];
00062    std::vector<double> x(n+1); 
00063    std::vector<double> y(n+1); 
00064    for (int i = 0; i < n+1; ++i) {  
00065       x[i] = i + 0.5 * std::sin(i+0.0); 
00066       y[i] = i + std::cos(i * i+0.0); 
00067       xorig[i] = x[i];
00068       yorig[i] = y[i];
00069    } 
00070 
00071    TCanvas *c1 = new TCanvas("c1","Original (red), Linear (upper left), Polynomial (upper right), Spline , Spline periodic, Akima (lower left) and Akima Periodic (lower right) Interpolation",10,10,1000,800);
00072    c1->Divide(2,3);
00073    c1->cd(1);
00074 
00075    grorig = new TGraph(n+1,xorig,yorig);
00076    grorig->SetMarkerColor(kRed);
00077    grorig->SetMarkerStyle(20);
00078    grorig->GetYaxis()->SetRange(0,40);
00079    grorig->Draw("AP"); 
00080 
00081    //ROOT::Math::Interpolator itp1(x, y, ROOT::Math::Interpolation::kLINEAR); 
00082    ROOT::Math::Interpolator itp1(x.size(), ROOT::Math::Interpolation::kLINEAR); 
00083    itp1.SetData(x,y);
00084    interpolate(itp1);
00085 
00086 
00087    c1->cd(2);
00088    grorig->Draw("AP"); 
00089 
00090    ROOT::Math::Interpolator itp2(x, y, ROOT::Math::Interpolation::kPOLYNOMIAL); 
00091    interpolate(itp2);
00092 
00093 
00094    c1->cd(3);
00095    grorig->Draw("AP"); 
00096 
00097    //std::cout << "Cubic Spline Interpolation: " << std::endl;
00098    ROOT::Math::Interpolator itp3( 2*x.size(), ROOT::Math::Interpolation::kCSPLINE); 
00099    itp3.SetData(x.size(), &x[0], &y[0]);
00100    interpolate(itp3);
00101 
00102    c1->cd(4);
00103    grorig->Draw("AP"); 
00104 
00105    //std::cout << "Akima  Interpolation: " << std::endl;
00106    ROOT::Math::Interpolator itp4(x, y, ROOT::Math::Interpolation::kAKIMA); 
00107    interpolate(itp4);
00108 
00109    c1->cd(5);
00110    grorig->Draw("AP"); 
00111 
00112    //std::cout << "Cubic Spline Periodic Interpolation: " << std::endl;
00113    ROOT::Math::Interpolator itp5(x, y, ROOT::Math::Interpolation::kCSPLINE_PERIODIC); 
00114    interpolate(itp5);
00115 
00116 
00117 
00118    c1->cd(6);
00119    grorig->Draw("AP"); 
00120 
00121    //std::cout << "Akima Periodic Interpolation: " << std::endl;
00122    ROOT::Math::Interpolator itp6(x, y, ROOT::Math::Interpolation::kAKIMA_PERIODIC); 
00123    interpolate(itp6);
00124 
00125    std::cout << "\n***********************************" << std::endl;
00126    std::cout << "Using default Interpolation type: " << std::endl;
00127 
00128    ROOT::Math::Interpolator itp7; 
00129    itp7.SetData(x,y);
00130    interpolate(itp7);
00131 
00132 
00133 }
00134 
00135 #ifndef __CINT__
00136 int main(int argc, char **argv)
00137 {
00138    using std::cerr;
00139    using std::cout;
00140    using std::endl;
00141 
00142    if ( argc > 1 && argc != 2 )
00143    {
00144       cerr << "Usage: " << argv[0] << " [-ng]\n";
00145       cerr << "  where:\n";
00146       cerr << "     -ng : no graphics mode";
00147       cerr << endl;
00148       exit(1);
00149    }
00150 
00151    if ( argc == 2 && strcmp( argv[1], "-ng") == 0 ) 
00152    {
00153       showGraphics = false;
00154    }
00155 
00156    TApplication* theApp = 0;
00157 
00158    if ( showGraphics )
00159       theApp = new TApplication("App",&argc,argv);
00160 
00161    testInterpolation();
00162 
00163    if ( showGraphics )
00164    {
00165       theApp->Run();
00166       delete theApp;
00167       theApp = 0;
00168    }
00169 
00170    return 0;
00171 }
00172 #endif

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