Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf615_simulation_based_inference.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit_main
3## \notebook
4## Use Simulation Based Inference (SBI) in RooFit.
5##
6## This tutorial shows how to use SBI in ROOT. As reference distribution we
7## choose a simple uniform distribution. The target distribution is chosen to
8## be gaussian with different mean values.
9## The classifier is trained to discriminate between the reference and target
10## distribution.
11## We see how the neural networks generalize to unknown mean values.
12##
13## We compare the approach of using the likelihood ratio trick to morphing.
14##
15## An introduction of SBI can be found in https://arxiv.org/pdf/2010.06439.
16##
17## A short recap:
18## The idea of SBI is to fit a surrogate model to the data, in order to really
19## learn the likelihood function instead of calculating it. Therefore, a classifier is trained to discriminate between
20## samples from a target distribution (here the Gaussian) $$x\sim p(x|\theta)$$ and a reference distribution (here the Uniform)
21## $$x\sim p_{ref}(x|\theta)$$.
22##
23## The output of the classifier $$\hat{s}(\theta)$$ is a monotonic function of the likelihood ration and can be turned into an estimate of the likelihood ratio
24## via $$\hat{r}(\theta)=\frac{1-\hat{s}(\theta)}{\hat{s}(\theta)}.$$
25## This is called the likelihood ratio trick.
26##
27## In the end we compare the negative logarithmic likelihoods of the learned, morphed and analytical likelihood with minuit and as a plot.
28##
29## \macro_image
30## \macro_code
31## \macro_output
32##
33## \date July 2024
34## \author Robin Syring
35
36import numpy as np
37import ROOT
38from sklearn.neural_network import MLPClassifier
39
40# The samples used for training the classifier in this tutorial / rescale for more accuracy
41n_samples = 10000
42
43# Kills warning messages
45
46
47# Morphing as a baseline
48def morphing(setting, workspace):
49 # Define binning for morphing
52
53 # Number of 'sampled' gaussians, if you change it, adjust the binning properly
54 n_grid = 5
55
56 for i in range(n_grid):
57 # Define the sampled gausians
58 mu_help = ROOT.RooRealVar(f"mu{i}", f"mu{i}", i)
59 help = ROOT.RooGaussian(f"g{i}", f"g{i}", x_var, mu_help, sigma)
60 workspace.Import(help, Silence=True)
61
62 # Fill the histograms
63 hist = workspace[f"g{i}"].generateBinned([x_var], n_samples)
64
65 # Make sure that every bin is filled and we don't get zero probability
66 for i_bin in range(hist.numEntries()):
67 hist.add(hist.get(i_bin), 1.0)
68
69 # Add the pdf to the workspace
70 workspace.Import(ROOT.RooHistPdf(f"histpdf{i}", f"histpdf{i}", [x_var], hist, 1), Silence=True)
71
72 # Add the pdf to the grid
73 grid.addPdf(workspace[f"histpdf{i}"], i)
74
75 # Create the morphing and add it to the workspace
76 morph_func = ROOT.RooMomentMorphFuncND("morph_func", "morph_func", [mu_var], [x_var], grid, setting)
77 morph = ROOT.RooWrapperPdf("morph", "morph", morph_func, True)
78 workspace.Import(morph, Silence=True)
79
80 # Uncomment to see input plots for the first dimension (you might need to increase the morphed samples)
81 # f1 = x_var.frame(Title="linear morphing;x;pdf", Range=(-4, 8))
82 # for i in range(n_grid):
83 # workspace[f"histpdf{i}"].plotOn(f1)
84 # workspace["morph"].plotOn(f1, LineColor="r")
85 # c0 = ROOT.TCanvas()
86 # f1.Draw()
87 # input() # Wait for user input to proceed
88
89
90# Class used in this case to demonstrate the use of SBI in Root
91class SBI:
92 # Initializing the class SBI
93 def __init__(self, workspace):
94 # Choose the hyperparameters for training the neural network
95 self.classifier = MLPClassifier(hidden_layer_sizes=(20, 20), max_iter=1000, random_state=42)
96 self.data_model = None
97 self.data_ref = None
98 self.X_train = None
99 self.y_train = None
100 self.workspace = workspace
101
102 # Defining the target / training data for different values of mean value mu
103 def model_data(self, model, x, mu, n_samples):
104 ws = self.workspace
105 data_test_model = []
106 samples_gaussian = ws[model].generate([ws[x], ws[mu]], n_samples).to_numpy()
107 self._training_mus = samples_gaussian[mu]
108 data_test_model.extend(samples_gaussian[x])
109
110 self.data_model = np.array(data_test_model).reshape(-1, 1)
111
112 # Generating samples for the reference distribution
113 def reference_data(self, model, x, n_samples):
114 ws = self.workspace
115 # Ensuring the normalization with generating as many reference data as target data
116 samples_uniform = ws[model].generate(ws[x], n_samples)
117 data_reference_model = np.array(
118 [samples_uniform.get(i).getRealValue("x") for i in range(samples_uniform.numEntries())]
119 )
120 self.data_ref = data_reference_model.reshape(-1, 1)
121
122 # Bringing the data in the right format for training
123 def preprocessing(self):
124 thetas_model = self._training_mus.reshape(-1, 1)
125 thetas_reference = self._training_mus.reshape(-1, 1)
126 thetas = np.concatenate((thetas_model, thetas_reference), axis=0)
127 X = np.concatenate([self.data_model, self.data_ref])
128 self.y_train = np.concatenate([np.ones(len(self.data_model)), np.zeros(len(self.data_ref))])
129 self.X_train = np.concatenate([X, thetas], axis=1)
130
131 # Train the classifier
132 def train_classifier(self):
133 self.classifier.fit(self.X_train, self.y_train)
134
135
136# Setting the training and toy data samples; the factor 5 to enable a fair comparison to morphing
137n_samples_train = n_samples * 5
138
139
140# Define the "observed" data in a workspace
141def build_ws(mu_observed, sigma):
142 # using a workspace for easier processing inside the class
143 ws = ROOT.RooWorkspace()
144 ws.factory(f"Gaussian::gauss(x[-5,15], mu[0,4], {sigma})")
145 ws.factory("Uniform::uniform(x)")
146 ws["mu"].setVal(mu_observed)
147 ws.Print("v")
148 obs_data = ws["gauss"].generate(ws["x"], 1000)
149 obs_data.SetName("obs_data")
150 ws.Import(obs_data, Silence=True)
151
152 return ws
153
154
155# The "observed" data
156mu_observed = 2.5
157sigma = 1.5
158workspace = build_ws(mu_observed, sigma)
159x_var = workspace["x"]
160mu_var = workspace["mu"]
161gauss = workspace["gauss"]
162uniform = workspace["uniform"]
163obs_data = workspace["obs_data"]
164
165# Training the model
166model = SBI(workspace)
167model.model_data("gauss", "x", "mu", n_samples_train)
168model.reference_data("uniform", "x", n_samples_train)
171sbi_model = model
172
173
174# Compute the likelihood ratio of the classifier for analysis purposes
175def learned_likelihood_ratio(x, mu):
176 n = max(len(x), len(mu))
177 X = np.zeros((n, 2))
178 X[:, 0] = x
179 X[:, 1] = mu
181 return prob / (1 - prob)
182
183
184# Compute the learned likelihood ratio
185llhr_learned = ROOT.RooFit.bindFunction("MyBinFunc", learned_likelihood_ratio, x_var, mu_var)
186
187# Compute the real likelihood ratio
188llhr_calc = ROOT.RooFormulaVar("llhr_calc", "x[0] / x[1]", [gauss, uniform])
189
190# Create the exact negative log likelihood functions for Gaussian model
191nll_gauss = gauss.createNLL(obs_data)
192ROOT.SetOwnership(nll_gauss, True)
193
194# Create the learned pdf and NLL sum based on the learned likelihood ratio
195pdf_learned = ROOT.RooWrapperPdf("learned_pdf", "learned_pdf", llhr_learned, True)
196
197nllr_learned = pdf_learned.createNLL(obs_data)
198ROOT.SetOwnership(nllr_learned, True)
199
200# Compute the morphed nll
202nll_morph = workspace["morph"].createNLL(obs_data)
203ROOT.SetOwnership(nll_morph, True)
204
205# Plot the negative logarithmic summed likelihood
206frame1 = mu_var.frame(Title="NLL of SBI vs. Morphing;mu;NLL", Range=(2.2, 2.8))
207nllr_learned.plotOn(frame1, LineColor="kP6Blue", ShiftToZero=True, Name="learned")
208nll_gauss.plotOn(frame1, LineColor="kP6Yellow", ShiftToZero=True, Name="gauss")
209ROOT.RooAbsReal.setEvalErrorLoggingMode("Ignore") # Silence some warnings
210nll_morph.plotOn(frame1, LineColor="kP6Red", ShiftToZero=True, Name="morphed")
212
213# Plot the likelihood functions
214frame2 = x_var.frame(Title="Likelihood ratio r(x|#mu=2.5);x;p_{gauss}/p_{uniform}")
215llhr_learned.plotOn(frame2, LineColor="kP6Blue", Name="learned_ratio")
216llhr_calc.plotOn(frame2, LineColor="kP6Yellow", Name="exact")
217
218# Write the plots into one canvas to show, or into separate canvases for saving.
219single_canvas = True
220
221c = ROOT.TCanvas("", "", 1200 if single_canvas else 600, 600)
222if single_canvas:
223 c.Divide(2)
224 c.cd(1)
226 frame1.GetYaxis().SetTitleOffset(1.8)
228
229legend1 = ROOT.TLegend(0.43, 0.63, 0.8, 0.87)
230legend1.SetFillColor("kWhite")
231legend1.SetLineColor("kWhite")
233legend1.AddEntry("learned", "learned (SBI)", "L")
234legend1.AddEntry("gauss", "true NLL", "L")
235legend1.AddEntry("morphed", "moment morphing", "L")
237
238if single_canvas:
239 c.cd(2)
241 frame2.GetYaxis().SetTitleOffset(1.8)
242else:
243 c.SaveAs("rf615_plot_1.png")
244 c = ROOT.TCanvas("", "", 600, 600)
245
247
248legend2 = ROOT.TLegend(0.53, 0.73, 0.87, 0.87)
249legend2.SetFillColor("kWhite")
250legend2.SetLineColor("kWhite")
252legend2.AddEntry("learned_ratio", "learned (SBI)", "L")
253legend2.AddEntry("exact", "true ratio", "L")
255
256if not single_canvas:
257 c.SaveAs("rf615_plot_2.png")
258
259# Compute the minimum via minuit and display the results
260for nll in [nll_gauss, nllr_learned, nll_morph]:
261 minimizer = ROOT.RooMinimizer(nll)
262 minimizer.setErrorLevel(0.5) # Adjust the error level in the minimization to work with likelihoods
264 minimizer.minimize("Minuit2")
265 result = minimizer.save()
266 ROOT.SetOwnership(result, True)
268
269del nll_morph
270del nllr_learned
271del nll_gauss
272del workspace
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 UChar_t len