approx

Macro to test interpolation function Approx

Author: Christian Stratowa, Vienna, Austria.
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Tuesday, November 29, 2022 at 11:16 AM.

In [1]:
TCanvas *vC1;
TGraph *grxy, *grin, *grout;

Definition of a helper function:

In [2]:
%%cpp -d
void DrawSmooth(Int_t pad, const char *title, const char *xt, const char *yt)
{
  vC1->cd(pad);
  TH1F *vFrame = gPad->DrawFrame(0,0,15,150);
  vFrame->SetTitle(title);
  vFrame->SetTitleSize(0.2);
  vFrame->SetXTitle(xt);
  vFrame->SetYTitle(yt);
  grxy->SetMarkerColor(kBlue);
  grxy->SetMarkerStyle(21);
  grxy->SetMarkerSize(0.5);
  grxy->Draw("P");
  grin->SetMarkerColor(kRed);
  grin->SetMarkerStyle(5);
  grin->SetMarkerSize(0.7);
  grin->Draw("P");
  grout->DrawClone("LP");
}

Test data (square)

In [3]:
Int_t n = 11;
Double_t x[] = {1,2,3,4,5,6,6,6,8,9,10};
Double_t y[] = {1,4,9,16,25,25,36,49,64,81,100};
grxy = new TGraph(n,x,y);

X values, for which y values should be interpolated

In [4]:
Int_t nout = 14;
Double_t xout[] =
   {1.2,1.7,2.5,3.2,4.4,5.2,5.7,6.5,7.6,8.3,9.7,10.4,11.3,13};

Create Canvas

In [5]:
vC1 = new TCanvas("vC1","square",200,10,700,700);
vC1->Divide(2,2);

Initialize graph with data

In [6]:
grin = new TGraph(n,x,y);

Interpolate at equidistant points (use mean for tied x-values)

In [7]:
TGraphSmooth *gs = new TGraphSmooth("normal");
grout = gs->Approx(grin,"linear");
DrawSmooth(1,"Approx: ties = mean","X-axis","Y-axis");

Re-initialize graph with data (since graph points were set to unique vales)

In [8]:
grin = new TGraph(n,x,y);

Interpolate at given points xout

In [9]:
grout = gs->Approx(grin,"linear", 14, xout, 0, 130);
DrawSmooth(2,"Approx: ties = mean","","");

Print output variables for given values xout

In [10]:
Int_t vNout = grout->GetN();
Double_t vXout, vYout;
for (Int_t k=0;k<vNout;k++) {
   grout->GetPoint(k, vXout, vYout);
   cout << "k= " << k << "  vXout[k]= " << vXout
        << "  vYout[k]= " << vYout << endl;
}
k= 0  vXout[k]= 1.2  vYout[k]= 1.6
k= 1  vXout[k]= 1.7  vYout[k]= 3.1
k= 2  vXout[k]= 2.5  vYout[k]= 6.5
k= 3  vXout[k]= 3.2  vYout[k]= 10.4
k= 4  vXout[k]= 4.4  vYout[k]= 19.6
k= 5  vXout[k]= 5.2  vYout[k]= 27.3333
k= 6  vXout[k]= 5.7  vYout[k]= 33.1667
k= 7  vXout[k]= 6.5  vYout[k]= 43.5
k= 8  vXout[k]= 7.6  vYout[k]= 58.5333
k= 9  vXout[k]= 8.3  vYout[k]= 69.1
k= 10  vXout[k]= 9.7  vYout[k]= 94.3
k= 11  vXout[k]= 10.4  vYout[k]= 130
k= 12  vXout[k]= 11.3  vYout[k]= 130
k= 13  vXout[k]= 13  vYout[k]= 130

Re-initialize graph with data

In [11]:
grin = new TGraph(n,x,y);

Interpolate at equidistant points (use min for tied x-values) grout = gs->Approx(grin,"linear", 50, 0, 0, 0, 1, 0, "min");

In [12]:
grout = gs->Approx(grin,"constant", 50, 0, 0, 0, 1, 0.5, "min");
DrawSmooth(3,"Approx: ties = min","","");

Re-initialize graph with data

In [13]:
grin = new TGraph(n,x,y);

Interpolate at equidistant points (use max for tied x-values)

In [14]:
grout = gs->Approx(grin,"linear", 14, xout, 0, 0, 2, 0, "max");
DrawSmooth(4,"Approx: ties = max","","");

Cleanup

In [15]:
delete gs;

Draw all canvases

In [16]:
%jsroot on
gROOT->GetListOfCanvases()->Draw()