Logo ROOT  
Reference Guide
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
23#include <TMath.h>
24#include <TCanvas.h>
25#include <TF1NormSum.h>
26#include <TF1.h>
27#include <TH1.h>
28
29
30using namespace std;
31
32
33void fitNormSum()
34{
35 const int nsig = 5.E4;
36 const int nbkg = 1.e6;
37 Int_t NEvents = nsig+nbkg;
38 Int_t NBins = 1e3;
39
40 double signal_mean = 3;
41 TF1 *f_cb = new TF1("MyCrystalBall","crystalball",-5.,5.);
42 TF1 *f_exp = new TF1("MyExponential","expo",-5.,5.);
43
44 // I.:
45 f_exp-> SetParameters(1.,-0.3);
46 f_cb -> SetParameters(1,signal_mean,0.3,2,1.5);
47
48 // CONSTRUCTION OF THE TF1NORMSUM OBJECT ........................................
49 // 1) :
50 TF1NormSum *fnorm_exp_cb = new TF1NormSum(f_cb,f_exp,nsig,nbkg);
51 // 4) :
52
53 TF1 * f_sum = new TF1("fsum", *fnorm_exp_cb, -5., 5., fnorm_exp_cb->GetNpar());
54 f_sum->Draw();
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 ..............................................................
64 TStopwatch w;
65 w.Start();
66 TH1D *h_sum = new TH1D("h_ExpCB", "Exponential Bkg + CrystalBall function", NBins, -5., 5.);
67 for (int i=0; i<NEvents; i++)
68 {
69 h_sum -> Fill(f_sum -> GetRandom());
70 }
71 printf("Time to generate %d events: ",NEvents);
72 w.Print();
73 //TH1F *h_orig = new TH1F(*h_sum);
74
75 // need to scale histogram with width since we are fitting a density
76 h_sum -> Sumw2();
77 h_sum -> Scale(1., "width");
78
79 //fit - use Minuit2 if available
81 new TCanvas("Fit","Fit",800,1000);
82 // do a least-square fit of the spectrum
83 auto result = h_sum -> Fit("fsum","SQ");
84 result->Print();
85 h_sum -> Draw();
86 printf("Time to fit using ROOT TF1Normsum: ");
87 w.Print();
88
89 // test if parameters are fine
90 std::vector<double> pref = {nsig, nbkg, signal_mean};
91 for (unsigned int i = 0; i< pref.size(); ++i) {
92 if (!TMath::AreEqualAbs(pref[i], f_sum->GetParameter(i), f_sum->GetParError(i)*10.) )
93 Error("testFitNormSum","Difference found in fitted %s - difference is %g sigma",f_sum->GetParName(i), (f_sum->GetParameter(i)-pref[i])/f_sum->GetParError(i));
94 }
95
97 // add parameters
98 auto t1 = new TLatex(-2.5, 300000, TString::Format("%s = %8.0f #pm %4.0f", "NSignal",f_sum->GetParameter(0), f_sum->GetParError(0) ) );
99 auto t2 = new TLatex(-2.5, 270000, TString::Format("%s = %8.0f #pm %4.0f", "Nbackgr",f_sum->GetParameter(1), f_sum->GetParError(1) ) );
100 t1->Draw();
101 t2->Draw();
102}
void Error(const char *location, const char *msgfmt,...)
void SetParameters(TFitEditor::FuncParams_t &pars, TF1 *func)
Restore the parameters from pars into the function.
Definition: TFitEditor.cxx:279
R__EXTERN TStyle * gStyle
Definition: TStyle.h:410
static void SetDefaultMinimizer(const char *type, const char *algo=0)
The Canvas class.
Definition: TCanvas.h:27
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.
Definition: TF1NormSum.cxx:292
Int_t GetNpar() const
Return the number of (non constant) parameters including the coefficients: for 2 functions: c1,...
Definition: TF1NormSum.cxx:364
1-Dim function class
Definition: TF1.h:210
virtual Double_t GetParError(Int_t ipar) const
Return value of parameter number ipar.
Definition: TF1.cxx:1913
virtual Int_t GetNpar() const
Definition: TF1.h:475
virtual const char * GetParName(Int_t ipar) const
Definition: TF1.h:523
virtual void SetParName(Int_t ipar, const char *name)
Set name of parameter number ipar.
Definition: TF1.cxx:3452
virtual void Draw(Option_t *option="")
Draw this function with its current attributes.
Definition: TF1.cxx:1320
virtual void SetParameters(const Double_t *params)
Definition: TF1.h:638
virtual Double_t GetParameter(Int_t ipar) const
Definition: TF1.h:506
1-D histogram with a double per channel (see TH1 documentation)}
Definition: TH1.h:614
To draw Mathematical Formula.
Definition: TLatex.h:18
Stopwatch class.
Definition: TStopwatch.h:28
void Start(Bool_t reset=kTRUE)
Start the stopwatch.
Definition: TStopwatch.cxx:58
void Print(Option_t *option="") const
Print the real and cpu time passed between the start and stop events.
Definition: TStopwatch.cxx:219
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:2311
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:1590
TFitResultPtr Fit(FitObject *h1, TF1 *f1, Foption_t &option, const ROOT::Math::MinimizerOptions &moption, const char *goption, ROOT::Fit::DataRange &range)
Definition: HFitImpl.cxx:134
Bool_t AreEqualAbs(Double_t af, Double_t bf, Double_t epsilon)
Definition: TMath.h:412
th1 Draw()
auto * t1
Definition: textangle.C:20