Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf711_lagrangianmorph.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -js
4## Morphing effective field theory distributions with RooLagrangianMorphFunc.
5## A morphing function as a function of one coefficient is setup and can be used
6## to obtain the distribution for any value of the coefficient.
7##
8## \macro_image
9## \macro_output
10## \macro_code
11##
12## \date January 2022
13## \authors Rahul Balasubramanian
14
15import ROOT
16
17ROOT.gStyle.SetOptStat(0)
18ROOT.PyConfig.IgnoreCommandLineOptions = True
19ROOT.gROOT.SetBatch(True)
20
21# Create functions
22# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
23observablename = "pTV"
24
25# Setup observable that is to be morphed
26obsvar = ROOT.RooRealVar(observablename, "p_{T}^{V}", 10, 600)
27
28# Setup two couplings that enters the morphing function
29# kSM -> SM coupling set to constant (1)
30# cHq3 -> EFT parameter with NewPhysics attribute set to true
31kSM = ROOT.RooRealVar("kSM", "sm modifier", 1.0)
32cHq3 = ROOT.RooRealVar("cHq3", "EFT modifier", 0.0, 1.0)
33cHq3.setAttribute("NewPhysics", True)
34
35# Inputs to setup config
36# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
37infilename = ROOT.gROOT.GetTutorialDir().Data() + "/roofit/input_histos_rf_lagrangianmorph.root"
38par = "cHq3"
39samplelist = ["SM_NPsq0", "cHq3_NPsq1", "cHq3_NPsq2"]
40
41# Set Config
42# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
43
44config = ROOT.RooLagrangianMorphFunc.Config()
45config.fileName = infilename
46config.observableName = observablename
47config.folderNames = samplelist
48config.couplings.add(cHq3)
49config.couplings.add(kSM)
50
51
52# Create morphing function
53# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
54
55morphfunc = ROOT.RooLagrangianMorphFunc("morphfunc", "morphed dist. of pTV", config)
56
57# Get morphed distribution at cHq3 = 0.01, 0.5
58# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
59morphfunc.setParameter("cHq3", 0.01)
60morph_hist_0p01 = morphfunc.createTH1("morph_cHq3=0.01")
61morphfunc.setParameter("cHq3", 0.25)
62morph_hist_0p25 = morphfunc.createTH1("morph_cHq3=0.25")
63morphfunc.setParameter("cHq3", 0.5)
64morph_hist_0p5 = morphfunc.createTH1("morph_cHq3=0.5")
65morph_datahist_0p01 = ROOT.RooDataHist("morph_dh_cHq3=0.01", "", [obsvar], morph_hist_0p01)
66morph_datahist_0p25 = ROOT.RooDataHist("morph_dh_cHq3=0.25", "", [obsvar], morph_hist_0p25)
67morph_datahist_0p5 = ROOT.RooDataHist("morph_dh_cHq3=0.5", "", [obsvar], morph_hist_0p5)
68
69# Extract input templates for plotting
70# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
71
72input_hists = {sample: ROOT.TFile.Open(infilename).Get(sample).FindObject(observablename) for sample in samplelist}
73input_datahists = {
74 sample: ROOT.RooDataHist("dh_" + sample, "dh_" + sample, [obsvar], input_hists[sample]) for sample in samplelist
75}
76
77# Plot input templates
78# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
79
80frame0 = obsvar.frame(Title="Input templates for p_{T}^{V}")
81for sample, color in zip(samplelist, "krb"):
82 input_datahists[sample].plotOn(frame0, Name=sample, LineColor=color, MarkerColor=color, MarkerSize=1)
83
84# Plot morphed templates for cHq3=0.01,0.25,0.5
85# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
86
87frame1 = obsvar.frame(Title="Morphed templates for selected values")
88plot_args = dict(
89 DrawOption="C",
90 DataError=None,
91 XErrorSize=0,
92)
93morph_datahist_0p01.plotOn(frame1, Name="morph_dh_cHq3=0.01", LineColor="kGreen", **plot_args)
94morph_datahist_0p25.plotOn(frame1, Name="morph_dh_cHq3=0.25", LineColor="kGreen+1", **plot_args)
95morph_datahist_0p5.plotOn(frame1, Name="morph_dh_cHq3=0.5", LineColor="kGreen+2", **plot_args)
96
97# Create wrapped pdf to generate 2D dataset of cHq3 as a function of pTV
98# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
99
100model = ROOT.RooWrapperPdf("wrap_pdf", "wrap_pdf", morphfunc)
101data = model.generate({cHq3, obsvar}, 1000000)
102hh_data = ROOT.RooAbsData.createHistogram(data, "x,y", obsvar, Binning=20, YVar=dict(var=cHq3, Binning=50))
103hh_data.SetTitle("Morphing prediction")
104
105# Draw plots on canvas
106# -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -- -
107
108c1 = ROOT.TCanvas("fig3", "fig3", 1200, 400)
109c1.Divide(3, 1)
110
111c1.cd(1)
112ROOT.gPad.SetLeftMargin(0.15)
113ROOT.gPad.SetRightMargin(0.05)
114
115frame0.Draw()
116leg1 = ROOT.TLegend(0.55, 0.65, 0.94, 0.87)
117leg1.SetTextSize(0.04)
118leg1.SetFillColor(ROOT.kWhite)
119leg1.SetLineColor(ROOT.kWhite)
120leg1.AddEntry("SM_NPsq0", "SM", "LP")
121leg1.AddEntry(0, "", "")
122leg1.AddEntry("cHq3_NPsq1", "c_{Hq^{(3)}}=1.0 at #Lambda^{-2}", "LP")
123leg1.AddEntry(0, "", "")
124leg1.AddEntry("cHq3_NPsq2", "c_{Hq^{(3)}}=1.0 at #Lambda^{-4}", "LP")
125leg1.Draw()
126
127c1.cd(2)
128ROOT.gPad.SetLeftMargin(0.15)
129ROOT.gPad.SetRightMargin(0.05)
130
131frame1.Draw()
132
133leg2 = ROOT.TLegend(0.62, 0.65, 0.94, 0.87)
134leg2.SetTextSize(0.04)
135leg2.SetFillColor(ROOT.kWhite)
136leg2.SetLineColor(ROOT.kWhite)
137
138leg2.AddEntry("morph_dh_cHq3=0.01", "c_{Hq^{(3)}}=0.01", "L")
139leg2.AddEntry(0, "", "")
140leg2.AddEntry("morph_dh_cHq3=0.025", "c_{Hq^{(3)}}=0.025", "L")
141leg2.AddEntry(0, "", "")
142leg2.AddEntry("morph_dh_cHq3=0.5", "c_{Hq^{(3)}}=0.5", "L")
143leg2.AddEntry(0, "", "")
144leg2.Draw()
145
146c1.cd(3)
147ROOT.gPad.SetLeftMargin(0.12)
148ROOT.gPad.SetRightMargin(0.18)
149ROOT.gStyle.SetNumberContours(255)
150ROOT.gStyle.SetPalette(ROOT.kGreyScale)
151ROOT.gStyle.SetOptStat(0)
152ROOT.TColor.InvertPalette()
153ROOT.gPad.SetLogz()
154hh_data.GetYaxis().SetTitle("c_{Hq^{(3)}}")
155hh_data.GetYaxis().SetRangeUser(0, 0.5)
156hh_data.GetZaxis().SetTitleOffset(1.8)
157hh_data.Draw("COLZ")
158c1.SaveAs("rf711_lagrangianmorph.png")