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
30
31RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist &dh, const RooAbsArg &x,
33 : RooAbsReal(name, title), _x("x", "x", this), _p("p", "p", this), _dh(dh), _relParam(paramRelative)
34{
35 // Populate x with observables
36 _x.add(x);
37
38 if (paramSource) {
39 // Now populate p with existing parameters
40 _p.add(paramSource->_p);
41 return;
42 }
43
44 // Now populate p with parameters
45 RooArgSet allVars;
46 for (Int_t i = 0; i < _dh.numEntries(); i++) {
47 const double wi = _dh.weight(i);
48 const char *vname = Form("%s_gamma_bin_%i", GetName(), i);
49 RooRealVar *var = new RooRealVar(vname, vname, 0, 1000);
50 var->setVal(_relParam ? 1 : wi);
51 var->setError(_relParam ? 1 / sqrt(wi) : sqrt(wi));
52 var->setConstant(true);
53 allVars.add(*var);
54 _p.add(*var);
55 }
56 addOwnedComponents(allVars);
57}
58
59////////////////////////////////////////////////////////////////////////////////
60
63 _x("x",this,other._x),
64 _p("p",this,other._p),
65 _dh(other._dh),
66 _relParam(other._relParam)
67{
68}
69
70////////////////////////////////////////////////////////////////////////////////
71
73{
74 Int_t idx = ((RooDataHist&)_dh).getIndex(_x,true) ;
75 double ret = (static_cast<RooAbsReal*>(_p.at(idx)))->getVal() ;
76 return _relParam ? ret * getNominal(idx) : ret;
77}
78
79////////////////////////////////////////////////////////////////////////////////
80
82{
83 return (static_cast<RooAbsReal&>(_p[ibin])).getVal() ;
84}
85
86////////////////////////////////////////////////////////////////////////////////
87
89{
90 (static_cast<RooRealVar&>(_p[ibin])).setVal(newVal) ;
91}
92
93////////////////////////////////////////////////////////////////////////////////
94
96{
97 return _dh.weight(ibin) ;
98}
99
100////////////////////////////////////////////////////////////////////////////////
101
103{
104 _dh.get(ibin) ;
105 return _dh.weightError() ;
106}
107
108////////////////////////////////////////////////////////////////////////////////
109/// Return sampling hint for making curves of (projections) of this function
110/// as the recursive division strategy of RooCurve cannot deal efficiently
111/// with the vertical lines that occur in a non-interpolated histogram
112
113std::list<double>* RooParamHistFunc::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
114{
115 // Check that observable is in dataset, if not no hint is generated
116 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
117 if (!lvarg) {
118 return nullptr ;
119 }
120
121 // Retrieve position of all bin boundaries
122 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
123 double* boundaries = binning->array() ;
124
125 std::list<double>* hint = new std::list<double> ;
126
127 // Widen range slightly
128 xlo = xlo - 0.01*(xhi-xlo) ;
129 xhi = xhi + 0.01*(xhi-xlo) ;
130
131 double delta = (xhi-xlo)*1e-8 ;
132
133 // Construct array with pairs of points positioned epsilon to the left and
134 // right of the bin boundaries
135 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
136 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
137 hint->push_back(boundaries[i]-delta) ;
138 hint->push_back(boundaries[i]+delta) ;
139 }
140 }
141
142 return hint ;
143}
144
145////////////////////////////////////////////////////////////////////////////////
146/// Return sampling hint for making curves of (projections) of this function
147/// as the recursive division strategy of RooCurve cannot deal efficiently
148/// with the vertical lines that occur in a non-interpolated histogram
149
150std::list<double>* RooParamHistFunc::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
151{
152 // Check that observable is in dataset, if not no hint is generated
153 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
154 if (!lvarg) {
155 return nullptr ;
156 }
157
158 // Retrieve position of all bin boundaries
159 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
160 double* boundaries = binning->array() ;
161
162 std::list<double>* hint = new std::list<double> ;
163
164 // Construct array with pairs of points positioned epsilon to the left and
165 // right of the bin boundaries
166 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
167 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
168 hint->push_back(boundaries[i]) ;
169 }
170 }
171
172 return hint ;
173}
174
175////////////////////////////////////////////////////////////////////////////////
176/// Advertise that all integrals can be handled internally.
177
179 const RooArgSet* /*normSet*/, const char* /*rangeName*/) const
180{
181 // Simplest scenario, integrate over all dependents
182 std::unique_ptr<RooAbsCollection> allVarsCommon{allVars.selectCommon(_x)};
183 bool intAllObs = (allVarsCommon->size()==_x.size()) ;
184 if (intAllObs && matchArgs(allVars,analVars,_x)) {
185 return 1 ;
186 }
187
188 return 0 ;
189}
190
191////////////////////////////////////////////////////////////////////////////////
192/// Implement analytical integrations by doing appropriate weighting from component integrals
193/// functions to integrators of components
194
195double RooParamHistFunc::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* rangeName) const
196{
197 // Supports only the scenario of integration over all dependents
198 R__ASSERT(code==1) ;
199
200 // The logic for summing over the histogram is borrowed from RooHistPdf with some differences:
201 //
202 // - a lambda function is used to inject the parameters for bin scaling into the RooDataHist::sum method
203 //
204 // - for simplicity, there is no check for the possibility of full-range integration with another overload of
205 // RooDataHist::sum
206 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
207 for (const auto obs : _x) {
209 }
210
211 auto getBinScale = [&](int iBin){ return static_cast<const RooAbsReal&>(_p[iBin]).getVal(); };
212
214 return const_cast<RooDataHist&>(_dh).sum(_x, sliceSet, true, false, ranges, getBinScale);
215}
#define e(i)
Definition RSha256.hxx:103
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#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:148
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2495
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
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
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.
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:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
bool matchArgs(const RooArgSet &allDeps, RooArgSet &analDeps, const RooArgProxy &a, const Proxies &... proxies) const
Definition RooAbsReal.h:425
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:61
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
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:2338