Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf616_morphing.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook
4## Use Morphing in RooFit.
5##
6## This tutorial shows how to use template morphing inside RooFit. As input we have several
7## Gaussian distributions. The output is one gaussian, with a specific mean value.
8## Since likelihoods are often used within the framework of morphing, we provide a
9## way to estimate the negative log likelihood (nll).
10##
11## Based on example of Kyle Cranmer https://gist.github.com/cranmer/46fff8d22015e5a26619.
12##
13## \macro_image
14## \macro_code
15## \macro_output
16##
17## \date August 2024
18## \author Robin Syring
19
20
21import ROOT
22
23# Number of samples to fill the histograms
24n_samples = 1000
25
26
27# Kills warning messages
28ROOT.RooMsgService.instance().setGlobalKillBelow(ROOT.RooFit.WARNING)
29
30
31# morphing as a baseline
32def morphing(setting):
33 # set up a frame for plotting
34 frame1 = x_var.frame()
35
36 # define binning for morphing
37 bin_mu_x = ROOT.RooBinning(4, 0.0, 4.0)
38 grid = ROOT.RooMomentMorphFuncND.Grid(bin_mu_x)
39 x_var.setBins(50)
40
41 # number of 'sampled' Gaussians, if you change it, adjust the binning properly
42 for i in range(5):
43 # Create the sampled Gaussian
44 workspace.factory(f"Gaussian::g{i}(x, mu{i}[{i}], {sigma})".format(i=i))
45
46 # Fill the histograms
47 hist = workspace[f"g{i}"].generateBinned([x_var], n_samples * 100)
48 # Make sure that every bin is filled and we don't get zero probability
49 for i_bin in range(hist.numEntries()):
50 hist.add(hist.get(i_bin), 1.0)
51
52 # Add the pdf to the workspace, the inOrder of 1 is necessary for calculation of the nll
53 # Adjust it to 0 to see binning
54 workspace.Import(ROOT.RooHistPdf(f"histpdf{i}", f"histpdf{i}", [x_var], hist, intOrder=1))
55
56 # Add the pdf to the grid and to the plot
57 grid.addPdf(workspace[f"histpdf{i}"], int(i))
58 workspace[f"histpdf{i}"].plotOn(frame1)
59
60 # Create the morphing and add it to the workspace
61 morph_func = ROOT.RooMomentMorphFuncND("morph_func", "morph_func", [mu_var], [x_var], grid, setting)
62
63 # Normalizes the morphed object to be a pdf, set it false to prevent warning messages and gain computational speed up
64 morph_func.setPdfMode()
65
66 # Creating the morphed pdf
67 morph = ROOT.RooWrapperPdf("morph", "morph", morph_func, True)
68 workspace.Import(morph)
69 workspace["morph"].plotOn(frame1, LineColor="r")
70
71 return frame1
72
73
74# Define the "observed" data in a workspade
75def build_ws(mu_observed, sigma):
76 ws = ROOT.RooWorkspace()
77 ws.factory(f"Gaussian::gauss(x[-5,15], mu[{mu_observed},0,4], {sigma})".format(mu_observed=mu_observed))
78
79 return ws
80
81
82# The "observed" data
83mu_observed = 2.5
84sigma = 1.5
85workspace = build_ws(mu_observed, sigma)
86x_var = workspace["x"]
87mu_var = workspace["mu"]
88gauss = workspace.pdf("gauss")
89obs_data = gauss.generate(x_var, n_samples)
90
91
92# Create the exact negative log likelihood functions for Gaussian model
93nll_gauss = gauss.createNLL(obs_data)
94
95# Compute the morphed nll
96frame1 = morphing(ROOT.RooMomentMorphFuncND.Linear)
97
98# TODO: Fix RooAddPdf::fixCoefNormalization(nset) warnings with new CPU backend
99nll_morph = workspace["morph"].createNLL(obs_data, EvalBackend="legacy")
100
101# Plot the negative logarithmic summed likelihood
102frame2 = mu_var.frame(Title="Negative log Likelihood")
103nll_gauss.plotOn(frame2, LineColor="b", ShiftToZero=True, Name="gauss")
104nll_morph.plotOn(frame2, LineColor="r", ShiftToZero=True, Name="morphed")
105
106c = ROOT.TCanvas("rf616_morphing", "rf616_morphing", 800, 400)
107c.Divide(2)
108c.cd(1)
109ROOT.gPad.SetLeftMargin(0.15)
110frame1.GetYaxis().SetTitleOffset(1.8)
111frame1.Draw()
112c.cd(2)
113ROOT.gPad.SetLeftMargin(0.15)
114frame2.GetYaxis().SetTitleOffset(1.8)
115frame2.Draw()
116
117
118# Compute the minimum via minuit and display the results
119for nll in [nll_gauss, nll_morph]:
120 minimizer = ROOT.RooMinimizer(nll)
121 minimizer.setPrintLevel(-1)
122 minimizer.minimize("Minuit2")
123 result = minimizer.save()
124 result.Print()
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void char Point_t Rectangle_t WindowAttributes_t Float_t Float_t Float_t Int_t Int_t UInt_t UInt_t Rectangle_t Int_t Int_t Window_t TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t format