ROOT  6.07/01
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
FittingDemo.C
Go to the documentation of this file.
1 //Example for fitting signal/background.
2 // This example can be executed with:
3 // root > .x FittingDemo.C (using the CINT interpreter)
4 // root > .x FittingDemo.C+ (using the native complier via ACLIC)
5 //Author: Rene Brun
6 
7 #include "TH1.h"
8 #include "TMath.h"
9 #include "TF1.h"
10 #include "TLegend.h"
11 #include "TCanvas.h"
12 
13 // Quadratic background function
15  return par[0] + par[1]*x[0] + par[2]*x[0]*x[0];
16 }
17 
18 
19 // Lorenzian Peak function
21  return (0.5*par[0]*par[1]/TMath::Pi()) /
22  TMath::Max( 1.e-10,(x[0]-par[2])*(x[0]-par[2])
23  + .25*par[1]*par[1]);
24 }
25 
26 // Sum of background and peak function
28  return background(x,par) + lorentzianPeak(x,&par[3]);
29 }
30 
31 void FittingDemo() {
32  //Bevington Exercise by Peter Malzacher, modified by Rene Brun
33 
34  const int nBins = 60;
35 
36  Double_t data[nBins] = { 6, 1,10,12, 6,13,23,22,15,21,
37  23,26,36,25,27,35,40,44,66,81,
38  75,57,48,45,46,41,35,36,53,32,
39  40,37,38,31,36,44,42,37,32,32,
40  43,44,35,33,33,39,29,41,32,44,
41  26,39,29,35,32,21,21,15,25,15};
42  TCanvas *c1 = new TCanvas("c1","Fitting Demo",10,10,700,500);
43  c1->SetFillColor(33);
44  c1->SetFrameFillColor(41);
45  c1->SetGrid();
46 
47  TH1F *histo = new TH1F("histo",
48  "Lorentzian Peak on Quadratic Background",60,0,3);
49  histo->SetMarkerStyle(21);
50  histo->SetMarkerSize(0.8);
51  histo->SetStats(0);
52 
53  for(int i=0; i < nBins; i++) histo->SetBinContent(i+1,data[i]);
54 
55  // create a TF1 with the range from 0 to 3 and 6 parameters
56  TF1 *fitFcn = new TF1("fitFcn",fitFunction,0,3,6);
57  fitFcn->SetNpx(500);
58  fitFcn->SetLineWidth(4);
59  fitFcn->SetLineColor(kMagenta);
60 
61  // first try without starting values for the parameters
62  // This defaults to 1 for each param.
63  // this results in an ok fit for the polynomial function
64  // however the non-linear part (lorenzian) does not
65  // respond well.
66  fitFcn->SetParameters(1,1,1,1,1,1);
67  histo->Fit("fitFcn","0");
68 
69  // second try: set start values for some parameters
70  fitFcn->SetParameter(4,0.2); // width
71  fitFcn->SetParameter(5,1); // peak
72 
73  histo->Fit("fitFcn","V+","ep");
74 
75  // improve the picture:
76  TF1 *backFcn = new TF1("backFcn",background,0,3,3);
77  backFcn->SetLineColor(kRed);
78  TF1 *signalFcn = new TF1("signalFcn",lorentzianPeak,0,3,3);
79  signalFcn->SetLineColor(kBlue);
80  signalFcn->SetNpx(500);
81  Double_t par[6];
82 
83  // writes the fit results into the par array
84  fitFcn->GetParameters(par);
85 
86  backFcn->SetParameters(par);
87  backFcn->Draw("same");
88 
89  signalFcn->SetParameters(&par[3]);
90  signalFcn->Draw("same");
91 
92  // draw the legend
93  TLegend *legend=new TLegend(0.6,0.65,0.88,0.85);
94  legend->SetTextFont(72);
95  legend->SetTextSize(0.04);
96  legend->AddEntry(histo,"Data","lpe");
97  legend->AddEntry(backFcn,"Background fit","l");
98  legend->AddEntry(signalFcn,"Signal fit","l");
99  legend->AddEntry(fitFcn,"Global Fit","l");
100  legend->Draw();
101 
102 }
virtual void SetLineWidth(Width_t lwidth)
Definition: TAttLine.h:57
double par[1]
Definition: unuranDistr.cxx:38
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:432
void FittingDemo()
Definition: FittingDemo.C:31
virtual void SetNpx(Int_t npx=100)
Set the number of points used to draw the function.
Definition: TF1.cxx:3116
This class displays a legend box (TPaveText) containing several legend entries.
Definition: TLegend.h:35
TCanvas * c1
Definition: legend1.C:2
Definition: Rtypes.h:61
virtual void Draw(Option_t *option="")
Draw this legend with its current attributes.
Definition: TLegend.cxx:373
TH1 * histo
1-D histogram with a float per channel (see TH1 documentation)}
Definition: TH1.h:570
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1052
TF1 * fitFcn
TLegend * legend
Definition: pirndm.C:35
virtual void SetTextFont(Font_t tfont=62)
Definition: TAttText.h:59
Double_t x[n]
Definition: legend1.C:17
Double_t fitFunction(Double_t *x, Double_t *par)
Definition: FittingDemo.C:27
virtual void SetGrid(Int_t valuex=1, Int_t valuey=1)
Definition: TPad.h:326
virtual void SetLineColor(Color_t lcolor)
Definition: TAttLine.h:54
virtual void SetFillColor(Color_t fcolor)
Definition: TAttFill.h:50
virtual void SetBinContent(Int_t bin, Double_t content)
Set bin content see convention for numbering bins in TH1::GetBin In case the bin number is greater th...
Definition: TH1.cxx:8543
virtual void SetMarkerStyle(Style_t mstyle=1)
Definition: TAttMarker.h:53
virtual void SetMarkerSize(Size_t msize=1)
Definition: TAttMarker.h:54
Double_t Pi()
Definition: TMath.h:44
The Canvas class.
Definition: TCanvas.h:48
Double_t lorentzianPeak(Double_t *x, Double_t *par)
Definition: FittingDemo.C:20
double Double_t
Definition: RtypesCore.h:55
TLegendEntry * AddEntry(const TObject *obj, const char *label="", Option_t *option="lpf")
Add a new entry to this legend.
Definition: TLegend.cxx:280
virtual Double_t * GetParameters() const
Definition: TF1.h:358
Short_t Max(Short_t a, Short_t b)
Definition: TMathBase.h:202
1-Dim function class
Definition: TF1.h:149
Double_t background(Double_t *x, Double_t *par)
Definition: FittingDemo.C:14
Definition: Rtypes.h:61
void SetFrameFillColor(Color_t color=1)
Definition: TAttPad.h:83
virtual void SetParameter(Int_t param, Double_t value)
Definition: TF1.h:424
virtual void SetTextSize(Float_t tsize=1)
Definition: TAttText.h:60
virtual TFitResultPtr Fit(const char *formula, Option_t *option="", Option_t *goption="", Double_t xmin=0, Double_t xmax=0)
Fit histogram with function fname.
Definition: TH1.cxx:3607
virtual void SetStats(Bool_t stats=kTRUE)
Set statistics option on/off.
Definition: TH1.cxx:8320