Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooParamHistFunc.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * *
4 * This code was autogenerated by RooClassFactory *
5 *****************************************************************************/
6
7/** \class RooParamHistFunc
8 * \ingroup Roofit
9 * A histogram function that assigns scale parameters to every bin. Instead of the bare bin contents,
10 * it therefore yields:
11 * \f[
12 * \gamma_{i} * \mathrm{bin}_i
13 * \f]
14 *
15 * The \f$ \gamma_i \f$ can therefore be used to parametrise statistical uncertainties of the histogram
16 * template. In conjunction with a constraint term, this can be used to implement the Barlow-Beeston method.
17 * The constraint can be implemented using RooHistConstraint.
18 *
19 * See also the tutorial rf709_BarlowBeeston.C
20 */
21
22#include "Riostream.h"
23#include "RooParamHistFunc.h"
24#include "RooAbsCategory.h"
25#include "RooRealVar.h"
26#include "RooHelpers.h"
27#include <math.h>
28#include "TMath.h"
29
30
31using namespace std ;
32
34
35////////////////////////////////////////////////////////////////////////////////
36/// Populate x with observables
37
38RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, bool paramRelative) :
39 RooAbsReal(name,title),
40 _x("x","x",this),
41 _p("p","p",this),
42 _dh(dh),
43 _relParam(paramRelative)
44{
45 _x.add(*_dh.get()) ;
46
47 // Now populate p with parameters
48 RooArgSet allVars ;
49 for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
50 _dh.get(i) ;
51
52 const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
53 RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
54 var->setVal(_relParam ? 1 : _dh.weight()) ;
55 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
56 var->setConstant(true) ;
57 allVars.add(*var) ;
58 _p.add(*var) ;
59 }
60 addOwnedComponents(allVars) ;
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Populate x with observables
65
66RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, const RooAbsArg& /*x*/, RooDataHist& dh, bool paramRelative) :
67 RooAbsReal(name,title),
68 _x("x","x",this),
69 _p("p","p",this),
70 _dh(dh),
71 _relParam(paramRelative)
72{
73 _x.add(*_dh.get()) ;
74
75 // Now populate p with parameters
76 RooArgSet allVars ;
77 for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
78 _dh.get(i) ;
79 const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
80 RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
81 var->setVal(_relParam ? 1 : _dh.weight()) ;
82 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
83 var->setConstant(true) ;
84 allVars.add(*var) ;
85 _p.add(*var) ;
86 }
87 addOwnedComponents(allVars) ;
88}
89
90////////////////////////////////////////////////////////////////////////////////
91
92RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, const RooParamHistFunc& paramSource, bool paramRelative) :
93 RooAbsReal(name,title),
94 _x("x","x",this),
95 _p("p","p",this),
96 _dh(dh),
97 _relParam(paramRelative)
98{
99 // Populate x with observables
100 _x.add(*_dh.get()) ;
101
102 // Now populate p with existing parameters
103 _p.add(paramSource._p) ;
104}
105
106////////////////////////////////////////////////////////////////////////////////
107
109 RooAbsReal(other,name),
110 _x("x",this,other._x),
111 _p("p",this,other._p),
112 _dh(other._dh),
113 _relParam(other._relParam)
114{
115}
116
117////////////////////////////////////////////////////////////////////////////////
118
120{
121 Int_t idx = ((RooDataHist&)_dh).getIndex(_x,true) ;
122 double ret = ((RooAbsReal*)_p.at(idx))->getVal() ;
123 if (_relParam) {
124 double nom = getNominal(idx) ;
125 ret *= nom ;
126 }
127 return ret ;
128}
129
130////////////////////////////////////////////////////////////////////////////////
131
133{
134 return ((RooAbsReal&)_p[ibin]).getVal() ;
135}
136
137////////////////////////////////////////////////////////////////////////////////
138
139void RooParamHistFunc::setActual(Int_t ibin, double newVal)
140{
141 ((RooRealVar&)_p[ibin]).setVal(newVal) ;
142}
143
144////////////////////////////////////////////////////////////////////////////////
145
147{
148 _dh.get(ibin) ;
149 return _dh.weight() ;
150}
151
152////////////////////////////////////////////////////////////////////////////////
153
155{
156 _dh.get(ibin) ;
157 return _dh.weightError() ;
158}
159
160////////////////////////////////////////////////////////////////////////////////
161/// Return sampling hint for making curves of (projections) of this function
162/// as the recursive division strategy of RooCurve cannot deal efficiently
163/// with the vertical lines that occur in a non-interpolated histogram
164
165list<double>* RooParamHistFunc::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
166{
167 // Check that observable is in dataset, if not no hint is generated
168 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
169 if (!lvarg) {
170 return 0 ;
171 }
172
173 // Retrieve position of all bin boundaries
174 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
175 double* boundaries = binning->array() ;
176
177 list<double>* hint = new list<double> ;
178
179 // Widen range slightly
180 xlo = xlo - 0.01*(xhi-xlo) ;
181 xhi = xhi + 0.01*(xhi-xlo) ;
182
183 double delta = (xhi-xlo)*1e-8 ;
184
185 // Construct array with pairs of points positioned epsilon to the left and
186 // right of the bin boundaries
187 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
188 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
189 hint->push_back(boundaries[i]-delta) ;
190 hint->push_back(boundaries[i]+delta) ;
191 }
192 }
193
194 return hint ;
195}
196
197////////////////////////////////////////////////////////////////////////////////
198/// Return sampling hint for making curves of (projections) of this function
199/// as the recursive division strategy of RooCurve cannot deal efficiently
200/// with the vertical lines that occur in a non-interpolated histogram
201
202std::list<double>* RooParamHistFunc::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
203{
204 // Check that observable is in dataset, if not no hint is generated
205 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
206 if (!lvarg) {
207 return 0 ;
208 }
209
210 // Retrieve position of all bin boundaries
211 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
212 double* boundaries = binning->array() ;
213
214 list<double>* hint = new list<double> ;
215
216 // Construct array with pairs of points positioned epsilon to the left and
217 // right of the bin boundaries
218 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
219 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
220 hint->push_back(boundaries[i]) ;
221 }
222 }
223
224 return hint ;
225}
226
227////////////////////////////////////////////////////////////////////////////////
228/// Advertise that all integrals can be handled internally.
229
231 const RooArgSet* /*normSet*/, const char* /*rangeName*/) const
232{
233 // Simplest scenario, integrate over all dependents
234 RooAbsCollection *allVarsCommon = allVars.selectCommon(_x) ;
235 bool intAllObs = (allVarsCommon->getSize()==_x.getSize()) ;
236 delete allVarsCommon ;
237 if (intAllObs && matchArgs(allVars,analVars,_x)) {
238 return 1 ;
239 }
240
241 return 0 ;
242}
243
244////////////////////////////////////////////////////////////////////////////////
245/// Implement analytical integrations by doing appropriate weighting from component integrals
246/// functions to integrators of components
247
248double RooParamHistFunc::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* rangeName) const
249{
250 // Supports only the scenario of integration over all dependents
251 R__ASSERT(code==1) ;
252
253 // The logic for summing over the histogram is borrowed from RooHistPdf with some differences:
254 //
255 // - a lambda function is used to inject the parameters for bin scaling into the RooDataHist::sum method
256 //
257 // - for simplicity, there is no check for the possibility of full-range integration with another overload of
258 // RooDataHist::sum
259 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
260 for (const auto obs : _x) {
261 ranges[obs] = RooHelpers::getRangeOrBinningInterval(obs, rangeName);
262 }
263
264 auto getBinScale = [&](int iBin){ return static_cast<const RooAbsReal&>(_p[iBin]).getVal(); };
265
266 RooArgSet sliceSet{};
267 return const_cast<RooDataHist&>(_dh).sum(_x, sliceSet, true, false, ranges, getBinScale);
268}
#define e(i)
Definition RSha256.hxx:103
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:117
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
Formats a string in a circular formatting buffer.
Definition TString.cxx:2467
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:74
bool addOwnedComponents(const RooAbsCollection &comps)
Take ownership of the contents of 'comps'.
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
virtual Int_t numBoundaries() const =0
virtual double * array() const =0
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
Return the number of elements in the collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
bool selectCommon(const RooAbsCollection &refColl, RooAbsCollection &outColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
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
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
void setConstant(bool value=true)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
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:55
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...
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
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:76
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 > * 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.
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...
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
void setVal(double value) override
Set value of variable to 'value'.
void setError(double value)
Definition RooRealVar.h:64
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
std::pair< double, double > getRangeOrBinningInterval(RooAbsArg const *arg, const char *rangeName)
Get the lower and upper bound of parameter range if arg can be casted to RooAbsRealLValue.
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345