Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fitLinear.C File Reference

Detailed Description

View in nbviewer Open in SWAN Example of fitting with a linear function, using TLinearFitter This example is for a TGraphErrors, but it can also be used when fitting a histogram, a TGraph2D or a TMultiGraph

****************************************
Minimizer is Linear / Migrad
Chi2 = 36.5406
NDf = 36
p0 = -7.07142 +/- 0.0233493
p1 = -0.0194368 +/- 0.0354128
p2 = 2.03968 +/- 0.0136149
p3 = 1.00594 +/- 0.0139068
****************************************
Minimizer is Linear / Migrad
Chi2 = 46.7362
NDf = 38
p0 = 1.0005 +/- 0.0242765
p1 = 0.985942 +/- 0.0279149
****************************************
Minimizer is Linear / Migrad
Chi2 = 43.6161
NDf = 38
p0 = -2.04095 +/- 0.0220454
p1 = 1.01171 +/- 0.00904363
#include "TGraphErrors.h"
#include "TF1.h"
#include "TRandom.h"
#include "TCanvas.h"
#include "TLegend.h"
#include "TMath.h"
void makePoints(Int_t n, Double_t *x, Double_t *y, Double_t *e, Int_t p);
void fitLinear()
{
Int_t n = 40;
Double_t *x = new Double_t[n];
Double_t *y = new Double_t[n];
Double_t *e = new Double_t[n];
TCanvas *myc = new TCanvas("myc",
"Fitting 3 TGraphErrors with linear functions");
myc->SetGrid();
//Generate points along a 3rd degree polynomial:
makePoints(n, x, y, e, 3);
TGraphErrors *gre3 = new TGraphErrors(n, x, y, 0, e);
gre3->Draw("a*");
//Fit the graph with the predefined "pol3" function
gre3->Fit("pol3");
//Access the fit resuts
TF1 *f3 = gre3->GetFunction("pol3");
f3->SetLineWidth(1);
//Generate points along a sin(x)+sin(2x) function
makePoints(n, x, y, e, 2);
TGraphErrors *gre2=new TGraphErrors(n, x, y, 0, e);
gre2->Draw("*same");
//The fitting function can be predefined and passed to the Fit function
//The "++" mean that the linear fitter should be used, and the following
//formula is equivalent to "[0]*sin(x) + [1]*sin(2*x)"
//A function, defined this way, is in no way different from any other TF1,
//it can be evaluated, drawn, you can get its parameters, etc.
//The fit result (parameter values, parameter errors, chisquare, etc) are
//written into the fitting function.
TF1 *f2 = new TF1("f2", "sin(x) ++ sin(2*x)", -2, 2);
gre2->Fit(f2);
f2 = gre2->GetFunction("f2");
f2->SetLineWidth(1);
//Generate points along a -2+exp(-x) function
makePoints(n, x, y, e, 4);
TGraphErrors *gre4=new TGraphErrors(n, x, y, 0, e);
gre4->Draw("*same");
//If you don't want to define the function, you can just pass the string
//with the the formula:
gre4->Fit("1 ++ exp(-x)");
//Access the fit results:
TF1 *f4 = gre4->GetFunction("1 ++ exp(-x)");
f4->SetName("f4");
f4->SetLineWidth(1);
TLegend *leg = new TLegend(0.3, 0.7, 0.65, 0.9);
leg->AddEntry(gre3, " -7 + 2*x*x + x*x*x", "p");
leg->AddEntry(gre2, "sin(x) + sin(2*x)", "p");
leg->AddEntry(gre4, "-2 + exp(-x)", "p");
leg->Draw();
}
void makePoints(Int_t n, Double_t *x, Double_t *y, Double_t *e, Int_t p)
{
Int_t i;
if (p==2) {
for (i=0; i<n; i++) {
x[i] = r.Uniform(-2, 2);
y[i]=TMath::Sin(x[i]) + TMath::Sin(2*x[i]) + r.Gaus()*0.1;
e[i] = 0.1;
}
}
if (p==3) {
for (i=0; i<n; i++) {
x[i] = r.Uniform(-2, 2);
y[i] = -7 + 2*x[i]*x[i] + x[i]*x[i]*x[i]+ r.Gaus()*0.1;
e[i] = 0.1;
}
}
if (p==4) {
for (i=0; i<n; i++) {
x[i] = r.Uniform(-2, 2);
y[i]=-2 + TMath::Exp(-x[i]) + r.Gaus()*0.1;
e[i] = 0.1;
}
}
}
ROOT::R::TRInterface & r
Definition Object.C:4
#define e(i)
Definition RSha256.hxx:103
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
virtual void SetLineWidth(Width_t lwidth)
Set the line width.
Definition TAttLine.h:43
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
virtual void SetMarkerColor(Color_t mcolor=1)
Set the marker color.
Definition TAttMarker.h:38
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:213
A TGraphErrors is a TGraph with error bars.
TF1 * GetFunction(const char *name) const
Return pointer to function with name.
Definition TGraph.cxx:1473
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Axis_t xmin=0, Axis_t xmax=0)
Fit this graph with function with name fname.
Definition TGraph.cxx:1073
virtual void Draw(Option_t *chopt="")
Draw this graph with its current attributes.
Definition TGraph.cxx:769
This class displays a legend box (TPaveText) containing several legend entries.
Definition TLegend.h:23
virtual void SetName(const char *name)
Set the name of the TNamed.
Definition TNamed.cxx:140
void SetGrid(Int_t valuex=1, Int_t valuey=1) override
Definition TPad.h:326
This is the base class for the ROOT Random number generators.
Definition TRandom.h:27
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
leg
Definition legend1.C:34
Double_t Exp(Double_t x)
Definition TMath.h:677
Double_t Sin(Double_t)
Definition TMath.h:589
Author
Anna Kreshuk

Definition in file fitLinear.C.