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