Tutorial for normalized sum of two functions Here: a background exponential and a crystalball function Parameters can be set:
Sum can be constructed by:
Author: Lorenzo Moneta
This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Monday, March 27, 2023 at 09:47 AM.
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 ........................................
TF1NormSum *fnorm_exp_cb = new TF1NormSum(f_cb, f_exp, nsig, nbkg);
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 ..............................................................
TStopwatch w;
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();
Time to generate 1050000 events: Real time 0:00:00, CP time 0.300
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
ROOT::Math::MinimizerOptions::SetDefaultMinimizer("Minuit2");
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();
**************************************** 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.300
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));
}
gStyle->SetOptStat(0);
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();
Draw all canvases
gROOT->GetListOfCanvases()->Draw()