Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf712_lagrangianmorphfit.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -js
4## Performing a simple fit with RooLagrangianMorphFunc
5##
6## \macro_image
7## \macro_output
8## \macro_code
9##
10## \date January 2022
11## \authors Rahul Balasubramanian
12
13import ROOT
14
15ROOT.gStyle.SetOptStat(0)
16ROOT.PyConfig.IgnoreCommandLineOptions = True
17ROOT.gROOT.SetBatch(True)
18
19# Create functions
20# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
21observablename = "pTV"
22obsvar = ROOT.RooRealVar(observablename, "observable of pTV", 10, 600)
23
24# Setup three EFT coefficent and constant SM modifier
25kSM = ROOT.RooRealVar("kSM", "sm modifier", 1.0)
26cHq3 = ROOT.RooRealVar("cHq3", "EFT modifier", -10.0, 10.0)
27cHq3.setAttribute("NewPhysics", True)
28cHl3 = ROOT.RooRealVar("cHl3", "EFT modifier", -10.0, 10.0)
29cHl3.setAttribute("NewPhysics", True)
30cHDD = ROOT.RooRealVar("cHDD", "EFT modifier", -10.0, 10.0)
31cHDD.setAttribute("NewPhysics", True)
32
33# Inputs to setup config
34# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
35infilename = ROOT.gROOT.GetTutorialDir().Data() + "/roofit/input_histos_rf_lagrangianmorph.root"
36par = "cHq3"
37samplelist = [
38 "SM_NPsq0",
39 "cHq3_NPsq1",
40 "cHq3_NPsq2",
41 "cHl3_NPsq1",
42 "cHl3_NPsq2",
43 "cHDD_NPsq1",
44 "cHDD_NPsq2",
45 "cHl3_cHDD_NPsq2",
46 "cHq3_cHDD_NPsq2",
47 "cHl3_cHq3_NPsq2",
48]
49
50# Set Config
51# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
52
53config = ROOT.RooLagrangianMorphFunc.Config()
54config.fileName = infilename
55config.observableName = observablename
56config.folderNames = samplelist
57config.couplings.add(cHq3)
58config.couplings.add(cHDD)
59config.couplings.add(cHl3)
60config.couplings.add(kSM)
61
62
63# Create morphing function
64# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
65
66morphfunc = ROOT.RooLagrangianMorphFunc("morphfunc", "morphed dist. of pTV", config)
67
68# Create pseudo data histogram to fit at cHq3 = 0.01, cHl3 = 1.0, cHDD = 0.2
69# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
70morphfunc.setParameter("cHq3", 0.01)
71morphfunc.setParameter("cHl3", 1.0)
72morphfunc.setParameter("cHDD", 0.2)
73
74pseudo_hist = morphfunc.createTH1("pseudo_hist")
75pseudo_dh = ROOT.RooDataHist("pseudo_dh", "pseudo_dh", [obsvar], pseudo_hist)
76
77# reset parameters to zeros before fit
78morphfunc.setParameter("cHq3", 0.0)
79morphfunc.setParameter("cHl3", 0.0)
80morphfunc.setParameter("cHDD", 0.0)
81
82# set error to set initial step size in fit
83cHq3.setError(0.1)
84cHl3.setError(0.1)
85cHDD.setError(0.1)
86
87# Wrap pdf on morphfunc and fit to data histogram
88# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
89
90# wrapper pdf to normalise morphing function to a morphing pdf
91model = ROOT.RooWrapperPdf("wrap_pdf", "wrap_pdf", morphfunc)
92fitres = model.fitTo(pseudo_dh, SumW2Error=True, Optimize=False, Save=True)
93# run the fit
94# Get the correlation matrix
95hcorr = fitres.correlationHist()
96
97# Extract postfit distribution and plot with initial histogram
98# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
99
100postfit_hist = morphfunc.createTH1("morphing_postfit_hist")
101postfit_dh = ROOT.RooDataHist("morphing_postfit_dh", "morphing_postfit_dh", [obsvar], postfit_hist)
102
103frame0 = obsvar.frame(Title="Input templates for p_{T}^{V}")
104postfit_dh.plotOn(
105 frame0,
106 Name="postfit_dist",
107 DrawOption="C",
108 LineColor="b",
109 DataError=None,
110 XErrorSize=0,
111)
112pseudo_dh.plotOn(frame0, Name="input")
113
114# Draw plots on canvas
115# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
116
117c1 = ROOT.TCanvas("fig3", "fig3", 800, 400)
118c1.Divide(2, 1)
119
120c1.cd(1)
121ROOT.gPad.SetLeftMargin(0.15)
122ROOT.gPad.SetRightMargin(0.05)
123
124model.paramOn(frame0, ROOT.RooFit.Layout(0.50, 0.75, 0.9))
125frame0.GetXaxis().SetTitle("p_{T}^{V}")
126frame0.Draw()
127
128c1.cd(2)
129ROOT.gPad.SetLeftMargin(0.15)
130ROOT.gPad.SetRightMargin(0.15)
131ROOT.gStyle.SetPaintTextFormat("4.1f")
132ROOT.gStyle.SetOptStat(0)
133hcorr.SetMarkerSize(3.0)
134hcorr.SetTitle("correlation matrix")
135hcorr.GetYaxis().SetTitleOffset(1.4)
136hcorr.GetYaxis().SetLabelSize(0.1)
137hcorr.GetXaxis().SetLabelSize(0.1)
138hcorr.GetYaxis().SetBinLabel(1, "c_{HDD}")
139hcorr.GetYaxis().SetBinLabel(2, "c_{Hl^{(3)}}")
140hcorr.GetYaxis().SetBinLabel(3, "c_{Hq^{(3)}}")
141hcorr.GetXaxis().SetBinLabel(3, "c_{HDD}")
142hcorr.GetXaxis().SetBinLabel(2, "c_{Hl^{(3)}}")
143hcorr.GetXaxis().SetBinLabel(1, "c_{Hq^{(3)}}")
144hcorr.GetYaxis().SetTitleOffset(1.4)
145hcorr.Draw("colz text")
146
147c1.SaveAs("rf712_lagrangianmorphfit.png")