Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fitNormSum.C
Go to the documentation of this file.
1/// \file
2/// \ingroup tutorial_fit
3/// \notebook
4/// Tutorial for normalized sum of two functions
5/// Here: a background exponential and a crystalball function
6/// Parameters can be set:
7/// 1. with the TF1 object before adding the function (for 3) and 4))
8/// 2. with the TF1NormSum object (first two are the coefficients, then the non constant parameters)
9/// 3. with the TF1 object after adding the function
10///
11/// Sum can be constructed by:
12/// 1. by a string containing the names of the functions and/or the coefficient in front
13/// 2. by a string containg formulas like expo, gaus...
14/// 3. by the list of functions and coefficients (which are 1 by default)
15/// 4. by a std::vector for functions and coefficients
16///
17/// \macro_image
18/// \macro_output
19/// \macro_code
20///
21/// \author Lorenzo Moneta
22
24#include <TCanvas.h>
25#include <TF1.h>
26#include <TF1NormSum.h>
27#include <TFitResult.h>
28#include <TH1.h>
29#include <TLatex.h>
30#include <TMath.h>
31#include <TStopwatch.h>
32#include <TStyle.h>
33
34void fitNormSum()
35{
36 const int nsig = 5.E4;
37 const int nbkg = 1.e6;
38 int nEvents = nsig + nbkg;
39 int nBins = 1e3;
40
41 double signal_mean = 3;
42 TF1 *f_cb = new TF1("MyCrystalBall", "crystalball", -5., 5.);
43 TF1 *f_exp = new TF1("MyExponential", "expo", -5., 5.);
44
45 // I.:
46 f_exp->SetParameters(1., -0.3);
47 f_cb->SetParameters(1, signal_mean, 0.3, 2, 1.5);
48
49 // CONSTRUCTION OF THE TF1NORMSUM OBJECT ........................................
50 // 1) :
51 TF1NormSum *fnorm_exp_cb = new TF1NormSum(f_cb, f_exp, nsig, nbkg);
52 // 4) :
53
54 TF1 *f_sum = new TF1("fsum", *fnorm_exp_cb, -5., 5., fnorm_exp_cb->GetNpar());
55
56 // III.:
57 f_sum->SetParameters(fnorm_exp_cb->GetParameters().data());
58 f_sum->SetParName(1, "NBackground");
59 f_sum->SetParName(0, "NSignal");
60 for (int i = 2; i < f_sum->GetNpar(); ++i)
61 f_sum->SetParName(i, fnorm_exp_cb->GetParName(i));
62
63 // GENERATE HISTOGRAM TO FIT ..............................................................
65 w.Start();
66 TH1D *h_sum = new TH1D("h_ExpCB", "Exponential Bkg + CrystalBall function", nBins, -5., 5.);
67 h_sum->FillRandom("fsum", nEvents);
68 printf("Time to generate %d events: ", nEvents);
69 w.Print();
70
71 // need to scale histogram with width since we are fitting a density
72 h_sum->Sumw2();
73 h_sum->Scale(1., "width");
74
75 // fit - use Minuit2 if available
77 new TCanvas("Fit", "Fit", 800, 1000);
78 // do a least-square fit of the spectrum
79 auto result = h_sum->Fit("fsum", "SQ");
80 result->Print();
81 h_sum->Draw();
82 printf("Time to fit using ROOT TF1Normsum: ");
83 w.Print();
84
85 // test if parameters are fine
86 std::vector<double> pref = {nsig, nbkg, signal_mean};
87 for (unsigned int i = 0; i < pref.size(); ++i) {
88 if (!TMath::AreEqualAbs(pref[i], f_sum->GetParameter(i), f_sum->GetParError(i) * 10.))
89 Error("testFitNormSum", "Difference found in fitted %s - difference is %g sigma", f_sum->GetParName(i),
90 (f_sum->GetParameter(i) - pref[i]) / f_sum->GetParError(i));
91 }
92
94 // add parameters
95 auto t1 = new TLatex(
96 -2.5, 300000, TString::Format("%s = %8.0f #pm %4.0f", "NSignal", f_sum->GetParameter(0), f_sum->GetParError(0)));
97 auto t2 = new TLatex(
98 -2.5, 270000, TString::Format("%s = %8.0f #pm %4.0f", "Nbackgr", f_sum->GetParameter(1), f_sum->GetParError(1)));
99 t1->Draw();
100 t2->Draw();
101}
void Error(const char *location, const char *msgfmt,...)
Use this function in case an error occurred.
Definition TError.cxx:185
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
R__EXTERN TStyle * gStyle
Definition TStyle.h:433
static void SetDefaultMinimizer(const char *type, const char *algo=nullptr)
Set the default Minimizer type and corresponding algorithms.
The Canvas class.
Definition TCanvas.h:23
Class adding two functions: c1*f1+c2*f2.
Definition TF1NormSum.h:19
const char * GetParName(Int_t ipar) const
Definition TF1NormSum.h:66
std::vector< double > GetParameters() const
Return array of parameters.
Int_t GetNpar() const
Return the number of (non constant) parameters including the coefficients: for 2 functions: c1,...
1-Dim function class
Definition TF1.h:233
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition TF1.cxx:1932
virtual Int_t GetNpar() const
Definition TF1.h:507
virtual const char * GetParName(Int_t ipar) const
Definition TF1.h:555
virtual void SetParName(Int_t ipar, const char *name)
Set name of parameter number ipar.
Definition TF1.cxx:3450
virtual void SetParameters(const Double_t *params)
Definition TF1.h:670
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:538
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:669
virtual void FillRandom(const char *fname, Int_t ntimes=5000, TRandom *rng=nullptr)
Fill histogram following distribution in function fname.
Definition TH1.cxx:3519
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:3898
void Draw(Option_t *option="") override
Draw this histogram with options.
Definition TH1.cxx:3066
virtual void Scale(Double_t c1=1, Option_t *option="")
Multiply this histogram by a constant c1.
Definition TH1.cxx:6572
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition TH1.cxx:8988
To draw Mathematical Formula.
Definition TLatex.h:18
Stopwatch class.
Definition TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
static TString Format(const char *fmt,...)
Static method which formats a string using a printf style format descriptor and return a TString.
Definition TString.cxx:2378
void SetOptStat(Int_t stat=1)
The type of information printed in the histogram statistics box can be selected via the parameter mod...
Definition TStyle.cxx:1636
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Comparing floating points.
Definition TMath.h:418
auto * t1
Definition textangle.C:20