ROOT  6.06/09
Reference Guide
RooHistPdf.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitCore *
4  * @(#)root/roofit:$Id$
5  * Authors: *
6  * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7  * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8  * *
9  * Copyright (c) 2000-2005, Regents of the University of California *
10  * and Stanford University. All rights reserved. *
11  * *
12  * Redistribution and use in source and binary forms, *
13  * with or without modification, are permitted according to the terms *
14  * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15  *****************************************************************************/
16 
17 //////////////////////////////////////////////////////////////////////////////
18 //
19 // BEGIN_HTML
20 // RooHistPdf implements a probablity density function sampled from a
21 // multidimensional histogram. The histogram distribution is explicitly
22 // normalized by RooHistPdf and can have an arbitrary number of real or
23 // discrete dimensions.
24 // END_HTML
25 //
26 
27 #include "RooFit.h"
28 #include "Riostream.h"
29 
30 #include "RooHistPdf.h"
31 #include "RooDataHist.h"
32 #include "RooMsgService.h"
33 #include "RooRealVar.h"
34 #include "RooCategory.h"
35 #include "RooWorkspace.h"
36 #include "RooGlobalFunc.h"
37 
38 #include "TError.h"
39 
40 using namespace std;
41 
43 ;
44 
45 
46 
47 ////////////////////////////////////////////////////////////////////////////////
48 /// Default constructor
49 /// coverity[UNINIT_CTOR]
50 
51 RooHistPdf::RooHistPdf() : _dataHist(0), _totVolume(0), _unitNorm(kFALSE)
52 {
55 }
56 
57 
58 ////////////////////////////////////////////////////////////////////////////////
59 /// Constructor from a RooDataHist. RooDataHist dimensions
60 /// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
61 /// RooHistPdf neither owns or clone 'dhist' and the user must ensure the input histogram exists
62 /// for the entire life span of this PDF.
63 
64 RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgSet& vars,
65  const RooDataHist& dhist, Int_t intOrder) :
66  RooAbsPdf(name,title),
67  _pdfObsList("pdfObs","List of p.d.f. observables",this),
68  _dataHist((RooDataHist*)&dhist),
69  _codeReg(10),
70  _intOrder(intOrder),
71  _cdfBoundaries(kFALSE),
72  _totVolume(0),
73  _unitNorm(kFALSE)
74 {
75  _histObsList.addClone(vars) ;
76  _pdfObsList.add(vars) ;
77 
78  // Verify that vars and dhist.get() have identical contents
79  const RooArgSet* dvars = dhist.get() ;
80  if (vars.getSize()!=dvars->getSize()) {
81  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
82  << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
83  assert(0) ;
84  }
85  TIterator* iter = vars.createIterator() ;
86  RooAbsArg* arg ;
87  while((arg=(RooAbsArg*)iter->Next())) {
88  if (!dvars->find(arg->GetName())) {
89  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
90  << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
91  assert(0) ;
92  }
93  }
94  delete iter ;
95 
98 
99 
100  // Adjust ranges of _histObsList to those of _dataHist
101  RooFIter oiter = _histObsList.fwdIterator() ;
102  RooAbsArg* hobs ;
103  while ((hobs = oiter.next())) {
104  // Guaranteed to succeed, since checked above in ctor
105  RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
106  RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
107  if (dhreal){
108  ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
109  }
110  }
111 
112 }
113 
114 
115 
116 
117 ////////////////////////////////////////////////////////////////////////////////
118 /// Constructor from a RooDataHist. The first list of observables are the p.d.f.
119 /// observables, which may any RooAbsReal (function or variable). The second list
120 /// are the corresponding observables in the RooDataHist which must be of type
121 /// RooRealVar or RooCategory This constructor thus allows to apply a coordinate transformation
122 /// on the histogram data to be applied.
123 
124 RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs,
125  const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder) :
126  RooAbsPdf(name,title),
127  _pdfObsList("pdfObs","List of p.d.f. observables",this),
128  _dataHist((RooDataHist*)&dhist),
129  _codeReg(10),
130  _intOrder(intOrder),
131  _cdfBoundaries(kFALSE),
132  _totVolume(0),
133  _unitNorm(kFALSE)
134 {
135  _histObsList.addClone(histObs) ;
136  _pdfObsList.add(pdfObs) ;
137 
138  // Verify that vars and dhist.get() have identical contents
139  const RooArgSet* dvars = dhist.get() ;
140  if (histObs.getSize()!=dvars->getSize()) {
141  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
142  << ") ERROR histogram variable list and RooDataHist must contain the same variables." << endl ;
143  throw(string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
144  }
145  TIterator* iter = histObs.createIterator() ;
146  RooAbsArg* arg ;
147  while((arg=(RooAbsArg*)iter->Next())) {
148  if (!dvars->find(arg->GetName())) {
149  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
150  << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
151  throw(string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
152  }
153  if (!arg->isFundamental()) {
154  coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
155  << ") ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory." << endl ;
156  throw(string("RooHistPdf::ctor() ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory.")) ;
157  }
158  }
159  delete iter ;
160 
163 
164  // Adjust ranges of _histObsList to those of _dataHist
165  RooFIter oiter = _histObsList.fwdIterator() ;
166  RooAbsArg* hobs ;
167  while ((hobs = oiter.next())) {
168  // Guaranteed to succeed, since checked above in ctor
169  RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
170  RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
171  if (dhreal){
172  ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
173  }
174  }
175 }
176 
177 
178 
179 ////////////////////////////////////////////////////////////////////////////////
180 /// Copy constructor
181 
182 RooHistPdf::RooHistPdf(const RooHistPdf& other, const char* name) :
183  RooAbsPdf(other,name),
184  _pdfObsList("pdfObs",this,other._pdfObsList),
185  _dataHist(other._dataHist),
186  _codeReg(other._codeReg),
187  _intOrder(other._intOrder),
188  _cdfBoundaries(other._cdfBoundaries),
189  _totVolume(other._totVolume),
190  _unitNorm(other._unitNorm)
191 {
193 
196 }
197 
198 
199 
200 
201 ////////////////////////////////////////////////////////////////////////////////
202 /// Destructor
203 
205 {
206  delete _histObsIter ;
207  delete _pdfObsIter ;
208 }
209 
210 
211 
212 
213 
214 ////////////////////////////////////////////////////////////////////////////////
215 /// Return the current value: The value of the bin enclosing the current coordinates
216 /// of the observables, normalized by the histograms contents. Interpolation
217 /// is applied if the RooHistPdf is configured to do that
218 
220 {
221  // Transfer values from
222  if (_pdfObsList.getSize()>0) {
223  _histObsIter->Reset() ;
224  _pdfObsIter->Reset() ;
225  RooAbsArg* harg, *parg ;
226  while((harg=(RooAbsArg*)_histObsIter->Next())) {
227  parg = (RooAbsArg*)_pdfObsIter->Next() ;
228  if (harg != parg) {
229  parg->syncCache() ;
230  harg->copyCache(parg,kTRUE) ;
231  if (!harg->inRange(0)) {
232  return 0 ;
233  }
234  }
235  }
236  }
237 
239  //cout << "RooHistPdf::evaluate(" << GetName() << ") ret = " << ret << endl ;
240  if (ret<0) {
241  ret=0 ;
242  }
243  return ret ;
244 }
245 
246 
247 ////////////////////////////////////////////////////////////////////////////////
248 /// Return the total volume spanned by the observables of the RooHistPdf
249 
251 {
252  // Return previously calculated value, if any
253  if (_totVolume>0) {
254  return _totVolume ;
255  }
256  _totVolume = 1. ;
258  RooAbsArg* arg ;
259  while((arg=(RooAbsArg*)iter->Next())) {
260  RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
261  if (real) {
262  _totVolume *= (real->getMax()-real->getMin()) ;
263  } else {
264  RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
265  if (cat) {
266  _totVolume *= cat->numTypes() ;
267  }
268  }
269  }
270  delete iter ;
271  return _totVolume ;
272 }
273 
274 namespace {
275  bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
276  {
277  const RooAbsRealLValue *_x = dynamic_cast<const RooAbsRealLValue*>(&x);
278  const RooAbsRealLValue *_y = dynamic_cast<const RooAbsRealLValue*>(&y);
279  if (!_x || !_y) return false;
280  if (!range || !strlen(range) || !_x->hasRange(range) ||
281  _x->getBinningPtr(range)->isParameterized()) {
282  // parameterized ranges may be full range now, but that might change,
283  // so return false
284  if (range && strlen(range) && _x->getBinningPtr(range)->isParameterized())
285  return false;
286  return (_x->getMin() == _y->getMin() && _x->getMax() == _y->getMax());
287  }
288  return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
289  }
290 }
291 
292 
293 ////////////////////////////////////////////////////////////////////////////////
294 /// Determine integration scenario. If no interpolation is used,
295 /// RooHistPdf can perform all integrals over its dependents
296 /// analytically via partial or complete summation of the input
297 /// histogram. If interpolation is used on the integral over
298 /// all histogram observables is supported
299 
300 Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
301 {
302  // First make list of pdf observables to histogram observables
303  // and select only those for which the integral is over the full range
304 
307  Int_t code = 0, frcode = 0, n = 0;
308  for (RooAbsArg *pa = 0, *ha = 0; (pa = it.next()) && (ha = jt.next()); ++n) {
309  if (allVars.find(*pa)) {
310  code |= 2 << n;
311  analVars.add(*pa);
312  if (fullRange(*pa, *ha, rangeName)) {
313  frcode |= 2 << n;
314  }
315  }
316  }
317 
318  if (code == frcode) {
319  // integrate over full range of all observables - use bit 0 to indicate
320  // full range integration over all observables
321  code |= 1;
322  }
323  // Disable partial analytical integrals if interpolation is used, and we
324  // integrate over sub-ranges, but leave them enabled when we integrate over
325  // the full range of one or several variables
326  if (_intOrder > 1 && !(code & 1)) {
327  analVars.removeAll();
328  return 0;
329  }
330  return (code >= 2) ? code : 0;
331 }
332 
333 
334 
335 ////////////////////////////////////////////////////////////////////////////////
336 /// Return integral identified by 'code'. The actual integration
337 /// is deferred to RooDataHist::sum() which implements partial
338 /// or complete summation over the histograms contents
339 
340 Double_t RooHistPdf::analyticalIntegral(Int_t code, const char* rangeName) const
341 {
342  // Simplest scenario, full-range integration over all dependents
343  if (((2 << _histObsList.getSize()) - 1) == code) {
344  return _dataHist->sum(kFALSE);
345  }
346 
347  // Partial integration scenario, retrieve set of variables, calculate partial
348  // sum, figure out integration ranges (if needed)
349  RooArgSet intSet;
350  std::map<const RooAbsArg*, std::pair<Double_t, Double_t> > ranges;
353  Int_t n(0);
354  for (RooAbsArg *pa = 0, *ha = 0; (pa = it.next()) && (ha = jt.next()); ++n) {
355  if (code & (2 << n)) {
356  intSet.add(*ha);
357  }
358  if (!(code & 1)) {
359  RooAbsRealLValue* rlv = dynamic_cast<RooAbsRealLValue*>(pa);
360  if (rlv) {
361  const RooAbsBinning* binning = rlv->getBinningPtr(rangeName);
362  if (rangeName && rlv->hasRange(rangeName)) {
363  ranges[ha] = std::make_pair(
364  rlv->getMin(rangeName), rlv->getMax(rangeName));
365  } else if (binning) {
366  if (!binning->isParameterized()) {
367  ranges[ha] = std::make_pair(
368  binning->lowBound(), binning->highBound());
369  } else {
370  ranges[ha] = std::make_pair(
371  binning->lowBoundFunc()->getVal(), binning->highBoundFunc()->getVal());
372  }
373  }
374  }
375  }
376  // WVE must sync hist slice list values to pdf slice list
377  // Transfer values from
378  if (ha != pa) {
379  pa->syncCache();
380  ha->copyCache(pa,kTRUE);
381  }
382  }
383 
384  Double_t ret = (code & 1) ?
386  _dataHist->sum(intSet,_histObsList,kFALSE,kTRUE, ranges);
387 
388  // cout << "intSet = " << intSet << endl ;
389  // cout << "slice position = " << endl ;
390  // _histObsList.Print("v") ;
391  // cout << "RooHistPdf::ai(" << GetName() << ") code = " << code << " ret = " << ret << endl ;
392 
393  return ret ;
394 }
395 
396 
397 
398 ////////////////////////////////////////////////////////////////////////////////
399 /// Return sampling hint for making curves of (projections) of this function
400 /// as the recursive division strategy of RooCurve cannot deal efficiently
401 /// with the vertical lines that occur in a non-interpolated histogram
402 
403 list<Double_t>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
404 {
405  // No hints are required when interpolation is used
406  if (_intOrder>0) {
407  return 0 ;
408  }
409 
410  // Check that observable is in dataset, if not no hint is generated
411  _histObsIter->Reset() ;
412  _pdfObsIter->Reset() ;
413  RooAbsArg *pdfObs, *histObs, *dhObs(0) ;
414  while ((pdfObs = (RooAbsArg*)_pdfObsIter->Next()) && !dhObs) {
415  histObs = (RooAbsArg*) _histObsIter->Next() ;
416  if (TString(obs.GetName())==pdfObs->GetName()) {
417  dhObs = _dataHist->get()->find(histObs->GetName()) ;
418  }
419  }
420 
421  if (!dhObs) {
422  return 0 ;
423  }
424  RooAbsLValue* lval = dynamic_cast<RooAbsLValue*>(dhObs) ;
425  if (!lval) {
426  return 0 ;
427  }
428 
429  // Retrieve position of all bin boundaries
430 
431  const RooAbsBinning* binning = lval->getBinningPtr(0) ;
432  Double_t* boundaries = binning->array() ;
433 
434  list<Double_t>* hint = new list<Double_t> ;
435 
436  // Widen range slighty
437  xlo = xlo - 0.01*(xhi-xlo) ;
438  xhi = xhi + 0.01*(xhi-xlo) ;
439 
440  Double_t delta = (xhi-xlo)*1e-8 ;
441 
442  // Construct array with pairs of points positioned epsilon to the left and
443  // right of the bin boundaries
444  for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
445  if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
446  hint->push_back(boundaries[i]-delta) ;
447  hint->push_back(boundaries[i]+delta) ;
448  }
449  }
450 
451  return hint ;
452 }
453 
454 
455 
456 ////////////////////////////////////////////////////////////////////////////////
457 /// Return sampling hint for making curves of (projections) of this function
458 /// as the recursive division strategy of RooCurve cannot deal efficiently
459 /// with the vertical lines that occur in a non-interpolated histogram
460 
461 std::list<Double_t>* RooHistPdf::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
462 {
463  // No hints are required when interpolation is used
464  if (_intOrder>0) {
465  return 0 ;
466  }
467 
468  // Check that observable is in dataset, if not no hint is generated
469  RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
470  if (!lvarg) {
471  return 0 ;
472  }
473 
474  // Retrieve position of all bin boundaries
475  const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
476  Double_t* boundaries = binning->array() ;
477 
478  list<Double_t>* hint = new list<Double_t> ;
479 
480  // Construct array with pairs of points positioned epsilon to the left and
481  // right of the bin boundaries
482  for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
483  if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
484  hint->push_back(boundaries[i]) ;
485  }
486  }
487 
488  return hint ;
489 }
490 
491 
492 
493 
494 ////////////////////////////////////////////////////////////////////////////////
495 /// Only handle case of maximum in all variables
496 
498 {
499  RooAbsCollection* common = _pdfObsList.selectCommon(vars) ;
500  if (common->getSize()==_pdfObsList.getSize()) {
501  delete common ;
502  return 1;
503  }
504  delete common ;
505  return 0 ;
506 }
507 
508 
509 ////////////////////////////////////////////////////////////////////////////////
510 
512 {
513  R__ASSERT(code==1) ;
514 
515  Double_t max(-1) ;
516  for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
517  _dataHist->get(i) ;
518  Double_t wgt = _dataHist->weight() ;
519  if (wgt>max) max=wgt ;
520  }
521 
522  return max*1.05 ;
523 }
524 
525 
526 
527 
528 ////////////////////////////////////////////////////////////////////////////////
529 
531 {
532  if (fabs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return kFALSE ;
533  if (dh1.numEntries() != dh2.numEntries()) return kFALSE ;
534  for (int i=0 ; i < dh1.numEntries() ; i++) {
535  dh1.get(i) ;
536  dh2.get(i) ;
537  if (fabs(dh1.weight()-dh2.weight())>1e-8) return kFALSE ;
538  }
539  return kTRUE ;
540 }
541 
542 
543 
544 ////////////////////////////////////////////////////////////////////////////////
545 /// Check if our datahist is already in the workspace
546 
548 {
549  std::list<RooAbsData*> allData = ws.allData() ;
550  std::list<RooAbsData*>::const_iterator iter ;
551  for (iter = allData.begin() ; iter != allData.end() ; ++iter) {
552  // If your dataset is already in this workspace nothing needs to be done
553  if (*iter == _dataHist) {
554  return kFALSE ;
555  }
556  }
557 
558  // Check if dataset with given name already exists
559  RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName()) ;
560 
561  if (wsdata) {
562 
563  // Yes it exists - now check if it is identical to our internal histogram
564  if (wsdata->InheritsFrom(RooDataHist::Class())) {
565 
566  // Check if histograms are identical
567  if (areIdentical((RooDataHist&)*wsdata,*_dataHist)) {
568 
569  // Exists and is of correct type, and identical -- adjust internal pointer to WS copy
570  _dataHist = (RooDataHist*) wsdata ;
571  } else {
572 
573  // not identical, clone rename and import
574  TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
575  Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
576  if (flag) {
577  coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
578  return kTRUE ;
579  }
580  _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
581  }
582 
583  } else {
584 
585  // Exists and is NOT of correct type: clone rename and import
586  TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
587  Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
588  if (flag) {
589  coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
590  return kTRUE ;
591  }
592  _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
593 
594  }
595  return kFALSE ;
596  }
597 
598  // We need to import our datahist into the workspace
600 
601  // Redirect our internal pointer to the copy in the workspace
603  return kFALSE ;
604 }
605 
606 
607 ////////////////////////////////////////////////////////////////////////////////
608 /// Stream an object of class RooHistPdf.
609 
610 void RooHistPdf::Streamer(TBuffer &R__b)
611 {
612  if (R__b.IsReading()) {
613  R__b.ReadClassBuffer(RooHistPdf::Class(),this);
614  // WVE - interim solution - fix proxies here
615  //_proxyList.Clear() ;
616  //registerProxy(_pdfObsList) ;
617  } else {
618  R__b.WriteClassBuffer(RooHistPdf::Class(),this);
619  }
620 }
621 
#define coutE(a)
Definition: RooMsgService.h:35
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual RooAbsReal * lowBoundFunc() const
Definition: RooAbsBinning.h:83
virtual void Reset()=0
Bool_t IsReading() const
Definition: TBuffer.h:81
virtual Bool_t InheritsFrom(const char *classname) const
Returns kTRUE if object inherits from class "classname".
Definition: TObject.cxx:487
RooAbsCollection * selectCommon(const RooAbsCollection &refColl) const
Create a subset of the current collection, consisting only of those elements that are contained as we...
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const
RooFIter fwdIterator() const
#define assert(cond)
Definition: unittest.h:542
virtual RooAbsReal * highBoundFunc() const
Definition: RooAbsBinning.h:87
virtual Bool_t isFundamental() const
Definition: RooAbsArg.h:157
std::list< RooAbsData * > allData() const
Return list of all dataset in the workspace.
Buffer base class used for serializing objects.
Definition: TBuffer.h:40
#define R__ASSERT(e)
Definition: TError.h:98
Basic string class.
Definition: TString.h:137
virtual Double_t getMin(const char *name=0) const
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
const Bool_t kFALSE
Definition: Rtypes.h:92
RooAbsData * embeddedData(const char *name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found...
virtual ~RooHistPdf()
Destructor.
Definition: RooHistPdf.cxx:204
STL namespace.
Iterator abstract base class.
Definition: TIterator.h:32
virtual Int_t numBoundaries() const =0
virtual const RooArgSet * get() const
Definition: RooDataHist.h:77
virtual void syncCache(const RooArgSet *nset=0)=0
const char * Data() const
Definition: TString.h:349
Double_t x[n]
Definition: legend1.C:17
Double_t evaluate() const
Return the current value: The value of the bin enclosing the current coordinates of the observables...
Definition: RooHistPdf.cxx:219
void Class()
Definition: Class.C:29
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
Definition: RooHistPdf.cxx:511
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
TIterator * _pdfObsIter
Definition: RooHistPdf.h:101
Int_t numTypes(const char *=0) const
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const =0
Double_t totVolume() const
Return the total volume spanned by the observables of the RooHistPdf.
Definition: RooHistPdf.cxx:250
std::map< std::string, std::string >::const_iterator iter
Definition: TAlienJob.cxx:54
TIterator * createIterator(Bool_t dir=kIterForward) const
Double_t getVal(const RooArgSet *set=0) const
Definition: RooAbsReal.h:64
RooDataHist * _dataHist
Definition: RooHistPdf.h:99
virtual Double_t highBound() const =0
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...
Definition: RooHistPdf.cxx:461
Double_t sum(Bool_t correctForBinSize, Bool_t inverseCorr=kFALSE) const
Return the sum of the weights of all hist bins.
RooHistPdf()
Default constructor coverity[UNINIT_CTOR].
Definition: RooHistPdf.cxx:51
virtual Bool_t hasRange(const char *name) const
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
RooAbsArg * find(const char *name) const
Find object with given name in list.
virtual Double_t weight() const
Definition: RooDataHist.h:96
ClassImp(RooHistPdf)
virtual void copyCache(const RooAbsArg *source, Bool_t valueOnly=kFALSE, Bool_t setValDirty=kTRUE)=0
Bool_t _unitNorm
Total volume of space (product of ranges of observables)
Definition: RooHistPdf.h:106
Double_t _totVolume
Definition: RooHistPdf.h:105
char * Form(const char *fmt,...)
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
virtual Bool_t inRange(const char *) const
Definition: RooAbsArg.h:289
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Determine integration scenario.
Definition: RooHistPdf.cxx:300
virtual Double_t lowBound() const =0
RooCmdArg Rename(const char *suffix)
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Return integral identified by 'code'.
Definition: RooHistPdf.cxx:340
Bool_t areIdentical(const RooDataHist &dh1, const RooDataHist &dh2)
Definition: RooHistPdf.cxx:530
RooAbsArg * next()
Bool_t _cdfBoundaries
Definition: RooHistPdf.h:104
Int_t _intOrder
Auxiliary class keeping tracking of analytical integration code.
Definition: RooHistPdf.h:103
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
virtual Bool_t isParameterized() const
Definition: RooAbsBinning.h:79
double Double_t
Definition: RtypesCore.h:55
Double_t y[n]
Definition: legend1.C:17
RooCmdArg Embedded(Bool_t flag=kTRUE)
virtual Int_t numEntries() const
Return the number of bins.
virtual Int_t getMaxVal(const RooArgSet &vars) const
Only handle case of maximum in all variables.
Definition: RooHistPdf.cxx:497
static Vc_ALWAYS_INLINE int_v max(const int_v &x, const int_v &y)
Definition: vector.h:440
#define name(a, b)
Definition: linkTestLib0.cpp:5
void setRange(const char *name, Double_t min, Double_t max)
Set range named 'name to [min,max].
Definition: RooRealVar.cxx:446
Bool_t importWorkspaceHook(RooWorkspace &ws)
Check if our datahist is already in the workspace.
Definition: RooHistPdf.cxx:547
virtual Double_t getMax(const char *name=0) const
RooSetProxy _pdfObsList
Definition: RooHistPdf.h:98
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...
Definition: RooHistPdf.cxx:403
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Bool_t import(const RooAbsArg &arg, const RooCmdArg &arg1=RooCmdArg(), const RooCmdArg &arg2=RooCmdArg(), const RooCmdArg &arg3=RooCmdArg(), const RooCmdArg &arg4=RooCmdArg(), const RooCmdArg &arg5=RooCmdArg(), const RooCmdArg &arg6=RooCmdArg(), const RooCmdArg &arg7=RooCmdArg(), const RooCmdArg &arg8=RooCmdArg(), const RooCmdArg &arg9=RooCmdArg())
Import a RooAbsArg object, e.g.
virtual TObject * Next()=0
virtual RooAbsArg * addClone(const RooAbsArg &var, Bool_t silent=kFALSE)
Add clone of specified element to an owning set.
Definition: RooArgSet.cxx:475
virtual Double_t * array() const =0
RooArgSet _histObsList
Definition: RooHistPdf.h:97
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: Rtypes.h:91
const Int_t n
Definition: legend1.C:16
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Add element to non-owning set.
Definition: RooArgSet.cxx:448
TIterator * _histObsIter
Definition: RooHistPdf.h:100
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE)
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
virtual Double_t sumEntries() const