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