Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooParamHistFunc.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2023, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11/** \class RooParamHistFunc
12 * \ingroup Roofit
13 * A histogram function that assigns scale parameters to every bin. Instead of the bare bin contents,
14 * it therefore yields:
15 * \f[
16 * \gamma_{i} * \mathrm{bin}_i
17 * \f]
18 *
19 * The \f$ \gamma_i \f$ can therefore be used to parametrise statistical uncertainties of the histogram
20 * template. In conjunction with a constraint term, this can be used to implement the Barlow-Beeston method.
21 * The constraint can be implemented using RooHistConstraint.
22 *
23 * See also the tutorial rf709_BarlowBeeston.C
24 */
25
26#include <RooParamHistFunc.h>
27#include <RooRealVar.h>
28#include <RooFitImplHelpers.h>
29
31
32RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist &dh, const RooAbsArg &x,
33 const RooParamHistFunc *paramSource, bool paramRelative)
34 : RooAbsReal(name, title), _x("x", "x", this), _p("p", "p", this), _dh(dh), _relParam(paramRelative)
35{
36 // Populate x with observables
37 _x.add(x);
38
39 if (paramSource) {
40 // Now populate p with existing parameters
41 _p.add(paramSource->_p);
42 return;
43 }
44
45 // Now populate p with parameters
46 RooArgSet allVars;
47 for (Int_t i = 0; i < _dh.numEntries(); i++) {
48 _dh.get(i);
49 const char *vname = Form("%s_gamma_bin_%i", GetName(), i);
50 RooRealVar *var = new RooRealVar(vname, vname, 0, 1000);
51 var->setVal(_relParam ? 1 : _dh.weight());
52 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
53 var->setConstant(true);
54 allVars.add(*var);
55 _p.add(*var);
56 }
57 addOwnedComponents(allVars);
58}
59
60////////////////////////////////////////////////////////////////////////////////
61
63 RooAbsReal(other,name),
64 _x("x",this,other._x),
65 _p("p",this,other._p),
66 _dh(other._dh),
67 _relParam(other._relParam)
68{
69}
70
71////////////////////////////////////////////////////////////////////////////////
72
74{
75 Int_t idx = ((RooDataHist&)_dh).getIndex(_x,true) ;
76 double ret = (static_cast<RooAbsReal*>(_p.at(idx)))->getVal() ;
77 return _relParam ? ret * getNominal(idx) : ret;
78}
79
80////////////////////////////////////////////////////////////////////////////////
81
83{
84 return (static_cast<RooAbsReal&>(_p[ibin])).getVal() ;
85}
86
87////////////////////////////////////////////////////////////////////////////////
88
89void RooParamHistFunc::setActual(Int_t ibin, double newVal)
90{
91 (static_cast<RooRealVar&>(_p[ibin])).setVal(newVal) ;
92}
93
94////////////////////////////////////////////////////////////////////////////////
95
97{
98 _dh.get(ibin) ;
99 return _dh.weight() ;
100}
101
102////////////////////////////////////////////////////////////////////////////////
103
105{
106 _dh.get(ibin) ;
107 return _dh.weightError() ;
108}
109
110////////////////////////////////////////////////////////////////////////////////
111/// Return sampling hint for making curves of (projections) of this function
112/// as the recursive division strategy of RooCurve cannot deal efficiently
113/// with the vertical lines that occur in a non-interpolated histogram
114
115std::list<double>* RooParamHistFunc::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
116{
117 // Check that observable is in dataset, if not no hint is generated
118 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
119 if (!lvarg) {
120 return nullptr ;
121 }
122
123 // Retrieve position of all bin boundaries
124 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
125 double* boundaries = binning->array() ;
126
127 std::list<double>* hint = new std::list<double> ;
128
129 // Widen range slightly
130 xlo = xlo - 0.01*(xhi-xlo) ;
131 xhi = xhi + 0.01*(xhi-xlo) ;
132
133 double delta = (xhi-xlo)*1e-8 ;
134
135 // Construct array with pairs of points positioned epsilon to the left and
136 // right of the bin boundaries
137 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
138 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
139 hint->push_back(boundaries[i]-delta) ;
140 hint->push_back(boundaries[i]+delta) ;
141 }
142 }
143
144 return hint ;
145}
146
147////////////////////////////////////////////////////////////////////////////////
148/// Return sampling hint for making curves of (projections) of this function
149/// as the recursive division strategy of RooCurve cannot deal efficiently
150/// with the vertical lines that occur in a non-interpolated histogram
151
152std::list<double>* RooParamHistFunc::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
153{
154 // Check that observable is in dataset, if not no hint is generated
155 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
156 if (!lvarg) {
157 return nullptr ;
158 }
159
160 // Retrieve position of all bin boundaries
161 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
162 double* boundaries = binning->array() ;
163
164 std::list<double>* hint = new std::list<double> ;
165
166 // Construct array with pairs of points positioned epsilon to the left and
167 // right of the bin boundaries
168 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
169 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
170 hint->push_back(boundaries[i]) ;
171 }
172 }
173
174 return hint ;
175}
176
177////////////////////////////////////////////////////////////////////////////////
178/// Advertise that all integrals can be handled internally.
179
181 const RooArgSet* /*normSet*/, const char* /*rangeName*/) const
182{
183 // Simplest scenario, integrate over all dependents
184 std::unique_ptr<RooAbsCollection> allVarsCommon{allVars.selectCommon(_x)};
185 bool intAllObs = (allVarsCommon->size()==_x.size()) ;
186 if (intAllObs && matchArgs(allVars,analVars,_x)) {
187 return 1 ;
188 }
189
190 return 0 ;
191}
192
193////////////////////////////////////////////////////////////////////////////////
194/// Implement analytical integrations by doing appropriate weighting from component integrals
195/// functions to integrators of components
196
197double RooParamHistFunc::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* rangeName) const
198{
199 // Supports only the scenario of integration over all dependents
200 R__ASSERT(code==1) ;
201
202 // The logic for summing over the histogram is borrowed from RooHistPdf with some differences:
203 //
204 // - a lambda function is used to inject the parameters for bin scaling into the RooDataHist::sum method
205 //
206 // - for simplicity, there is no check for the possibility of full-range integration with another overload of
207 // RooDataHist::sum
208 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
209 for (const auto obs : _x) {
210 ranges[obs] = RooHelpers::getRangeOrBinningInterval(obs, rangeName);
211 }
212
213 auto getBinScale = [&](int iBin){ return static_cast<const RooAbsReal&>(_p[iBin]).getVal(); };
214
215 RooArgSet sliceSet{};
216 return const_cast<RooDataHist&>(_dh).sum(_x, sliceSet, true, false, ranges, getBinScale);
217}
#define e(i)
Definition RSha256.hxx:103
#define ClassImp(name)
Definition Rtypes.h:382
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2489
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
Abstract base class for RooRealVar binning definitions.
virtual Int_t numBoundaries() const =0
virtual double * array() const =0
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract base class for objects that are lvalues, i.e.
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const =0
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
void setConstant(bool value=true)
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
double weight(std::size_t i) const
Return weight of i-th bin.
void weightError(double &lo, double &hi, ErrorType etype=Poisson) const override
Return the asymmetric errors on the current weight.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
A histogram function that assigns scale parameters to every bin.
void setActual(Int_t ibin, double newVal)
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
double getNominalError(Int_t ibin) const
double getActual(Int_t ibin)
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
double getNominal(Int_t ibin) const
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Advertise that all integrals can be handled internally.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
void setError(double value)
Definition RooRealVar.h:60
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
std::pair< double, double > getRangeOrBinningInterval(RooAbsArg const *arg, const char *rangeName)
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345