Logo ROOT   6.10/00
Reference Guide
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
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 
20 using namespace std ;
21 
23 
24 ////////////////////////////////////////////////////////////////////////////////
25 /// Populate x with observables
26 
27 RooParamHistFunc::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->setConstant(kTRUE) ;
44  allVars.add(*var) ;
45  _p.add(*var) ;
46  }
47  addOwnedComponents(allVars) ;
48 }
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 /// Populate x with observables
52 
53 RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, const RooAbsArg& /*x*/, RooDataHist& dh, Bool_t paramRelative) :
54  RooAbsReal(name,title),
55  _x("x","x",this),
56  _p("p","p",this),
57  _dh(dh),
58  _relParam(paramRelative)
59 {
60  _x.add(*_dh.get()) ;
61 
62  // Now populate p with parameters
63  RooArgSet allVars ;
64  for (Int_t i=0 ; i<_dh.numEntries() ; i++) {
65  _dh.get(i) ;
66  const char* vname = Form("%s_gamma_bin_%i",GetName(),i) ;
67  RooRealVar* var = new RooRealVar(vname,vname,0,1000) ;
68  var->setVal(_relParam ? 1 : _dh.weight()) ;
69  var->setConstant(kTRUE) ;
70  allVars.add(*var) ;
71  _p.add(*var) ;
72  }
73  addOwnedComponents(allVars) ;
74 }
75 
76 ////////////////////////////////////////////////////////////////////////////////
77 
78 RooParamHistFunc::RooParamHistFunc(const char *name, const char *title, RooDataHist& dh, const RooParamHistFunc& paramSource, Bool_t paramRelative) :
79  RooAbsReal(name,title),
80  _x("x","x",this),
81  _p("p","p",this),
82  _dh(dh),
83  _relParam(paramRelative)
84 {
85  // Populate x with observables
86  _x.add(*_dh.get()) ;
87 
88  // Now populate p with existing parameters
89  _p.add(paramSource._p) ;
90 }
91 
92 ////////////////////////////////////////////////////////////////////////////////
93 
95  RooAbsReal(other,name),
96  _x("x",this,other._x),
97  _p("p",this,other._p),
98  _dh(other._dh),
99  _relParam(other._relParam)
100 {
101 }
102 
103 ////////////////////////////////////////////////////////////////////////////////
104 
106 {
107  Int_t idx = ((RooDataHist&)_dh).getIndex(_x,kTRUE) ;
108  Double_t ret = ((RooAbsReal*)_p.at(idx))->getVal() ;
109  if (_relParam) {
110  Double_t nom = getNominal(idx) ;
111  ret *= nom ;
112  }
113  return ret ;
114 }
115 
116 ////////////////////////////////////////////////////////////////////////////////
117 
119 {
120  return ((RooAbsReal&)_p[ibin]).getVal() ;
121 }
122 
123 ////////////////////////////////////////////////////////////////////////////////
124 
126 {
127  ((RooRealVar&)_p[ibin]).setVal(newVal) ;
128 }
129 
130 ////////////////////////////////////////////////////////////////////////////////
131 
133 {
134  _dh.get(ibin) ;
135  return _dh.weight() ;
136 }
137 
138 ////////////////////////////////////////////////////////////////////////////////
139 
141 {
142  _dh.get(ibin) ;
143  return _dh.weightError() ;
144 }
145 
146 ////////////////////////////////////////////////////////////////////////////////
147 /// Return sampling hint for making curves of (projections) of this function
148 /// as the recursive division strategy of RooCurve cannot deal efficiently
149 /// with the vertical lines that occur in a non-interpolated histogram
150 
152 {
153  // Check that observable is in dataset, if not no hint is generated
154  RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
155  if (!lvarg) {
156  return 0 ;
157  }
158 
159  // Retrieve position of all bin boundaries
160  const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
161  Double_t* boundaries = binning->array() ;
162 
163  list<Double_t>* hint = new list<Double_t> ;
164 
165  // Widen range slightly
166  xlo = xlo - 0.01*(xhi-xlo) ;
167  xhi = xhi + 0.01*(xhi-xlo) ;
168 
169  Double_t delta = (xhi-xlo)*1e-8 ;
170 
171  // Construct array with pairs of points positioned epsilon to the left and
172  // right of the bin boundaries
173  for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
174  if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
175  hint->push_back(boundaries[i]-delta) ;
176  hint->push_back(boundaries[i]+delta) ;
177  }
178  }
179 
180  return hint ;
181 }
182 
183 ////////////////////////////////////////////////////////////////////////////////
184 /// Return sampling hint for making curves of (projections) of this function
185 /// as the recursive division strategy of RooCurve cannot deal efficiently
186 /// with the vertical lines that occur in a non-interpolated histogram
187 
188 std::list<Double_t>* RooParamHistFunc::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
189 {
190  // Check that observable is in dataset, if not no hint is generated
191  RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dh.get()->find(obs.GetName())) ;
192  if (!lvarg) {
193  return 0 ;
194  }
195 
196  // Retrieve position of all bin boundaries
197  const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
198  Double_t* boundaries = binning->array() ;
199 
200  list<Double_t>* hint = new list<Double_t> ;
201 
202  // Construct array with pairs of points positioned epsilon to the left and
203  // right of the bin boundaries
204  for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
205  if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
206  hint->push_back(boundaries[i]) ;
207  }
208  }
209 
210  return hint ;
211 }
212 
213 ////////////////////////////////////////////////////////////////////////////////
214 /// Advertise that all integrals can be handled internally.
215 
217  const RooArgSet* /*normSet*/, const char* /*rangeName*/) const
218 {
219  // Simplest scenario, integrate over all dependents
220  RooAbsCollection *allVarsCommon = allVars.selectCommon(_x) ;
221  Bool_t intAllObs = (allVarsCommon->getSize()==_x.getSize()) ;
222  delete allVarsCommon ;
223  if (intAllObs && matchArgs(allVars,analVars,_x)) {
224  return 1 ;
225  }
226 
227  return 0 ;
228 }
229 
230 ////////////////////////////////////////////////////////////////////////////////
231 /// Implement analytical integrations by doing appropriate weighting from component integrals
232 /// functions to integrators of components
233 
234 Double_t RooParamHistFunc::analyticalIntegralWN(Int_t code, const RooArgSet* /*normSet2*/,const char* /*rangeName*/) const
235 {
236  R__ASSERT(code==1) ;
237 
238  RooFIter iter = _p.fwdIterator() ;
239  RooAbsReal* p ;
240  Double_t ret(0) ;
241  Int_t i(0) ;
242  while((p=(RooAbsReal*)iter.next())) {
243  Double_t bin = p->getVal() ;
244  if (_relParam) bin *= getNominal(i++) ;
245  ret += bin ;
246  }
247 
248  // WVE fix this!!! Assume uniform binning for now!
249  RooFIter xiter = _x.fwdIterator() ;
250  RooAbsArg* obs ;
251  Double_t binV(1) ;
252  while((obs=xiter.next())) {
253  RooRealVar* xx = (RooRealVar*) obs ;
254  binV *= (xx->getMax()-xx->getMin())/xx->numBins() ;
255  }
256 
257  return ret*binV ;
258 }
std::string GetName(const std::string &scope_name)
Definition: Cppyy.cxx:145
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...
virtual Int_t numBins(const char *rangeName=0) const
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:86
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
RooFIter fwdIterator() const
Double_t getActual(Int_t ibin)
#define R__ASSERT(e)
Definition: TError.h:96
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 Double_t getMin(const char *name=0) const
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t evaluate() const
Bool_t addOwnedComponents(const RooArgSet &comps)
Take ownership of the contents of &#39;comps&#39;.
Definition: RooAbsArg.cxx:2282
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=0) const
Advertise that all integrals can be handled internally.
Float_t delta
virtual Int_t numBoundaries() const =0
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
RooDataSet is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:40
Double_t x[n]
Definition: legend1.C:17
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const =0
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
virtual void weightError(Double_t &lo, Double_t &hi, ErrorType etype=Poisson) const
Return the error on current weight.
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Reimplementation of standard RooArgList::add()
RooRealVar represents a fundamental (non-derived) real valued object.
Definition: RooRealVar.h:36
virtual void setVal(Double_t value)
Set value of variable to &#39;value&#39;.
Definition: RooRealVar.cxx:205
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Double_t weight() const
Definition: RooDataHist.h:96
void setConstant(Bool_t value=kTRUE)
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:47
RooAbsArg * next()
RooAbsBinning is the abstract base class for RooRealVar binning definitions This class defines the in...
Definition: RooAbsBinning.h:26
RooAbsArg * at(Int_t idx) const
Definition: RooArgList.h:84
#define ClassImp(name)
Definition: Rtypes.h:336
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
Definition: TRolke.cxx:630
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...
virtual Int_t numEntries() const
Return the number of bins.
Double_t getNominalError(Int_t ibin) const
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects...
virtual Double_t getMax(const char *name=0) const
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
void setActual(Int_t ibin, Double_t newVal)
Double_t getNominal(Int_t ibin) const
virtual Double_t * array() const =0
Abstract base class for objects that are lvalues, i.e.
Definition: RooAbsLValue.h:26
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Int_t getSize() const
RooAbsArg is the common abstract base class for objects that represent a value (of arbitrary type) an...
Definition: RooAbsArg.h:66
const Bool_t kTRUE
Definition: RtypesCore.h:91
char name[80]
Definition: TGX11.cxx:109