Logo ROOT  
Reference Guide
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 conjuction 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 "RooAbsReal.h"
25#include "RooAbsCategory.h"
26#include "RooRealVar.h"
27#include <math.h>
28#include "TMath.h"
29
30using namespace std ;
31
33
34////////////////////////////////////////////////////////////////////////////////
35/// Populate x with observables
36
37RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, Bool_t paramRelative) :
38 RooAbsReal(name,title),
39 _x("x","x",this),
40 _p("p","p",this),
41 _dh(dh),
42 _relParam(paramRelative)
43{
44 _x.add(*_dh.get()) ;
45
46 // Now populate p with parameters
47 RooArgSet allVars ;
48 for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
49 _dh.get(i) ;
50
51 const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
52 RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
53 var->setVal(_relParam ? 1 : _dh.weight()) ;
54 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
55 var->setConstant(kTRUE) ;
56 allVars.add(*var) ;
57 _p.add(*var) ;
58 }
59 addOwnedComponents(allVars) ;
60}
61
62////////////////////////////////////////////////////////////////////////////////
63/// Populate x with observables
64
65RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, const RooAbsArg& /*x*/, RooDataHist& dh, Bool_t paramRelative) :
66 RooAbsReal(name,title),
67 _x("x","x",this),
68 _p("p","p",this),
69 _dh(dh),
70 _relParam(paramRelative)
71{
72 _x.add(*_dh.get()) ;
73
74 // Now populate p with parameters
75 RooArgSet allVars ;
76 for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
77 _dh.get(i) ;
78 const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
79 RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
80 var->setVal(_relParam ? 1 : _dh.weight()) ;
81 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
82 var->setConstant(kTRUE) ;
83 allVars.add(*var) ;
84 _p.add(*var) ;
85 }
86 addOwnedComponents(allVars) ;
87}
88
89////////////////////////////////////////////////////////////////////////////////
90
91RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, const RooParamHistFunc& paramSource, Bool_t paramRelative) :
92 RooAbsReal(name,title),
93 _x("x","x",this),
94 _p("p","p",this),
95 _dh(dh),
96 _relParam(paramRelative)
97{
98 // Populate x with observables
99 _x.add(*_dh.get()) ;
100
101 // Now populate p with existing parameters
102 _p.add(paramSource._p) ;
103}
104
105////////////////////////////////////////////////////////////////////////////////
106
108 RooAbsReal(other,name),
109 _x("x",this,other._x),
110 _p("p",this,other._p),
111 _dh(other._dh),
112 _relParam(other._relParam)
113{
114}
115
116////////////////////////////////////////////////////////////////////////////////
117
119{
120 Int_t idx = ((RooDataHist&)_dh).getIndex(_x,kTRUE) ;
121 Double_t ret = ((RooAbsReal*)_p.at(idx))->getVal() ;
122 if (_relParam) {
123 Double_t nom = getNominal(idx) ;
124 ret *= nom ;
125 }
126 return ret ;
127}
128
129////////////////////////////////////////////////////////////////////////////////
130
132{
133 return ((RooAbsReal&)_p[ibin]).getVal() ;
134}
135
136////////////////////////////////////////////////////////////////////////////////
137
139{
140 ((RooRealVar&)_p[ibin]).setVal(newVal) ;
141}
142
143////////////////////////////////////////////////////////////////////////////////
144
146{
147 _dh.get(ibin) ;
148 return _dh.weight() ;
149}
150
151////////////////////////////////////////////////////////////////////////////////
152
154{
155 _dh.get(ibin) ;
156 return _dh.weightError() ;
157}
158
159////////////////////////////////////////////////////////////////////////////////
160/// Return sampling hint for making curves of (projections) of this function
161/// as the recursive division strategy of RooCurve cannot deal efficiently
162/// with the vertical lines that occur in a non-interpolated histogram
163
165{
166 // Check that observable is in dataset, if not no hint is generated
167 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
168 if (!lvarg) {
169 return 0 ;
170 }
171
172 // Retrieve position of all bin boundaries
173 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
174 Double_t* boundaries = binning->array() ;
175
176 list<Double_t>* hint = new list<Double_t> ;
177
178 // Widen range slightly
179 xlo = xlo - 0.01*(xhi-xlo) ;
180 xhi = xhi + 0.01*(xhi-xlo) ;
181
182 Double_t delta = (xhi-xlo)*1e-8 ;
183
184 // Construct array with pairs of points positioned epsilon to the left and
185 // right of the bin boundaries
186 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
187 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
188 hint->push_back(boundaries[i]-delta) ;
189 hint->push_back(boundaries[i]+delta) ;
190 }
191 }
192
193 return hint ;
194}
195
196////////////////////////////////////////////////////////////////////////////////
197/// Return sampling hint for making curves of (projections) of this function
198/// as the recursive division strategy of RooCurve cannot deal efficiently
199/// with the vertical lines that occur in a non-interpolated histogram
200
201std::list<Double_t>* RooParamHistFunc::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
202{
203 // Check that observable is in dataset, if not no hint is generated
204 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
205 if (!lvarg) {
206 return 0 ;
207 }
208
209 // Retrieve position of all bin boundaries
210 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
211 Double_t* boundaries = binning->array() ;
212
213 list<Double_t>* hint = new list<Double_t> ;
214
215 // Construct array with pairs of points positioned epsilon to the left and
216 // right of the bin boundaries
217 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
218 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
219 hint->push_back(boundaries[i]) ;
220 }
221 }
222
223 return hint ;
224}
225
226////////////////////////////////////////////////////////////////////////////////
227/// Advertise that all integrals can be handled internally.
228
230 const RooArgSet* /*normSet*/, const char* /*rangeName*/) const
231{
232 // Simplest scenario, integrate over all dependents
233 RooAbsCollection *allVarsCommon = allVars.selectCommon(_x) ;
234 Bool_t intAllObs = (allVarsCommon->getSize()==_x.getSize()) ;
235 delete allVarsCommon ;
236 if (intAllObs && matchArgs(allVars,analVars,_x)) {
237 return 1 ;
238 }
239
240 return 0 ;
241}
242
243////////////////////////////////////////////////////////////////////////////////
244/// Implement analytical integrations by doing appropriate weighting from component integrals
245/// functions to integrators of components
246
247Double_t RooParamHistFunc::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
248{
249 R__ASSERT(code==1) ;
250
251 Double_t ret(0) ;
252 Int_t i(0) ;
253 for (const auto param : _p) {
254 auto p = static_cast<const RooAbsReal*>(param);
255 Double_t bin = p->getVal() ;
256 if (_relParam) bin *= getNominal(i++) ;
257 ret += bin ;
258 }
259
260 // WVE fix this!!! Assume uniform binning for now!
261 Double_t binV(1) ;
262 for (const auto obs : _x) {
263 auto xx = static_cast<const RooRealVar*>(obs);
264 binV *= (xx->getMax()-xx->getMin())/xx->numBins() ;
265 }
266
267 return ret*binV ;
268}
#define e(i)
Definition: RSha256.hxx:103
const Bool_t kTRUE
Definition: RtypesCore.h:89
#define ClassImp(name)
Definition: Rtypes.h:361
#define R__ASSERT(e)
Definition: TError.h:96
char name[80]
Definition: TGX11.cxx:109
double sqrt(double)
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:73
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2107
RooAbsBinning is the abstract base class for RooRealVar binning definitions This class defines the in...
Definition: RooAbsBinning.h:26
virtual Double_t * array() const =0
virtual Int_t numBoundaries() const =0
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
Int_t getSize() const
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for objects that are lvalues, i.e.
Definition: RooAbsLValue.h:26
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_t value=kTRUE)
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:60
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition: RooAbsReal.h:90
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition: RooArgList.h:74
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition: RooArgSet.h:28
virtual Bool_t add(const RooAbsCollection &col, Bool_t silent=kFALSE)
Add a collection of arguments to this collection by calling add() for each element in the source coll...
Definition: RooArgSet.h:88
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
virtual Double_t weight() const
Definition: RooDataHist.h:106
virtual void weightError(Double_t &lo, Double_t &hi, ErrorType etype=Poisson) const
Return the error on current weight.
virtual Int_t numEntries() const
Return the number of bins.
virtual const RooArgSet * get() const
Definition: RooDataHist.h:79
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
A histogram function that assigns scale parameters to every bin.
void setActual(Int_t ibin, Double_t newVal)
Double_t evaluate() const
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Double_t getNominal(Int_t ibin) const
Double_t getNominalError(Int_t ibin) const
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Return sampling hint for making curves of (projections) of this function as the recursive division st...
Double_t getActual(Int_t ibin)
Double_t analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName=0) const
Implement analytical integrations by doing appropriate weighting from component integrals functions t...
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
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=0) const
Advertise that all integrals can be handled internally.
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:35
void setError(Double_t value)
Definition: RooRealVar.h:64
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:261
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47