Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
fitNormSum.C File Reference

Detailed Description

View in nbviewer Open in SWAN
Tutorial for normalized sum of two functions Here: a background exponential and a crystalball function Parameters can be set:

  1. with the TF1 object before adding the function (for 3) and 4))
  2. with the TF1NormSum object (first two are the coefficients, then the non constant parameters)
  3. with the TF1 object after adding the function

Sum can be constructed by:

  1. by a string containing the names of the functions and/or the coefficient in front
  2. by a string containg formulas like expo, gaus...
  3. by the list of functions and coefficients (which are 1 by default)
  4. by a std::vector for functions and coefficients
Time to generate 1050000 events: Real time 0:00:00, CP time 0.160
****************************************
Minimizer is Minuit2 / Migrad
Chi2 = 1018.73
NDf = 993
Edm = 9.65559e-06
NCalls = 233
NSignal = 50082 +/- 1231.21
NBackground = 998899 +/- 1569.86
Mean = 2.99896 +/- 0.0022426
Sigma = 0.297871 +/- 0.00230279
Alpha = 2.12493 +/- 0.1368
N = 1.1562 +/- 0.468136
Slope = -0.300341 +/- 0.000644187
Time to fit using ROOT TF1Normsum: Real time 0:00:00, CP time 0.160
#include <TCanvas.h>
#include <TF1.h>
#include <TF1NormSum.h>
#include <TFitResult.h>
#include <TH1.h>
#include <TLatex.h>
#include <TMath.h>
#include <TStopwatch.h>
#include <TStyle.h>
void fitNormSum()
{
const int nsig = 5.E4;
const int nbkg = 1.e6;
int nEvents = nsig + nbkg;
int nBins = 1e3;
double signal_mean = 3;
TF1 *f_cb = new TF1("MyCrystalBall", "crystalball", -5., 5.);
TF1 *f_exp = new TF1("MyExponential", "expo", -5., 5.);
// I.:
f_exp->SetParameters(1., -0.3);
f_cb->SetParameters(1, signal_mean, 0.3, 2, 1.5);
// CONSTRUCTION OF THE TF1NORMSUM OBJECT ........................................
// 1) :
TF1NormSum *fnorm_exp_cb = new TF1NormSum(f_cb, f_exp, nsig, nbkg);
// 4) :
TF1 *f_sum = new TF1("fsum", *fnorm_exp_cb, -5., 5., fnorm_exp_cb->GetNpar());
// III.:
f_sum->SetParameters(fnorm_exp_cb->GetParameters().data());
f_sum->SetParName(1, "NBackground");
f_sum->SetParName(0, "NSignal");
for (int i = 2; i < f_sum->GetNpar(); ++i)
f_sum->SetParName(i, fnorm_exp_cb->GetParName(i));
// GENERATE HISTOGRAM TO FIT ..............................................................
w.Start();
TH1D *h_sum = new TH1D("h_ExpCB", "Exponential Bkg + CrystalBall function", nBins, -5., 5.);
h_sum->FillRandom("fsum", nEvents);
printf("Time to generate %d events: ", nEvents);
w.Print();
// need to scale histogram with width since we are fitting a density
h_sum->Sumw2();
h_sum->Scale(1., "width");
// fit - use Minuit2 if available
new TCanvas("Fit", "Fit", 800, 1000);
// do a least-square fit of the spectrum
auto result = h_sum->Fit("fsum", "SQ");
result->Print();
h_sum->Draw();
printf("Time to fit using ROOT TF1Normsum: ");
w.Print();
// test if parameters are fine
std::vector<double> pref = {nsig, nbkg, signal_mean};
for (unsigned int i = 0; i < pref.size(); ++i) {
if (!TMath::AreEqualAbs(pref[i], f_sum->GetParameter(i), f_sum->GetParError(i) * 10.))
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));
}
// add parameters
auto t1 = new TLatex(
-2.5, 300000, TString::Format("%s = %8.0f #pm %4.0f", "NSignal", f_sum->GetParameter(0), f_sum->GetParError(0)));
auto t2 = new TLatex(
-2.5, 270000, TString::Format("%s = %8.0f #pm %4.0f", "Nbackgr", f_sum->GetParameter(1), f_sum->GetParError(1)));
t1->Draw();
t2->Draw();
}
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:436
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:1930
virtual Int_t GetNpar() const
Definition TF1.h:509
virtual const char * GetParName(Int_t ipar) const
Definition TF1.h:557
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:677
virtual Double_t GetParameter(Int_t ipar) const
Definition TF1.h:540
1-D histogram with a double per channel (see TH1 documentation)
Definition TH1.h:670
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:6604
virtual void Sumw2(Bool_t flag=kTRUE)
Create structure to store sum of squares of weights.
Definition TH1.cxx:9020
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:1640
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
Author
Lorenzo Moneta

Definition in file fitNormSum.C.