High Level Factory: creation of a simple model
[#1] INFO:Minimization -- p.d.f. provides expected number of events, including extended term in likelihood.
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions will be evaluated in cache-and-track mode: (gauss,argus)
Minuit2Minimizer: Minimize with max-calls 2500 convergence for edm < 1 strategy 1
Minuit2Minimizer : Valid minimum - status = 0
FVAL = -18178.1839194288186
Edm = 7.6129471840253625e-05
Nfcn = 112
argpar = -22.9054 +/- 3.42143 (limited)
mean = 5.27987 +/- 0.000215769 (limited)
nbkg = 1612.74 +/- 44.6723 (limited)
nsig = 387.37 +/- 27.7674 (limited)
width = 0.00300965 +/- 0.000199023 (limited)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:Plotting -- RooAbsPdf::plotOn(sum) directly selected PDF components: (argus)
[#1] INFO:Plotting -- RooAbsPdf::plotOn(sum) indirectly selected PDF components: ()
import ROOT
card_name = "HLFavtoryexample.rs"
with open(card_name, "w") as f:
f.write("// The simplest card\n\n")
f.write("gauss = Gaussian(mes[5.20,5.30],mean[5.28,5.2,5.3],width[0.0027,0.001,1]);\n")
f.write("argus = ArgusBG(mes,5.291,argpar[-20,-100,-1]);\n")
f.write("sum = SUM(nsig[200,0,10000]*gauss,nbkg[800,0,10000]*argus);\n\n")
hlf = ROOT.RooStats.HLFactory("HLFavtoryexample", card_name, False)
w = hlf.GetWs()
mes = w.arg("mes")
sumpdf = w["sum"]
argus = w["argus"]
data = sumpdf.generate(mes, 2000)
sumpdf.fitTo(data)
mesframe = mes.frame()
data.plotOn(mesframe)
sumpdf.plotOn(mesframe)
sumpdf.plotOn(mesframe, Components=argus, LineStyle="--")
ROOT.gROOT.SetStyle("Plain")
c = ROOT.TCanvas()
mesframe.Draw()
c.SaveAs("rs601_HLFactoryexample.png")
- Date
- July 2022
- Authors
- Artem Busorgin, Danilo Piparo (C++ version)
Definition in file rs601_HLFactoryexample.py.