Special pdf's: using a product of an (acceptance) efficiency and a pdf as pdf
import ROOT
t = ROOT.RooRealVar("t", "t", 0, 5)
tau = ROOT.RooRealVar("tau", "tau", -1.54, -4, -0.1)
model = ROOT.RooExponential("model", "model", t, tau)
eff = ROOT.RooFormulaVar("eff", "0.5*(TMath::Erf((t-1)/0.5)+1)", [t])
modelEff = ROOT.RooEffProd("modelEff", "model with efficiency", model, eff)
frame1 = t.frame(Title="Efficiency")
eff.plotOn(frame1, LineColor="r")
frame2 = t.frame(Title="Pdf with and without efficiency")
model.plotOn(frame2, LineStyle="--")
modelEff.plotOn(frame2)
data = modelEff.generate({t}, 10000)
modelEff.fitTo(data, PrintLevel=-1)
frame3 = t.frame(Title="Fitted pdf with efficiency")
data.plotOn(frame3)
modelEff.plotOn(frame3)
c = ROOT.TCanvas("rf703_effpdfprod", "rf703_effpdfprod", 1200, 400)
c.Divide(3)
c.cd(1)
ROOT.gPad.SetLeftMargin(0.15)
frame1.GetYaxis().SetTitleOffset(1.4)
frame1.Draw()
c.cd(2)
ROOT.gPad.SetLeftMargin(0.15)
frame2.GetYaxis().SetTitleOffset(1.6)
frame2.Draw()
c.cd(3)
ROOT.gPad.SetLeftMargin(0.15)
frame3.GetYaxis().SetTitleOffset(1.6)
frame3.Draw()
c.SaveAs("rf703_effpdfprod.png")
[#1] INFO:NumericIntegration -- RooRealIntegral::init(modelEff_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:NumericIntegration -- RooRealIntegral::init(modelEff_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: activating const optimization
[#1] INFO:Minimization -- The following expressions have been identified as constant and will be precalculated and cached: (eff)
[#1] INFO:Minimization -- RooAbsMinimizerFcn::setOptimizeConst: deactivating const optimization
[#1] INFO:NumericIntegration -- RooRealIntegral::init(modelEff_Int[t]) using numeric integrator RooIntegrator1D to calculate Int(t)
- Date
- February 2018
- Authors
- Clemens Lange, Wouter Verkerke (C++ version)
Definition in file rf703_effpdfprod.py.