Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
rf613_global_observables.py
Go to the documentation of this file.
1## \file
2## \ingroup tutorial_roofit
3## \notebook -js
4## This tutorial explains the concept of global observables in RooFit, and
5## showcases how their values can be stored either in the model or in the
6## dataset.
7##
8## ### Introduction
9##
10## Note: in this tutorial, we are multiplying the likelihood with an additional
11## likelihood to constrain the parameters with auxiliary measurements. This is
12## different from the `rf604_constraints` tutorial, where the likelihood is
13## multiplied with a Bayesian prior to constrain the parameters.
14##
15##
16## With RooFit, you usually optimize some model parameters `p` to maximize the
17## likelihood `L` given the per-event or per-bin observations `x`:
18##
19## \f[ L( x | p ) \f]
20##
21## Often, the parameters are constrained with some prior likelihood `C`, which
22## doesn't depend on the observables `x`:
23##
24## \f[ L'( x | p ) = L( x | p ) * C( p ) \f]
25##
26## Usually, these constraint terms depend on some auxiliary measurements of
27## other observables `g`. The constraint term is then the likelihood of the
28## so-called global observables:
29##
30## \f[ L'( x | p ) = L( x | p ) * C( g | p ) \f]
31##
32## For example, think of a model where the true luminosity `lumi` is a
33## nuisance parameter that is constrained by an auxiliary measurement
34## `lumi_obs` with uncertainty `lumi_obs_sigma`:
35##
36## \f[ L'(data | mu, lumi) = L(data | mu, lumi) * \text{Gauss}(lumi_obs | lumi, lumi_obs_sigma) \f]
37##
38## As a Gaussian is symmetric under exchange of the observable and the mean
39## parameter, you can also sometimes find this equivalent but less conventional
40## formulation for Gaussian constraints:
41##
42## \f[ L'(data | mu, lumi) = L(data | mu, lumi) * \text{Gauss}(lumi | lumi_obs, lumi_obs_sigma) \f]
43##
44## If you wanted to constrain a parameter that represents event counts, you
45## would use a Poissonian constraint, e.g.:
46##
47## \f[ L'(data | mu, count) = L(data | mu, count) * \text{Poisson}(count_obs | count) \f]
48##
49## Unlike a Gaussian, a Poissonian is not symmetric under exchange of the
50## observable and the parameter, so here you need to be more careful to follow
51## the global observable prescription correctly.
52##
53## \macro_code
54## \macro_output
55##
56## \date January 2022
57## \author Jonas Rembser
58
59
60import ROOT
61
62# Silence info output for this tutorial
63ROOT.RooMsgService.instance().getStream(1).removeTopic(ROOT.RooFit.Minimization)
64ROOT.RooMsgService.instance().getStream(1).removeTopic(ROOT.RooFit.Fitting)
65
66# Setting up the model and creating toy dataset
67# ---------------------------------------------
68
69# l'(x | mu, sigma) = l(x | mu, sigma) * Gauss(mu_obs | mu, 0.2)
70
71# event observables
72x = ROOT.RooRealVar("x", "x", -10, 10)
73
74# parameters
75mu = ROOT.RooRealVar("mu", "mu", 0.0, -10, 10)
76sigma = ROOT.RooRealVar("sigma", "sigma", 1.0, 0.1, 2.0)
77
78# Gaussian model for event observables
79gauss = ROOT.RooGaussian("gauss", "gauss", x, mu, sigma)
80
81# global observables (which are not parameters so they are constant)
82mu_obs = ROOT.RooRealVar("mu_obs", "mu_obs", 1.0, -10, 10)
83mu_obs.setConstant()
84# note: alternatively, one can create a constant with default limits using `RooRealVar("mu_obs", "mu_obs", 1.0)`
85
86# constraint pdf
87constraint = ROOT.RooGaussian("constraint", "constraint", mu_obs, mu, 0.1)
88
89# full pdf including constraint pdf
90model = ROOT.RooProdPdf("model", "model", [gauss, constraint])
91
92# Generating toy data with randomized global observables
93# ------------------------------------------------------
94
95# For most toy-based statistical procedures, it is necessary to also
96# randomize the global observable when generating toy datasets.
97#
98# To that end, let's generate a single event from the model and take the
99# global observable value (the same is done in the RooStats:ToyMCSampler
100# class):
101
102dataGlob = model.generate({mu_obs}, 1)
103
104# Next, we temporarily set the value of `mu_obs` to the randomized value for
105# generating our toy dataset:
106mu_obs_orig_val = mu_obs.getVal()
107
108ROOT.RooArgSet(mu_obs).assign(dataGlob.get(0))
109
110# Actually generate the toy dataset. We don't generate too many events,
111# otherwise, the constraint will not have much weight in the fit and the result
112# looks like it's unaffected by it.
113data = model.generate({x}, 50)
114
115# When fitting the toy dataset, it is important to set the global
116# observables in the fit to the values that were used to generate the toy
117# dataset. To facilitate the bookkeeping of global observable values, you
118# can attach a snapshot with the current global observable values to the
119# dataset like this (new feature introduced in ROOT 6.26):
120
121data.setGlobalObservables({mu_obs})
122
123# reset original mu_obs value
124mu_obs.setVal(mu_obs_orig_val)
125
126# Fitting a model with global observables
127# ---------------------------------------
128
129# Create snapshot of original parameters to reset parameters after fitting
130modelParameters = model.getParameters(data.get())
131origParameters = modelParameters.snapshot()
132
133# When you fit a model that includes global observables, you need to
134# specify them in the call to RooAbsPdf::fitTo with the
135# RooFit::GlobalObservables command argument. By default, the global
136# observable values attached to the dataset will be prioritized over the
137# values in the model, so the following fit correctly uses the randomized
138# global observable values from the toy dataset:
139print("1. model.fitTo(*data, GlobalObservables(mu_obs))")
140print("------------------------------------------------")
141model.fitTo(data, GlobalObservables=mu_obs, PrintLevel=-1, Save=True).Print()
142modelParameters.assign(origParameters)
143
144# In our example, the set of global observables is attached to the toy
145# dataset. In this case, you can actually drop the GlobalObservables()
146# command argument, because the global observables are automatically
147# figured out from the data set (this fit result should be identical to the
148# previous one).
149print("2. model.fitTo(*data)")
150print("---------------------")
151model.fitTo(data, PrintLevel=-1, Save=True).Print()
152modelParameters.assign(origParameters)
153
154# If you want to explicitly ignore the global observables in the dataset,
155# you can do that by specifying GlobalObservablesSource("model"). Keep in
156# mind that now it's also again your responsibility to define the set of
157# global observables.
158print('3. model.fitTo(*data, GlobalObservables(mu_obs), GlobalObservablesSource("model"))')
159print("------------------------------------------------")
160model.fitTo(data, GlobalObservables=mu_obs, GlobalObservablesSource="model", PrintLevel=-1, Save=True).Print()
161modelParameters.assign(origParameters)
void Print(GNN_Data &d, std::string txt="")