Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
multifit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_fit
3/// \notebook -js
4/// Fitting multiple functions to different ranges of a 1-D histogram
5/// Example showing how to fit in a sub-range of an histogram
6/// An histogram is created and filled with the bin contents and errors
7/// defined in the table below.
8/// 3 gaussians are fitted in sub-ranges of this histogram.
9/// A new function (a sum of 3 gaussians) is fitted on another subrange
10/// Note that when fitting simple functions, such as gaussians, the initial
11/// values of parameters are automatically computed by ROOT.
12/// In the more complicated case of the sum of 3 gaussians, the initial values
13/// of parameters must be given. In this particular case, the initial values
14/// are taken from the result of the individual fits.
15///
16/// \macro_image
17/// \macro_output
18/// \macro_code
19///
20/// \author Rene Brun
21
22#include "TH1.h"
23#include "TF1.h"
24
25void multifit() {
26 const Int_t np = 49;
27 Float_t x[np] = {1.913521, 1.953769, 2.347435, 2.883654, 3.493567,
28 4.047560, 4.337210, 4.364347, 4.563004, 5.054247,
29 5.194183, 5.380521, 5.303213, 5.384578, 5.563983,
30 5.728500, 5.685752, 5.080029, 4.251809, 3.372246,
31 2.207432, 1.227541, 0.8597788,0.8220503,0.8046592,
32 0.7684097,0.7469761,0.8019787,0.8362375,0.8744895,
33 0.9143721,0.9462768,0.9285364,0.8954604,0.8410891,
34 0.7853871,0.7100883,0.6938808,0.7363682,0.7032954,
35 0.6029015,0.5600163,0.7477068,1.188785, 1.938228,
36 2.602717, 3.472962, 4.465014, 5.177035};
37
38 TH1F *h = new TH1F("h","Example of several fits in subranges",np,85,134);
39 h->SetMaximum(7);
40
41 for (int i=0;i<np;i++) {
42 h->SetBinContent(i+1,x[i]);
43 }
44
45 Double_t par[9];
46 TF1 *g1 = new TF1("g1","gaus",85,95);
47 TF1 *g2 = new TF1("g2","gaus",98,108);
48 TF1 *g3 = new TF1("g3","gaus",110,121);
49 TF1 *total = new TF1("total","gaus(0)+gaus(3)+gaus(6)",85,125);
50 total->SetLineColor(2);
51 h->Fit(g1,"R");
52 h->Fit(g2,"R+");
53 h->Fit(g3,"R+");
54 g1->GetParameters(&par[0]);
55 g2->GetParameters(&par[3]);
56 g3->GetParameters(&par[6]);
57 total->SetParameters(par);
58 h->Fit(total,"R+");
59}
60
#define h(i)
Definition RSha256.hxx:106
int Int_t
Definition RtypesCore.h:45
double Double_t
Definition RtypesCore.h:59
float Float_t
Definition RtypesCore.h:57
static unsigned int total
1-Dim function class
Definition TF1.h:213
virtual Double_t * GetParameters() const
Definition TF1.h:520
1-D histogram with a float per channel (see TH1 documentation)}
Definition TH1.h:575
Double_t x[n]
Definition legend1.C:17