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 propability 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,
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{
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,
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{
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)
155{
156 initializeOwnedDataHist(std::move(dhist));
157}
158
159
160////////////////////////////////////////////////////////////////////////////////
161/// Copy constructor
162
165 _pdfObsList("pdfObs",this,other._pdfObsList),
166 _dataHist(other._dataHist),
167 _codeReg(other._codeReg),
168 _intOrder(other._intOrder),
169 _cdfBoundaries(other._cdfBoundaries),
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]);
195}
196
197
198////////////////////////////////////////////////////////////////////////////////
199/// Return the current value: The value of the bin enclosing the current coordinates
200/// of the observables, normalized by the histograms contents. Interpolation
201/// is applied if the RooHistPdf is configured to do that.
202
204{
205 // Transfer values from
206 for (unsigned int i=0; i < _pdfObsList.size(); ++i) {
209
210 if (harg != parg) {
211 parg->syncCache() ;
212 harg->copyCache(parg,true) ;
213 if (!harg->inRange(nullptr)) {
214 return 0 ;
215 }
216 }
217 }
218
220
221 return std::max(ret, 0.0);
222}
223
224////////////////////////////////////////////////////////////////////////////////
225/// Return the total volume spanned by the observables of the RooHistPdf
226
228{
229 // Return previously calculated value, if any
230 if (_totVolume>0) {
231 return _totVolume ;
232 }
233 _totVolume = 1. ;
234
235 for (const auto arg : _histObsList) {
236 RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
237 if (real) {
238 _totVolume *= (real->getMax()-real->getMin()) ;
239 } else {
240 RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
241 if (cat) {
242 _totVolume *= cat->numTypes() ;
243 }
244 }
245 }
246
247 return _totVolume ;
248}
249
250namespace {
251
252bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
253{
254 const RooAbsRealLValue *_x = dynamic_cast<const RooAbsRealLValue*>(&x);
255 const RooAbsRealLValue *_y = dynamic_cast<const RooAbsRealLValue*>(&y);
256 if (!_x || !_y) return false;
257 if (!range || !strlen(range) || !_x->hasRange(range) ||
258 _x->getBinningPtr(range)->isParameterized()) {
259 // parameterized ranges may be full range now, but that might change,
260 // so return false
261 if (range && strlen(range) && _x->getBinningPtr(range)->isParameterized())
262 return false;
263 return (_x->getMin() == _y->getMin() && _x->getMax() == _y->getMax());
264 }
265 return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
266}
267
268bool okayForAnalytical(RooAbsArg const& obs, RooArgSet const& allVars)
269{
270 auto lobs = dynamic_cast<RooAbsRealLValue const*>(&obs);
271 if(lobs == nullptr) return false;
272
273 bool isOkayForAnalyticalInt = false;
274
275 for(RooAbsArg *var : allVars) {
276 if(obs.dependsOn(*var)) {
277 if(!lobs->isJacobianOK(*var)) return false;
279 }
280 }
281
283}
284
285} // namespace
286
287
290 const char* rangeName,
291 RooArgSet const& histObsList,
292 RooArgSet const& pdfObsList,
294{
295 // First make list of pdf observables to histogram observables
296 // and select only those for which the integral is over the full range
297
298 Int_t code = 0;
299 Int_t frcode = 0;
300 for (unsigned int n=0; n < pdfObsList.size() && n < histObsList.size(); ++n) {
301 const auto pa = pdfObsList[n];
302 const auto ha = histObsList[n];
303
304 if (okayForAnalytical(*pa, allVars)) {
305 code |= 2 << n;
306 analVars.add(*pa);
307 if (fullRange(*pa, *ha, rangeName)) {
308 frcode |= 2 << n;
309 }
310 }
311 }
312
313 if (code == frcode) {
314 // integrate over full range of all observables - use bit 0 to indicate
315 // full range integration over all observables
316 code |= 1;
317 }
318
319 // Disable partial analytical integrals if interpolation is used, and we
320 // integrate over sub-ranges, but leave them enabled when we integrate over
321 // the full range of one or several variables
322 if (intOrder > 1 && !(code & 1)) {
323 analVars.removeAll();
324 return 0;
325 }
326 return (code >= 2) ? code : 0;
327}
328
329
331 const char* rangeName,
332 RooArgSet const& histObsList,
333 RooArgSet const& pdfObsList,
334 RooDataHist& dataHist,
335 bool histFuncMode) {
336 // Simplest scenario, full-range integration over all dependents
337 if (((2 << histObsList.size()) - 1) == code) {
338 return dataHist.sum(histFuncMode);
339 }
340
341 // Partial integration scenario, retrieve set of variables, calculate partial
342 // sum, figure out integration ranges (if needed)
344 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
345 for (unsigned int n=0; n < pdfObsList.size() && n < histObsList.size(); ++n) {
346 const auto pa = pdfObsList[n];
347 const auto ha = histObsList[n];
348
349 if (code & (2 << n)) {
350 intSet.add(*ha);
351 }
352 if (!(code & 1)) {
353 ranges[ha] = RooHelpers::getRangeOrBinningInterval(pa, rangeName);
354 }
355 // WVE must sync hist slice list values to pdf slice list
356 // Transfer values from
357 if (ha != pa) {
358 pa->syncCache();
359 ha->copyCache(pa,true);
360 }
361 }
362
363 double ret = (code & 1) ? dataHist.sum(intSet,histObsList,true,!histFuncMode) :
365
366 return ret ;
367}
368
369////////////////////////////////////////////////////////////////////////////////
370/// Determine integration scenario. If no interpolation is used,
371/// RooHistPdf can perform all integrals over its dependents
372/// analytically via partial or complete summation of the input
373/// histogram. If interpolation is used on the integral over
374/// all histogram observables is supported
375
377{
378 return getAnalyticalIntegral(allVars, analVars, rangeName, _histObsList, _pdfObsList, _intOrder);
379}
380
381
382////////////////////////////////////////////////////////////////////////////////
383/// Return integral identified by 'code'. The actual integration
384/// is deferred to RooDataHist::sum() which implements partial
385/// or complete summation over the histograms contents.
386
387double RooHistPdf::analyticalIntegral(Int_t code, const char* rangeName) const
388{
389 return analyticalIntegral(code, rangeName, _histObsList, _pdfObsList, *_dataHist, false);
390}
391
392
394{
395 bool isOkayForAnalyticalInt = false;
396
397 for (RooAbsArg * obs : pdfObsList) {
398 if(obs->dependsOn(dep)) {
399 // If the observable doesn't depend linearly on the integration
400 // variable we will not do analytical integration.
401 auto lvalue = dynamic_cast<RooAbsRealLValue const*>(obs);
402 if(!(lvalue && lvalue->isJacobianOK(dep))) return false;
404 }
405 }
406
408}
409
410
415
416
417////////////////////////////////////////////////////////////////////////////////
418/// Return sampling hint for making curves of (projections) of this function
419/// as the recursive division strategy of RooCurve cannot deal efficiently
420/// with the vertical lines that occur in a non-interpolated histogram
421
422std::list<double>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
423{
425}
426
427
428std::list<double>* RooHistPdf::plotSamplingHint(RooDataHist const& dataHist,
429 RooArgSet const& pdfObsList,
430 RooArgSet const& histObsList,
431 int intOrder,
432 RooAbsRealLValue& obs,
433 double xlo,
434 double xhi)
435{
436 // No hints are required when interpolation is used
437 if (intOrder>0) {
438 return nullptr;
439 }
440
441 // Check that observable is in dataset, if not no hint is generated
442 RooAbsArg* dhObs = nullptr;
443 for (unsigned int i=0; i < pdfObsList.size(); ++i) {
446 if (std::string(obs.GetName())==pdfObs->GetName()) {
447 dhObs = dataHist.get()->find(histObs->GetName()) ;
448 break;
449 }
450 }
451
452 if (!dhObs) {
453 return nullptr;
454 }
455 RooAbsLValue* lval = dynamic_cast<RooAbsLValue*>(dhObs) ;
456 if (!lval) {
457 return nullptr;
458 }
459
460 // Retrieve position of all bin boundaries
461
462 const RooAbsBinning* binning = lval->getBinningPtr(nullptr);
463 std::span<const double> boundaries{binning->array(), static_cast<std::size_t>(binning->numBoundaries())};
464
465 // Use the helper function from RooCurve to make sure to get sampling hints
466 // that work with the RooFitPlotting.
467 return RooCurve::plotSamplingHintForBinBoundaries(boundaries, xlo, xhi);
468}
469
470
471////////////////////////////////////////////////////////////////////////////////
472/// Return sampling hint for making curves of (projections) of this function
473/// as the recursive division strategy of RooCurve cannot deal efficiently
474/// with the vertical lines that occur in a non-interpolated histogram
475
476std::list<double>* RooHistPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
477{
478 // No hints are required when interpolation is used
479 if (_intOrder>0) {
480 return nullptr;
481 }
482
483 // Check that observable is in dataset, if not no hint is generated
484 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
485 if (!lvarg) {
486 return nullptr ;
487 }
488
489 // Retrieve position of all bin boundaries
490 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
491 double* boundaries = binning->array() ;
492
493 auto hint = new std::list<double> ;
494
495 // Construct array with pairs of points positioned epsilon to the left and
496 // right of the bin boundaries
497 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
498 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
499 hint->push_back(boundaries[i]) ;
500 }
501 }
502
503 return hint ;
504}
505
506
507
508
509////////////////////////////////////////////////////////////////////////////////
510/// Only handle case of maximum in all variables
511
513{
514 std::unique_ptr<RooAbsCollection> common{_pdfObsList.selectCommon(vars)};
515 if (common->size()==_pdfObsList.size()) {
516 return 1;
517 }
518 return 0 ;
519}
520
521
522////////////////////////////////////////////////////////////////////////////////
523
524double RooHistPdf::maxVal(Int_t code) const
525{
526 R__ASSERT(code==1) ;
527
528 double max(-1) ;
529 for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
530 _dataHist->get(i) ;
531 double wgt = _dataHist->weight() ;
532 if (wgt>max) max=wgt ;
533 }
534
535 return max*1.05 ;
536}
537
538
539
540
541////////////////////////////////////////////////////////////////////////////////
542
544{
545 if (std::abs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return false ;
546 if (dh1.numEntries() != dh2.numEntries()) return false ;
547 for (int i=0 ; i < dh1.numEntries() ; i++) {
548 dh1.get(i) ;
549 dh2.get(i) ;
550 if (std::abs(dh1.weight()-dh2.weight())>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
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();
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();
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()) {
622 R__b.ReadClassBuffer(RooHistPdf::Class(),this);
623 // WVE - interim solution - fix proxies here
624 //_proxyList.Clear() ;
625 //registerProxy(_pdfObsList) ;
626 } else {
627 R__b.WriteClassBuffer(RooHistPdf::Class(),this);
628 }
629}
#define e(i)
Definition RSha256.hxx:103
#define coutE(a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
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.
friend void RooRefArray::Streamer(TBuffer &)
Abstract base class for RooRealVar binning definitions.
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)
Storage_t::size_type size() const
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
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.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
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
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
Object to represent discrete states.
Definition RooCategory.h:28
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
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:892
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
double sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
void weights(double *output, std::span< double const > xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of RooDataHist::weight() for one dimensional histograms with up to one dimension...
static TClass * Class()
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
Definition RooDataHist.h:61
double weight(std::size_t i) const
Return weight of i-th bin.
double weightFast(const RooArgSet &bin, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A faster version of RooDataHist::weight that assumes the passed arguments are aligned with the histog...
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
A propability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
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'.
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
void setRange(const char *name, double min, double max)
Set a fit or plotting range.
Persistable container for RooFit projects.
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
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
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
std::pair< double, double > getRangeOrBinningInterval(RooAbsArg const *arg, const char *rangeName)
static void output()