Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
multifit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_fit
3## \notebook
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## A histogram is created and filled with the bin contents and errors
7## defined in the table below.
8## Three 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## \authors Jonas Rembser, Rene Brun (C++ version)
21
22import ROOT
23
24import numpy as np
25
26n_x = 49
27
28# fmt: off
29x = np.array( [ 1.913521, 1.953769, 2.347435, 2.883654, 3.493567, 4.047560,
30 4.337210, 4.364347, 4.563004, 5.054247, 5.194183, 5.380521, 5.303213,
31 5.384578, 5.563983, 5.728500, 5.685752, 5.080029, 4.251809, 3.372246,
32 2.207432, 1.227541, 0.8597788, 0.8220503, 0.8046592, 0.7684097, 0.7469761,
33 0.8019787, 0.8362375, 0.8744895, 0.9143721, 0.9462768, 0.9285364,
34 0.8954604, 0.8410891, 0.7853871, 0.7100883, 0.6938808, 0.7363682,
35 0.7032954, 0.6029015, 0.5600163, 0.7477068, 1.188785, 1.938228, 2.602717,
36 3.472962, 4.465014, 5.177035, ], dtype=np.float32,)
37# fmt: on
38
39# The histogram are filled with bins defined in the array x.
40h = ROOT.TH1F("h", "Example of several fits in subranges", n_x, 85, 134)
41h.SetMaximum(7)
42
43for i, x_i in enumerate(x):
44 h.SetBinContent(i + 1, x[i])
45
46# Define the parameter array for the total function.
47par = np.zeros(9)
48
49# Three TF1 objects are created, one for each subrange.
50g1 = ROOT.TF1("g1", "gaus", 85, 95)
51g2 = ROOT.TF1("g2", "gaus", 98, 108)
52g3 = ROOT.TF1("g3", "gaus", 110, 121)
53
54# The total is the sum of the three, each has three parameters.
55total = ROOT.TF1("total", "gaus(0)+gaus(3)+gaus(6)", 85, 125)
56total.SetLineColor(2)
57
58# The canvas that the histograms and fit functions are drawn on.
59c = ROOT.TCanvas("multifit", "multifit", 800, 400)
60
61# Fit each function and add it to the list of functions. By default, TH1::Fit()
62# fits the function on the defined histogram range. You can specify the "R"
63# option in the second parameter of TH1::Fit() to restrict the fit to the range
64# specified in the TF1 constructor. Alternatively, you can also specify the
65# range in the call to TH1::Fit(), which we demonstrate here with the 3rd
66# Gaussian. The "+" option needs to be added to the later fits to not replace
67# existing fitted functions in the histogram.
68h.Fit(g1, "R")
69h.Fit(g2, "R+")
70h.Fit(g3, "+", "", 110, 121);
71
72# Get the parameters from the fit.
73g1.GetParameters(par[:3])
74g2.GetParameters(par[3:6])
75g3.GetParameters(par[6:])
76
77print(par)
78
79# Use the parameters on the sum.
80total.SetParameters(par)
81
82h.Draw()
83h.Fit(total, "R+")
84
85# Save the plot for later inspection.
86c.SaveAs("multifit.png")