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