This tutorial explains the concept of global observables in RooFit, and showcases how their values can be stored either in the model or in the dataset.
Introduction
Note: in this tutorial, we are multiplying the likelihood with an additional likelihood to constrain the parameters with auxiliary measurements. This is different from the rf604_constraints
tutorial, where the likelihood is multiplied with a Bayesian prior to constrain the parameters.
With RooFit, you usually optimize some model parameters p
to maximize the likelihood L
given the per-event or per-bin observations x
:
\[ L( x | p ) \]
Often, the parameters are constrained with some prior likelihood C
, which doesn't depend on the observables x
:
\[ L'( x | p ) = L( x | p ) * C( p ) \]
Usually, these constraint terms depend on some auxiliary measurements of other observables g
. The constraint term is then the likelihood of the so-called global observables:
\[ L'( x | p ) = L( x | p ) * C( g | p ) \]
For example, think of a model where the true luminosity lumi
is a nuisance parameter that is constrained by an auxiliary measurement lumi_obs
with uncertainty lumi_obs_sigma
:
\[ L'(data | mu, lumi) = L(data | mu, lumi) * \text{Gauss}(lumi_obs | lumi, lumi_obs_sigma) \]
As a Gaussian is symmetric under exchange of the observable and the mean parameter, you can also sometimes find this equivalent but less conventional formulation for Gaussian constraints:
\[ L'(data | mu, lumi) = L(data | mu, lumi) * \text{Gauss}(lumi | lumi_obs, lumi_obs_sigma) \]
If you wanted to constrain a parameter that represents event counts, you would use a Poissonian constraint, e.g.:
\[ L'(data | mu, count) = L(data | mu, count) * \text{Poisson}(count_obs | count) \]
Unlike a Gaussian, a Poissonian is not symmetric under exchange of the observable and the parameter, so here you need to be more careful to follow the global observable prescription correctly.
{
RooProdPdf model(
"model",
"model", {gauss, constraint});
std::unique_ptr<RooDataSet>
data{model.generate({
x}, 50)};
using FitRes = std::unique_ptr<RooFitResult>;
std::cout << "1. model.fitTo(*data, GlobalObservables(mu_obs))\n";
std::cout << "------------------------------------------------\n";
std::cout << "2. model.fitTo(*data)\n";
std::cout << "---------------------\n";
std::cout << "3. model.fitTo(*data, GlobalObservables(mu_obs), GlobalObservablesSource(\"model\"))\n";
std::cout << "------------------------------------------------\n";
}
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 data
RooArgSet is a container object that can hold multiple RooAbsArg objects.
static RooMsgService & instance()
Return reference to singleton instance.
Efficient implementation of a product of PDFs of the form.
Variable that can be changed from the outside.
RooCmdArg Save(bool flag=true)
RooCmdArg GlobalObservables(Args_t &&... argsOrArgSet)
RooCmdArg GlobalObservablesSource(const char *sourceName)
RooCmdArg PrintLevel(Int_t code)
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
1. model.fitTo(*data, GlobalObservables(mu_obs))
------------------------------------------------
RooFitResult: minimized FCN value: 68.2482, estimated distance to minimum: 9.80327e-07
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
mu 5.2717e-02 +/- 8.11e-02
sigma 9.7190e-01 +/- 9.73e-02
2. model.fitTo(*data)
---------------------
RooFitResult: minimized FCN value: 68.2482, estimated distance to minimum: 9.80327e-07
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
mu 5.2717e-02 +/- 8.11e-02
sigma 9.7190e-01 +/- 9.73e-02
3. model.fitTo(*data, GlobalObservables(mu_obs), GlobalObservablesSource("model"))
------------------------------------------------
RooFitResult: minimized FCN value: 83.7181, estimated distance to minimum: 6.67911e-07
covariance matrix quality: Full, accurate covariance matrix
Status : MINIMIZE=0 HESSE=0
Floating Parameter FinalValue +/- Error
-------------------- --------------------------
mu 7.4744e-01 +/- 9.68e-02
sigma 1.2451e+00 +/- 1.38e-01
- Date
- January 2022
- Author
- Jonas Rembser
Definition in file rf613_global_observables.C.