#!/usr/bin/env python # coding: utf-8 # # rf306_condpereventerrors # Multidimensional models: complete example with use of conditional pdf with per-event errors # # # # # **Author:** Clemens Lange, Wouter Verkerke (C++ version) # This notebook tutorial was automatically generated with ROOTBOOK-izer from the macro found in the ROOT repository on Thursday, July 03, 2025 at 03:45 PM. # In[1]: import ROOT # B-physics pdf with per-event Gaussian resolution # ---------------------------------------------------------------------------------------------- # Observables # In[2]: dt = ROOT.RooRealVar("dt", "dt", -10, 10) dterr = ROOT.RooRealVar("dterr", "per-event error on dt", 0.01, 10) # Build a gaussian resolution model scaled by the per-error = # gauss(dt,bias,sigma*dterr) # In[3]: bias = ROOT.RooRealVar("bias", "bias", 0, -10, 10) sigma = ROOT.RooRealVar("sigma", "per-event error scale factor", 1, 0.1, 10) gm = ROOT.RooGaussModel("gm1", "gauss model scaled bt per-event error", dt, bias, sigma, dterr) # Construct decay(dt) (x) gauss1(dt|dterr) # In[4]: tau = ROOT.RooRealVar("tau", "tau", 1.548) decay_gm = ROOT.RooDecay("decay_gm", "decay", dt, tau, gm, type="DoubleSided") # Construct fake 'external' data with per-event error # ------------------------------------------------------------------------------------------------------ # Use landau pdf to get somewhat realistic distribution with long tail # In[5]: pdfDtErr = ROOT.RooLandau("pdfDtErr", "pdfDtErr", dterr, 1.0, 0.25) expDataDterr = pdfDtErr.generate({dterr}, 10000) # Sample data from conditional decay_gm(dt|dterr) # --------------------------------------------------------------------------------------------- # Specify external dataset with dterr values to use decay_dm as # conditional pdf # In[6]: data = decay_gm.generate({dt}, ProtoData=expDataDterr) # Fit conditional decay_dm(dt|dterr) # --------------------------------------------------------------------- # Specify dterr as conditional observable # In[7]: decay_gm.fitTo(data, ConditionalObservables={dterr}, PrintLevel=-1) # Plot conditional decay_dm(dt|dterr) # --------------------------------------------------------------------- # Make two-dimensional plot of conditional pdf in (dt,dterr) # In[8]: hh_decay = decay_gm.createHistogram("hh_decay", dt, Binning=50, YVar=dict(var=dterr, Binning=50)) hh_decay.SetLineColor("kBlue") # Plot decay_gm(dt|dterr) at various values of dterr # In[9]: frame = dt.frame(Title="Slices of decay(dt|dterr) at various dterr") for ibin in range(0, 100, 20): dterr.setBin(ibin) decay_gm.plotOn(frame, Normalization=5.0) # Make projection of data an dt # In[10]: frame2 = dt.frame(Title="Projection of decay(dt|dterr) on dt") data.plotOn(frame2) # Make projection of decay(dt|dterr) on dt. # # Instead of integrating out dterr, a weighted average of curves # at values dterr_i as given in the external dataset. # (The kTRUE argument bins the data before projection to speed up the process) # In[11]: decay_gm.plotOn(frame2, ProjWData=(expDataDterr, True)) # Draw all frames on canvas # In[12]: c = ROOT.TCanvas("rf306_condpereventerrors", "rf306_condperventerrors", 1200, 400) c.Divide(3) c.cd(1) ROOT.gPad.SetLeftMargin(0.20) hh_decay.GetZaxis().SetTitleOffset(2.5) hh_decay.Draw("surf") c.cd(2) ROOT.gPad.SetLeftMargin(0.15) frame.GetYaxis().SetTitleOffset(1.6) frame.Draw() c.cd(3) ROOT.gPad.SetLeftMargin(0.15) frame2.GetYaxis().SetTitleOffset(1.6) frame2.Draw() c.SaveAs("rf306_condpereventerrors.png") # Draw all canvases # In[13]: from ROOT import gROOT gROOT.GetListOfCanvases().Draw()