TCanvas *vC1; TGraph *grxy, *grin, *grout; void approx() { //********************************************** // Macro to test interpolation function Approx // Author: Christian Stratowa, Vienna, Austria. // Created: 26 Aug 2001 //********************************************** // test data (square) 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 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 vC1 = new TCanvas("vC1","square",200,10,700,700); vC1->Divide(2,2); // initialize graph with data grin = new TGraph(n,x,y); // interpolate at equidistant points (use mean for tied x-values) 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) grin = new TGraph(n,x,y); // interpolate at given points xout grout = gs->Approx(grin,"linear", 14, xout, 0, 130); DrawSmooth(2,"Approx: ties = mean","",""); // print output variables for given values xout 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; } // re-initialize graph with data 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"); grout = gs->Approx(grin,"constant", 50, 0, 0, 0, 1, 0.5, "min"); DrawSmooth(3,"Approx: ties = min","",""); // re-initialize graph with data grin = new TGraph(n,x,y); // interpolate at equidistant points (use max for tied x-values) grout = gs->Approx(grin,"linear", 14, xout, 0, 0, 2, 0, "max"); DrawSmooth(4,"Approx: ties = max","",""); // cleanup delete gs; } 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"); } approx.C:1 approx.C:2 approx.C:3 approx.C:4 approx.C:5 approx.C:6 approx.C:7 approx.C:8 approx.C:9 approx.C:10 approx.C:11 approx.C:12 approx.C:13 approx.C:14 approx.C:15 approx.C:16 approx.C:17 approx.C:18 approx.C:19 approx.C:20 approx.C:21 approx.C:22 approx.C:23 approx.C:24 approx.C:25 approx.C:26 approx.C:27 approx.C:28 approx.C:29 approx.C:30 approx.C:31 approx.C:32 approx.C:33 approx.C:34 approx.C:35 approx.C:36 approx.C:37 approx.C:38 approx.C:39 approx.C:40 approx.C:41 approx.C:42 approx.C:43 approx.C:44 approx.C:45 approx.C:46 approx.C:47 approx.C:48 approx.C:49 approx.C:50 approx.C:51 approx.C:52 approx.C:53 approx.C:54 approx.C:55 approx.C:56 approx.C:57 approx.C:58 approx.C:59 approx.C:60 approx.C:61 approx.C:62 approx.C:63 approx.C:64 approx.C:65 approx.C:66 approx.C:67 approx.C:68 approx.C:69 approx.C:70 approx.C:71 approx.C:72 approx.C:73 approx.C:74 approx.C:75 approx.C:76 approx.C:77 approx.C:78 approx.C:79 approx.C:80 approx.C:81 approx.C:82 approx.C:83 approx.C:84 approx.C:85 approx.C:86 approx.C:87 approx.C:88 |
|