Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
combinedFit.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_fit
3/// \notebook
4/// Combined (simultaneous) fit of two histogram with separate functions
5/// and some common parameters
6///
7/// See http://root.cern.ch/phpBB3//viewtopic.php?f=3&t=11740#p50908
8/// for a modified version working with Fumili or GSLMultiFit
9///
10/// N.B. this macro must be compiled with ACliC
11///
12/// \macro_image
13/// \macro_output
14/// \macro_code
15///
16/// \author Lorenzo Moneta
17
18#include <Fit/Fitter.h>
19#include <Fit/BinData.h>
20#include <Fit/Chi2FCN.h>
21#include <TH1.h>
23#include <HFitInterface.h>
24#include <TCanvas.h>
25#include <TStyle.h>
26
27// definition of shared parameter
28// background function
29int iparB[2] = {
30 0, // exp amplitude in B histo
31 2 // exp common parameter
32};
33
34// signal + background function
35int iparSB[5] = {
36 1, // exp amplitude in S+B histo
37 2, // exp common parameter
38 3, // Gaussian amplitude
39 4, // Gaussian mean
40 5 // Gaussian sigma
41};
42
43// Create the GlobalCHi2 structure
44
45struct GlobalChi2 {
46 GlobalChi2(ROOT::Math::IMultiGenFunction &f1, ROOT::Math::IMultiGenFunction &f2) : fChi2_1(&f1), fChi2_2(&f2) {}
47
48 // parameter vector is first background (in common 1 and 2)
49 // and then is signal (only in 2)
50 double operator()(const double *par) const
51 {
52 double p1[2];
53 for (int i = 0; i < 2; ++i)
54 p1[i] = par[iparB[i]];
55
56 double p2[5];
57 for (int i = 0; i < 5; ++i)
58 p2[i] = par[iparSB[i]];
59
60 return (*fChi2_1)(p1) + (*fChi2_2)(p2);
61 }
62
63 const ROOT::Math::IMultiGenFunction *fChi2_1;
64 const ROOT::Math::IMultiGenFunction *fChi2_2;
65};
66
67void combinedFit()
68{
69
70 TH1D *hB = new TH1D("hB", "histo B", 100, 0, 100);
71 TH1D *hSB = new TH1D("hSB", "histo S+B", 100, 0, 100);
72
73 TF1 *fB = new TF1("fB", "expo", 0, 100);
74 fB->SetParameters(1, -0.05);
75 hB->FillRandom("fB");
76
77 TF1 *fS = new TF1("fS", "gaus", 0, 100);
78 fS->SetParameters(1, 30, 5);
79
80 hSB->FillRandom("fB", 2000);
81 hSB->FillRandom("fS", 1000);
82
83 // perform now global fit
84
85 TF1 *fSB = new TF1("fSB", "expo + gaus(2)", 0, 100);
86
88 ROOT::Math::WrappedMultiTF1 wfSB(*fSB, 1);
89
92 // set the data range
93 rangeB.SetRange(10, 90);
94 ROOT::Fit::BinData dataB(opt, rangeB);
95 ROOT::Fit::FillData(dataB, hB);
96
98 rangeSB.SetRange(10, 50);
99 ROOT::Fit::BinData dataSB(opt, rangeSB);
100 ROOT::Fit::FillData(dataSB, hSB);
101
102 ROOT::Fit::Chi2Function chi2_B(dataB, wfB);
103 ROOT::Fit::Chi2Function chi2_SB(dataSB, wfSB);
104
105 GlobalChi2 globalChi2(chi2_B, chi2_SB);
106
107 ROOT::Fit::Fitter fitter;
108
109 const int Npar = 6;
110 double par0[Npar] = {5, 5, -0.1, 100, 30, 10};
111
112 // create before the parameter settings in order to fix or set range on them
113 fitter.Config().SetParamsSettings(6, par0);
114 // fix 5-th parameter
115 fitter.Config().ParSettings(4).Fix();
116 // set limits on the third and 4-th parameter
117 fitter.Config().ParSettings(2).SetLimits(-10, -1.E-4);
118 fitter.Config().ParSettings(3).SetLimits(0, 10000);
119 fitter.Config().ParSettings(3).SetStepSize(5);
120
122 fitter.Config().SetMinimizer("Minuit2", "Migrad");
123
124 // fit FCN function directly
125 // (specify optionally data size and flag to indicate that is a chi2 fit)
126 fitter.FitFCN(6, globalChi2, nullptr, dataB.Size() + dataSB.Size(), true);
128 result.Print(std::cout);
129
130 TCanvas *c1 = new TCanvas("Simfit", "Simultaneous fit of two histograms", 10, 10, 700, 700);
131 c1->Divide(1, 2);
132 c1->cd(1);
133 gStyle->SetOptFit(1111);
134
135 fB->SetFitResult(result, iparB);
136 fB->SetRange(rangeB().first, rangeB().second);
137 fB->SetLineColor(kBlue);
138 hB->GetListOfFunctions()->Add(fB);
139 hB->Draw();
140
141 c1->cd(2);
142 fSB->SetFitResult(result, iparSB);
143 fSB->SetRange(rangeSB().first, rangeSB().second);
144 fSB->SetLineColor(kRed);
145 hSB->GetListOfFunctions()->Add(fSB);
146 hSB->Draw();
147}
@ kRed
Definition Rtypes.h:66
@ kBlue
Definition Rtypes.h:66
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t result
TRObject operator()(const T1 &t1) const
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
Class describing the binned data sets : vectors of x coordinates, y values and optionally error on y ...
Definition BinData.h:52
Chi2FCN class for binned fits using the least square methods.
Definition Chi2FCN.h:46
class describing the range in the coordinates it supports multiple range in a coordinate.
Definition DataRange.h:35
void SetRange(unsigned int icoord, double xmin, double xmax)
set a range [xmin,xmax] for the new coordinate icoord If more range exists for other coordinates,...
void SetMinimizer(const char *type, const char *algo=nullptr)
set minimizer type
Definition FitConfig.h:181
void SetParamsSettings(unsigned int npar, const double *params, const double *vstep=nullptr)
set the parameter settings from number of parameters and a vector of values and optionally step value...
const ParameterSettings & ParSettings(unsigned int i) const
get the parameter settings for the i-th parameter (const method)
Definition FitConfig.h:76
ROOT::Math::MinimizerOptions & MinimizerOptions()
access to the minimizer control parameter (non const method)
Definition FitConfig.h:167
class containing the result of the fit and all the related information (fitted parameter values,...
Definition FitResult.h:47
Fitter class, entry point for performing all type of fits.
Definition Fitter.h:77
const FitResult & Result() const
get fit result
Definition Fitter.h:393
const FitConfig & Config() const
access to the fit configuration (const method)
Definition Fitter.h:421
bool FitFCN(unsigned int npar, Function &fcn, const double *params=nullptr, unsigned int dataSize=0, bool chi2fit=false)
Fit using the a generic FCN function as a C++ callable object implementing double () (const double *)...
Definition Fitter.h:649
void SetStepSize(double err)
set the step size
void SetLimits(double low, double up)
set a double side limit, if low == up the parameter is fixed if low > up the limits are removed The c...
void Fix()
fix the parameter
Documentation for the abstract class IBaseFunctionMultiDim.
Definition IFunction.h:61
void SetPrintLevel(int level)
set print level
Class to Wrap a ROOT Function class (like TF1) in a IParamMultiFunction interface of multi-dimensions...
virtual void SetLineColor(Color_t lcolor)
Set the line color.
Definition TAttLine.h:40
The Canvas class.
Definition TCanvas.h:23
1-Dim function class
Definition TF1.h:214
virtual void SetRange(Double_t xmin, Double_t xmax)
Initialize the upper and lower bounds to draw the function.
Definition TF1.cxx:3522
virtual void SetParameters(const Double_t *params)
Definition TF1.h:650
virtual void SetFitResult(const ROOT::Fit::FitResult &result, const Int_t *indpar=nullptr)
Set the result from the fit parameter values, errors, chi2, etc... Optionally a pointer to a vector (...
Definition TF1.cxx:3355
1-D histogram with a double per channel (see TH1 documentation)}
Definition TH1.h:620
virtual void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr)
Fill histogram following distribution in function fname.
Definition TH1.cxx:3520
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3067
TList * GetListOfFunctions() const
Definition TH1.h:242
void Add(TObject *obj) override
Definition TList.h:81
void SetOptFit(Int_t fit=1)
The type of information about fit parameters printed in the histogram statistics box can be selected ...
Definition TStyle.cxx:1589
return c1
Definition legend1.C:41
TF1 * f1
Definition legend1.C:11
void FillData(BinData &dv, const TH1 *hist, TF1 *func=nullptr)
fill the data vector from a TH1.
Definition first.py:1
DataOptions : simple structure holding the options on how the data are filled.
Definition DataOptions.h:28