Logo ROOT   6.16/01
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
10*/
11
12#include "Riostream.h"
13#include "RooParamHistFunc.h"
14#include "RooAbsReal.h"
15#include "RooAbsCategory.h"
16#include "RooRealVar.h"
17#include <math.h>
18#include "TMath.h"
19
20using namespace std ;
21
23
24////////////////////////////////////////////////////////////////////////////////
25/// Populate x with observables
26
27RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, Bool_t paramRelative) :
28 RooAbsReal(name,title),
29 _x("x","x",this),
30 _p("p","p",this),
31 _dh(dh),
32 _relParam(paramRelative)
33{
34 _x.add(*_dh.get()) ;
35
36 // Now populate p with parameters
37 RooArgSet allVars ;
38 for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
39 _dh.get(i) ;
40 const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
41 RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
42 var->setVal(_relParam ? 1 : _dh.weight()) ;
43 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
44 var->setConstant(kTRUE) ;
45 allVars.add(*var) ;
46 _p.add(*var) ;
47 }
48 addOwnedComponents(allVars) ;
49}
50
51////////////////////////////////////////////////////////////////////////////////
52/// Populate x with observables
53
54RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, const RooAbsArg& /*x*/, RooDataHist& dh, Bool_t paramRelative) :
55 RooAbsReal(name,title),
56 _x("x","x",this),
57 _p("p","p",this),
58 _dh(dh),
59 _relParam(paramRelative)
60{
61 _x.add(*_dh.get()) ;
62
63 // Now populate p with parameters
64 RooArgSet allVars ;
65 for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
66 _dh.get(i) ;
67 const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
68 RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
69 var->setVal(_relParam ? 1 : _dh.weight()) ;
70 var->setError(_relParam ? 1 / sqrt(_dh.weight()) : sqrt(_dh.weight()));
71 var->setConstant(kTRUE) ;
72 allVars.add(*var) ;
73 _p.add(*var) ;
74 }
75 addOwnedComponents(allVars) ;
76}
77
78////////////////////////////////////////////////////////////////////////////////
79
80RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, const RooParamHistFunc& paramSource, Bool_t paramRelative) :
81 RooAbsReal(name,title),
82 _x("x","x",this),
83 _p("p","p",this),
84 _dh(dh),
85 _relParam(paramRelative)
86{
87 // Populate x with observables
88 _x.add(*_dh.get()) ;
89
90 // Now populate p with existing parameters
91 _p.add(paramSource._p) ;
92}
93
94////////////////////////////////////////////////////////////////////////////////
95
97 RooAbsReal(other,name),
98 _x("x",this,other._x),
99 _p("p",this,other._p),
100 _dh(other._dh),
101 _relParam(other._relParam)
102{
103}
104
105////////////////////////////////////////////////////////////////////////////////
106
108{
109 Int_t idx = ((RooDataHist&)_dh).getIndex(_x,kTRUE) ;
110 Double_t ret = ((RooAbsReal*)_p.at(idx))->getVal() ;
111 if (_relParam) {
112 Double_t nom = getNominal(idx) ;
113 ret *= nom ;
114 }
115 return ret ;
116}
117
118////////////////////////////////////////////////////////////////////////////////
119
121{
122 return ((RooAbsReal&)_p[ibin]).getVal() ;
123}
124
125////////////////////////////////////////////////////////////////////////////////
126
128{
129 ((RooRealVar&)_p[ibin]).setVal(newVal) ;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
135{
136 _dh.get(ibin) ;
137 return _dh.weight() ;
138}
139
140////////////////////////////////////////////////////////////////////////////////
141
143{
144 _dh.get(ibin) ;
145 return _dh.weightError() ;
146}
147
148////////////////////////////////////////////////////////////////////////////////
149/// Return sampling hint for making curves of (projections) of this function
150/// as the recursive division strategy of RooCurve cannot deal efficiently
151/// with the vertical lines that occur in a non-interpolated histogram
152
154{
155 // Check that observable is in dataset, if not no hint is generated
156 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
157 if (!lvarg) {
158 return 0 ;
159 }
160
161 // Retrieve position of all bin boundaries
162 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
163 Double_t* boundaries = binning->array() ;
164
165 list<Double_t>* hint = new list<Double_t> ;
166
167 // Widen range slightly
168 xlo = xlo - 0.01*(xhi-xlo) ;
169 xhi = xhi + 0.01*(xhi-xlo) ;
170
171 Double_t delta = (xhi-xlo)*1e-8 ;
172
173 // Construct array with pairs of points positioned epsilon to the left and
174 // right of the bin boundaries
175 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
176 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
177 hint->push_back(boundaries[i]-delta) ;
178 hint->push_back(boundaries[i]+delta) ;
179 }
180 }
181
182 return hint ;
183}
184
185////////////////////////////////////////////////////////////////////////////////
186/// Return sampling hint for making curves of (projections) of this function
187/// as the recursive division strategy of RooCurve cannot deal efficiently
188/// with the vertical lines that occur in a non-interpolated histogram
189
190std::list<Double_t>* RooParamHistFunc::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
191{
192 // Check that observable is in dataset, if not no hint is generated
193 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
194 if (!lvarg) {
195 return 0 ;
196 }
197
198 // Retrieve position of all bin boundaries
199 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
200 Double_t* boundaries = binning->array() ;
201
202 list<Double_t>* hint = new list<Double_t> ;
203
204 // Construct array with pairs of points positioned epsilon to the left and
205 // right of the bin boundaries
206 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
207 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
208 hint->push_back(boundaries[i]) ;
209 }
210 }
211
212 return hint ;
213}
214
215////////////////////////////////////////////////////////////////////////////////
216/// Advertise that all integrals can be handled internally.
217
219 const RooArgSet* /*normSet*/, const char* /*rangeName*/) const
220{
221 // Simplest scenario, integrate over all dependents
222 RooAbsCollection *allVarsCommon = allVars.selectCommon(_x) ;
223 Bool_t intAllObs = (allVarsCommon->getSize()==_x.getSize()) ;
224 delete allVarsCommon ;
225 if (intAllObs && matchArgs(allVars,analVars,_x)) {
226 return 1 ;
227 }
228
229 return 0 ;
230}
231
232////////////////////////////////////////////////////////////////////////////////
233/// Implement analytical integrations by doing appropriate weighting from component integrals
234/// functions to integrators of components
235
236Double_t RooParamHistFunc::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
237{
238 R__ASSERT(code==1) ;
239
240 RooFIter iter = _p.fwdIterator() ;
241 RooAbsReal* p ;
242 Double_t ret(0) ;
243 Int_t i(0) ;
244 while((p=(RooAbsReal*)iter.next())) {
245 Double_t bin = p->getVal() ;
246 if (_relParam) bin *= getNominal(i++) ;
247 ret += bin ;
248 }
249
250 // WVE fix this!!! Assume uniform binning for now!
251 RooFIter xiter = _x.fwdIterator() ;
252 RooAbsArg* obs ;
253 Double_t binV(1) ;
254 while((obs=xiter.next())) {
255 RooRealVar* xx = (RooRealVar*) obs ;
256 binV *= (xx->getMax()-xx->getMin())/xx->numBins() ;
257 }
258
259 return ret*binV ;
260}
#define e(i)
Definition: RSha256.hxx:103
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
double Double_t
Definition: RtypesCore.h:55
const Bool_t kTRUE
Definition: RtypesCore.h:87
#define ClassImp(name)
Definition: Rtypes.h:363
#define R__ASSERT(e)
Definition: TError.h:96
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:66
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of 'comps'.
Definition: RooAbsArg.cxx:2274
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
RooFIter fwdIterator() 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...
virtual Double_t getMax(const char *name=0) const
virtual Int_t numBins(const char *rangeName=0) const
void setConstant(Bool_t value=kTRUE)
virtual Double_t getMin(const char *name=0) const
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t getVal(const RooArgSet *set=0) const
Evaluate object. Returns either cached value or triggers a recalculation.
Definition: RooAbsReal.h:64
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
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
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
virtual Double_t weight() const
Definition: RooDataHist.h:96
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:77
RooAbsArg * next()
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
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 fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
void setError(Double_t value)
Definition: RooRealVar.h:55
virtual void setVal(Double_t value)
Set value of variable to 'value'.
Definition: RooRealVar.cxx:204
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
STL namespace.