ROOT
Version v6.32
master
v6.34
v6.30
v6.28
v6.26
v6.24
v6.22
v6.20
v6.18
v6.16
v6.14
v6.12
v6.10
v6.08
v6.06
Reference Guide
►
ROOT
•
All
Classes
Namespaces
Files
Functions
Variables
Typedefs
Enumerations
Enumerator
Properties
Friends
Macros
Modules
Pages
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_code
8
## \macro_output
9
##
10
## \date January 2022
11
## \author Rahul Balasubramanian
12
13
import
ROOT
14
15
ROOT.gStyle.SetOptStat
(0)
16
ROOT.PyConfig.IgnoreCommandLineOptions
=
True
17
ROOT.gROOT.SetBatch
(
True
)
18
19
# Create functions
20
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
21
observablename =
"pTV"
22
obsvar =
ROOT.RooRealVar
(observablename,
"observable of pTV"
, 10, 600)
23
24
# Setup three EFT coefficient and constant SM modifier
25
kSM =
ROOT.RooRealVar
(
"kSM"
,
"sm modifier"
, 1.0)
26
cHq3 =
ROOT.RooRealVar
(
"cHq3"
,
"EFT modifier"
, -10.0, 10.0)
27
cHq3.setAttribute
(
"NewPhysics"
,
True
)
28
cHl3 =
ROOT.RooRealVar
(
"cHl3"
,
"EFT modifier"
, -10.0, 10.0)
29
cHl3.setAttribute
(
"NewPhysics"
,
True
)
30
cHDD =
ROOT.RooRealVar
(
"cHDD"
,
"EFT modifier"
, -10.0, 10.0)
31
cHDD.setAttribute
(
"NewPhysics"
,
True
)
32
33
# Inputs to setup config
34
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
35
infilename =
ROOT.gROOT.GetTutorialDir
().Data() +
"/roofit/input_histos_rf_lagrangianmorph.root"
36
par =
"cHq3"
37
samplelist = [
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
53
config =
ROOT.RooLagrangianMorphFunc.Config
()
54
config.fileName
= infilename
55
config.observableName
= observablename
56
config.folderNames
= samplelist
57
config.couplings.add
(cHq3)
58
config.couplings.add
(cHDD)
59
config.couplings.add
(cHl3)
60
config.couplings.add
(kSM)
61
62
63
# Create morphing function
64
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
65
66
morphfunc =
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
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
70
morphfunc.setParameter
(
"cHq3"
, 0.01)
71
morphfunc.setParameter
(
"cHl3"
, 1.0)
72
morphfunc.setParameter
(
"cHDD"
, 0.2)
73
74
pseudo_hist =
morphfunc.createTH1
(
"pseudo_hist"
)
75
pseudo_dh =
ROOT.RooDataHist
(
"pseudo_dh"
,
"pseudo_dh"
, [obsvar], pseudo_hist)
76
77
# reset parameters to zeros before fit
78
morphfunc.setParameter
(
"cHq3"
, 0.0)
79
morphfunc.setParameter
(
"cHl3"
, 0.0)
80
morphfunc.setParameter
(
"cHDD"
, 0.0)
81
82
# set error to set initial step size in fit
83
cHq3.setError
(0.1)
84
cHl3.setError
(0.1)
85
cHDD.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
91
model =
ROOT.RooWrapperPdf
(
"wrap_pdf"
,
"wrap_pdf"
, morphfunc)
92
fitres =
model.fitTo
(pseudo_dh, SumW2Error=
True
, Optimize=
False
, Save=
True
, PrintLevel=-1)
93
# run the fit
94
# Get the correlation matrix
95
hcorr =
fitres.correlationHist
()
96
97
# Extract postfit distribution and plot with initial histogram
98
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
99
100
postfit_hist =
morphfunc.createTH1
(
"morphing_postfit_hist"
)
101
postfit_dh =
ROOT.RooDataHist
(
"morphing_postfit_dh"
,
"morphing_postfit_dh"
, [obsvar], postfit_hist)
102
103
frame0 =
obsvar.frame
(Title=
"Input templates for p_{T}^{V}"
)
104
postfit_dh.plotOn
(
105
frame0,
106
Name=
"postfit_dist"
,
107
DrawOption=
"C"
,
108
LineColor=
"b"
,
109
DataError=
None
,
110
XErrorSize=0,
111
)
112
pseudo_dh.plotOn
(frame0, Name=
"input"
)
113
114
# Draw plots on canvas
115
# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
116
117
c1 =
ROOT.TCanvas
(
"fig3"
,
"fig3"
, 800, 400)
118
c1.Divide
(2, 1)
119
120
c1.cd
(1)
121
ROOT.gPad.SetLeftMargin
(0.15)
122
ROOT.gPad.SetRightMargin
(0.05)
123
124
model.paramOn
(frame0,
ROOT.RooFit.Layout
(0.50, 0.75, 0.9))
125
frame0.GetXaxis
().SetTitle(
"p_{T}^{V}"
)
126
frame0.Draw
()
127
128
c1.cd
(2)
129
ROOT.gPad.SetLeftMargin
(0.15)
130
ROOT.gPad.SetRightMargin
(0.15)
131
ROOT.gStyle.SetPaintTextFormat
(
"4.1f"
)
132
ROOT.gStyle.SetOptStat
(0)
133
hcorr.SetMarkerSize
(3.0)
134
hcorr.SetTitle
(
"correlation matrix"
)
135
hcorr.GetYaxis
().SetTitleOffset(1.4)
136
hcorr.GetYaxis
().SetLabelSize(0.1)
137
hcorr.GetXaxis
().SetLabelSize(0.1)
138
hcorr.GetYaxis
().SetBinLabel(1,
"c_{HDD}"
)
139
hcorr.GetYaxis
().SetBinLabel(2,
"c_{Hl^{(3)}}"
)
140
hcorr.GetYaxis
().SetBinLabel(3,
"c_{Hq^{(3)}}"
)
141
hcorr.GetXaxis
().SetBinLabel(3,
"c_{HDD}"
)
142
hcorr.GetXaxis
().SetBinLabel(2,
"c_{Hl^{(3)}}"
)
143
hcorr.GetXaxis
().SetBinLabel(1,
"c_{Hq^{(3)}}"
)
144
hcorr.GetYaxis
().SetTitleOffset(1.4)
145
hcorr.Draw
(
"colz text"
)
146
147
c1.SaveAs
(
"rf712_lagrangianmorphfit.png"
)
TRangeDynCast
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
Definition
TCollection.h:358
tutorials
roofit
rf712_lagrangianmorphfit.py
ROOT v6-32 - Reference Guide Generated on Wed Apr 2 2025 08:24:49 (GVA Time) using Doxygen 1.10.0