Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooHistFunc.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 RooHistFunc.cxx
19\class RooHistFunc
20\ingroup Roofitcore
21
22A real-valued function sampled from a
23multidimensional histogram. The histogram can have an arbitrary number of real or
24discrete dimensions and may have negative values.
25**/
26
27#include "RooHistFunc.h"
28#include "RooDataHist.h"
29#include "RooMsgService.h"
30#include "RooRealVar.h"
31#include "RooCategory.h"
32#include "RooWorkspace.h"
33#include "RooHistPdf.h"
34#include "RooFitImplHelpers.h"
35
36#include "TError.h"
37#include "TBuffer.h"
38
39#include <stdexcept>
40
41
42
43
44////////////////////////////////////////////////////////////////////////////////
45/// Constructor from a RooDataHist. The variable listed in 'vars' control the dimensionality of the
46/// function. Any additional dimensions present in 'dhist' will be projected out. RooDataHist dimensions
47/// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
48/// RooHistFunc neither owns or clone 'dhist' and the user must ensure the input histogram exists
49/// for the entire life span of this function.
50
51RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgSet& vars,
53 RooAbsReal(name,title),
54 _depList("depList","List of dependents",this),
55 _dataHist(const_cast<RooDataHist*>(&dhist)),
56 _codeReg(10),
57 _intOrder(intOrder)
58{
60 _depList.add(vars) ;
61
62 // Verify that vars and dhist.get() have identical contents
63 const RooArgSet* dvars = dhist.get() ;
64 if (vars.size()!=dvars->size()) {
65 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
66 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
67 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
68 }
69
70 for (const auto arg : vars) {
71 if (!dvars->find(arg->GetName())) {
72 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
73 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
74 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
75 }
76 }
77
79}
80
81
82
83////////////////////////////////////////////////////////////////////////////////
84/// Constructor from a RooDataHist. The variable listed in 'vars' control the dimensionality of the
85/// function. Any additional dimensions present in 'dhist' will be projected out. RooDataHist dimensions
86/// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
87/// RooHistFunc neither owns or clone 'dhist' and the user must ensure the input histogram exists
88/// for the entire life span of this function.
89
90RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgList& funcObs, const RooArgList& histObs,
92 RooAbsReal(name,title),
93 _depList("depList","List of dependents",this),
94 _dataHist(const_cast<RooDataHist*>(&dhist)),
95 _codeReg(10),
96 _intOrder(intOrder)
97{
100
101 // Verify that vars and dhist.get() have identical contents
102 const RooArgSet* dvars = dhist.get() ;
103 if (histObs.size()!=dvars->size()) {
104 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
105 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
106 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
107 }
108
109 for (const auto arg : histObs) {
110 if (!dvars->find(arg->GetName())) {
111 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
112 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
113 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
114 }
115 }
116
118}
119
120RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgSet &vars, std::unique_ptr<RooDataHist> dhist,
121 int intOrder)
122 : RooHistFunc{name, title, vars, *dhist, intOrder}
123{
124 initializeOwnedDataHist(std::move(dhist));
125}
126
127RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgList &pdfObs, const RooArgList &histObs,
128 std::unique_ptr<RooDataHist> dhist, int intOrder)
130{
131 initializeOwnedDataHist(std::move(dhist));
132}
133
134
135////////////////////////////////////////////////////////////////////////////////
136/// Copy constructor
137
140 _depList("depList",this,other._depList),
141 _dataHist(other._dataHist),
142 _codeReg(other._codeReg),
143 _intOrder(other._intOrder),
144 _cdfBoundaries(other._cdfBoundaries),
145 _totVolume(other._totVolume),
146 _unitNorm(other._unitNorm)
147{
149
150 _histObsList.addClone(other._histObsList) ;
151}
152
153
154
155////////////////////////////////////////////////////////////////////////////////
156
161
162
163
164
165////////////////////////////////////////////////////////////////////////////////
166/// Return the current value: The value of the bin enclosing the current coordinates
167/// of the dependents, normalized by the histograms contents. Interpolation
168/// is applied if the RooHistFunc is configured to do that
169
171{
172 // Transfer values from
173 if (!_depList.empty()) {
174 for (auto i = 0u; i < _histObsList.size(); ++i) {
175 const auto harg = _histObsList[i];
176 const auto parg = _depList[i];
177
178 if (harg != parg) {
179 parg->syncCache() ;
180 harg->copyCache(parg,true) ;
181 if (!harg->inRange(nullptr)) {
182 return 0 ;
183 }
184 }
185 }
186 }
187
189 return ret ;
190}
191
193{
194 std::span<double> output = ctx.output();
195 std::size_t nEvents = output.size();
196
197 if (_depList.size() == 1) {
198 auto xVals = ctx.at(_depList[0]);
200 return;
201 }
202
203 std::vector<std::span<const double>> inputValues;
204 for (const auto& obs : _depList) {
205 auto realObs = dynamic_cast<const RooAbsReal*>(obs);
206 if (realObs) {
207 inputValues.push_back(ctx.at(realObs));
208 } else {
209 inputValues.emplace_back();
210 }
211 }
212
213 for (std::size_t i = 0; i < nEvents; ++i) {
214 bool skip = false;
215
216 for (auto j = 0u; j < _histObsList.size(); ++j) {
217 const auto histObs = _histObsList[j];
218
219 if (i < inputValues[j].size()) {
220 histObs->setCachedValue(inputValues[j][i], false);
221 if (!histObs->inRange(nullptr)) {
222 skip = true;
223 break;
224 }
225 }
226 }
227
229 }
230}
231
232
233////////////////////////////////////////////////////////////////////////////////
234/// Only handle case of maximum in all variables
235
237{
238 std::unique_ptr<RooAbsCollection> common{_depList.selectCommon(vars)};
239 return common->size() == _depList.size() ? 1 : 0;
240}
241
242////////////////////////////////////////////////////////////////////////////////
243
244double RooHistFunc::maxVal(Int_t code) const
245{
246 R__ASSERT(code==1) ;
247
248 double max(-1) ;
249 for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
250 double wgt = _dataHist->weight(i) ;
251 if (wgt>max) max=wgt ;
252 }
253
254 return max*1.05 ;
255}
256
258 if (_ownedDataHist) return _ownedDataHist.get();
259 _ownedDataHist.reset(static_cast<RooDataHist*>(_dataHist->Clone(newname)));
261 return _dataHist;
262}
263
264////////////////////////////////////////////////////////////////////////////////
265/// Return the total volume spanned by the observables of the RooDataHist
266
268{
269 // Return previously calculated value, if any
270 if (_totVolume>0) {
271 return _totVolume ;
272 }
273 _totVolume = 1. ;
274 for (const auto arg : _depList) {
275 RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
276 if (real) {
277 _totVolume *= (real->getMax()-real->getMin()) ;
278 } else {
279 RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
280 if (cat) {
281 _totVolume *= cat->numTypes() ;
282 }
283 }
284 }
285
286 return _totVolume ;
287}
288
289
290////////////////////////////////////////////////////////////////////////////////
291/// Determine integration scenario. If no interpolation is used,
292/// RooHistFunc can perform all integrals over its dependents
293/// analytically via partial or complete summation of the input
294/// histogram. If interpolation is used, only the integral
295/// over all RooHistPdf observables is implemented.
296
301
302
303////////////////////////////////////////////////////////////////////////////////
304/// Return integral identified by 'code'. The actual integration
305/// is deferred to RooDataHist::sum() which implements partial
306/// or complete summation over the histograms contents
307
308double RooHistFunc::analyticalIntegral(Int_t code, const char* rangeName) const
309{
311}
312
317
318////////////////////////////////////////////////////////////////////////////////
319/// Return sampling hint for making curves of (projections) of this function
320/// as the recursive division strategy of RooCurve cannot deal efficiently
321/// with the vertical lines that occur in a non-interpolated histogram
322
323std::list<double>* RooHistFunc::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
324{
326}
327
328
329////////////////////////////////////////////////////////////////////////////////
330/// Return sampling hint for making curves of (projections) of this function
331/// as the recursive division strategy of RooCurve cannot deal efficiently
332/// with the vertical lines that occur in a non-interpolated histogram
333
334std::list<double>* RooHistFunc::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
335{
336 // No hints are required when interpolation is used
337 if (_intOrder>1) {
338 return nullptr ;
339 }
340
341 // Find histogram observable corresponding to pdf observable
342 RooAbsArg* hobs(nullptr) ;
343 for (auto i = 0u; i < _histObsList.size(); ++i) {
344 const auto harg = _histObsList[i];
345 const auto parg = _depList[i];
346 if (std::string(parg->GetName())==obs.GetName()) {
347 hobs=harg ;
348 }
349 }
350
351 // std::cout << "RooHistFunc::bb(" << GetName() << ") histObs = " << _histObsList << std::endl ;
352 // std::cout << "RooHistFunc::bb(" << GetName() << ") pdfObs = " << _depList << std::endl ;
353
354 RooAbsRealLValue* transform = nullptr;
355 if (!hobs) {
356
357 // Considering alternate: input observable is histogram observable and pdf observable is transformation in terms of it
358 RooAbsArg* pobs = nullptr;
359 for (auto i = 0u; i < _histObsList.size(); ++i) {
360 const auto harg = _histObsList[i];
361 const auto parg = _depList[i];
362 if (std::string(harg->GetName())==obs.GetName()) {
363 pobs=parg ;
364 hobs=harg ;
365 }
366 }
367
368 // Not found, or check that matching pdf observable is an l-value dependent on histogram observable fails
369 if (!hobs || !(pobs->dependsOn(obs) && dynamic_cast<RooAbsRealLValue*>(pobs))) {
370 std::cout << "RooHistFunc::binBoundaries(" << GetName() << ") obs = " << obs.GetName() << " hobs is not found, returning null" << std::endl ;
371 return nullptr ;
372 }
373
374 // Now we are in business - we are in a situation where the pdf observable LV(x), mapping to a histogram observable x
375 // We can return bin boundaries by mapping the histogram boundaties through the inverse of the LV(x) transformation
376 transform = dynamic_cast<RooAbsRealLValue*>(pobs) ;
377 }
378
379
380 // std::cout << "hobs = " << hobs->GetName() << std::endl ;
381 // std::cout << "transform = " << (transform?transform->GetName():"<none>") << std::endl ;
382
383 // Check that observable is in dataset, if not no hint is generated
384 RooAbsArg* xtmp = _dataHist->get()->find(hobs->GetName()) ;
385 if (!xtmp) {
386 std::cout << "RooHistFunc::binBoundaries(" << GetName() << ") hobs = " << hobs->GetName() << " is not found in dataset?" << std::endl ;
387 _dataHist->get()->Print("v") ;
388 return nullptr ;
389 }
390 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(hobs->GetName())) ;
391 if (!lvarg) {
392 std::cout << "RooHistFunc::binBoundaries(" << GetName() << ") hobs = " << hobs->GetName() << " but is not an LV, returning null" << std::endl ;
393 return nullptr ;
394 }
395
396 // Retrieve position of all bin boundaries
397 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
398 double* boundaries = binning->array() ;
399
400 auto hint = new std::list<double> ;
401
402 double delta = (xhi-xlo)*1e-8 ;
403
404 // Construct array with pairs of points positioned epsilon to the left and
405 // right of the bin boundaries
406 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
407 if (boundaries[i]>xlo-delta && boundaries[i]<xhi+delta) {
408
409 double boundary = boundaries[i] ;
410 if (transform) {
411 transform->setVal(boundary) ;
412 //cout << "transform bound " << boundary << " using " << transform->GetName() << " result " << obs.getVal() << std::endl ;
413 hint->push_back(obs.getVal()) ;
414 } else {
415 hint->push_back(boundary) ;
416 }
417 }
418 }
419
420 return hint ;
421}
422
423
424
425////////////////////////////////////////////////////////////////////////////////
426/// Check if our datahist is already in the workspace.
427/// In case of error, return true.
429{
430 // Check if dataset with given name already exists
432
433 if (wsdata) {
434 // If our data is exactly the same, we are done:
435 if (static_cast<RooDataHist*>(wsdata) == _dataHist)
436 return false;
437
438 // Yes it exists - now check if it is identical to our internal histogram
439 if (wsdata->InheritsFrom(RooDataHist::Class())) {
440
441 // Check if histograms are identical
442 if (areIdentical(static_cast<RooDataHist&>(*wsdata),*_dataHist)) {
443
444 // Exists and is of correct type, and identical -- adjust internal pointer to WS copy
445 _dataHist = static_cast<RooDataHist*>(wsdata) ;
446 } else {
447
448 // not identical, clone rename and import
449 auto uniqueName = std::string(_dataHist->GetName()) + "_" + GetName();
451 if (flag) {
452 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << std::endl ;
453 return true ;
454 }
455 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(uniqueName)) ;
456 }
457
458 } else {
459
460 // Exists and is NOT of correct type: clone rename and import
461 auto uniqueName = std::string(_dataHist->GetName()) + "_" + GetName();
463 if (flag) {
464 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << std::endl ;
465 return true ;
466 }
467 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(uniqueName));
468
469 }
470 return false ;
471 }
472
473 // We need to import our datahist into the workspace
475
476 // Redirect our internal pointer to the copy in the workspace
477 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(_dataHist->GetName())) ;
478 return false ;
479}
480
481
482////////////////////////////////////////////////////////////////////////////////
483
485{
486 if (std::abs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return false ;
487 if (dh1.numEntries() != dh2.numEntries()) return false ;
488 for (int i=0 ; i < dh1.numEntries() ; i++) {
489 if (std::abs(dh1.weight(i)-dh2.weight(i))>1e-8) return false ;
490 }
492 if (getColonSeparatedNameString(*dh1.get()) != getColonSeparatedNameString(*dh2.get())) return false ;
493 return true ;
494}
495
496
497
498////////////////////////////////////////////////////////////////////////////////
499/// Stream an object of class RooHistFunc.
500
502{
503 if (R__b.IsReading()) {
504 R__b.ReadClassBuffer(RooHistFunc::Class(),this);
505 // WVE - interim solution - fix proxies here
506 _proxyList.Clear() ;
508 } else {
509 R__b.WriteClassBuffer(RooHistFunc::Class(),this);
510 }
511}
512
513
514////////////////////////////////////////////////////////////////////////////////
515/// Schema evolution: if histObsList wasn't filled from persistence (v1)
516/// then fill it here. Can't be done in regular schema evolution in LinkDef
517/// as _depList content is not guaranteed to be initialized there
518
520{
521 RooAbsReal::ioStreamerPass2(); // call the baseclass method
522
523 if (_histObsList.empty()) {
525 }
526}
527
528
529////////////////////////////////////////////////////////////////////////////////
530/// Compute bin number corresponding to current coordinates.
531/// \return If a bin is not in the current range of the observables, return -1.
533 if (!_depList.empty()) {
534 for (auto i = 0u; i < _histObsList.size(); ++i) {
535 const auto harg = _histObsList[i];
536 const auto parg = _depList[i];
537
538 if (harg != parg) {
539 parg->syncCache() ;
540 harg->copyCache(parg,true) ;
541 if (!harg->inRange(nullptr)) {
542 return -1;
543 }
544 }
545 }
546 }
547
548 return _dataHist->getIndex(_histObsList, true);
549}
550
551
552////////////////////////////////////////////////////////////////////////////////
553/// Compute bin numbers corresponding to all coordinates in `evalData`.
554/// \return Vector of bin numbers. If a bin is not in the current range of the observables, return -1.
555std::vector<Int_t> RooHistFunc::getBins(RooFit::EvalContext & ctx) const {
556 std::vector<std::span<const double>> depData;
557 for (const auto dep : _depList) {
558 auto real = dynamic_cast<const RooAbsReal*>(dep);
559 if (real) {
560 depData.push_back(ctx.at(real));
561 } else {
562 depData.emplace_back();
563 }
564 }
565
566 const auto batchSize = std::max_element(depData.begin(), depData.end(),
567 [](const std::span<const double>& a, const std::span<const double>& b){ return a.size() < b.size(); })->size();
568 std::vector<Int_t> results;
569
570 for (std::size_t evt = 0; evt < batchSize; ++evt) {
571 if (!_depList.empty()) {
572 for (auto i = 0u; i < _histObsList.size(); ++i) {
573 const auto harg = _histObsList[i];
574
575 if (evt < depData[i].size())
576 harg->setCachedValue(depData[i][evt], false);
577
578 if (!harg->inRange(nullptr)) {
579 results.push_back(-1);
580 continue;
581 }
582 }
583 }
584
585 results.push_back(_dataHist->getIndex(_histObsList, true));
586 }
587
588 return results;
589}
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define e(i)
Definition RSha256.hxx:103
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define coutE(a)
#define TRACE_DESTROY
Definition RooTrace.h:24
#define TRACE_CREATE
Definition RooTrace.h:23
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
char name[80]
Definition TGX11.cxx:145
const_iterator begin() const
const_iterator end() const
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:76
RooRefArray _proxyList
Definition RooAbsArg.h:569
void registerProxy(RooArgProxy &proxy)
Register an RooArgProxy in the proxy list.
friend void RooRefArray::Streamer(TBuffer &)
virtual void ioStreamerPass2()
Method called by workspace container to finalize schema evolution issues that cannot be handled in a ...
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: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.
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
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...
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
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
Int_t getIndex(const RooAbsCollection &coord, bool fast=false) const
Calculate bin number of the given coordinates.
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 real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:31
bool _cdfBoundaries
Use boundary conditions for CDFs.
RooDataHist * _dataHist
Unowned pointer to underlying histogram.
double _totVolume
! Total volume of space (product of ranges of observables)
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
bool forceAnalyticalInt(const RooAbsArg &dep) const override
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Return integral identified by 'code'.
~RooHistFunc() override
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...
std::unique_ptr< RooDataHist > _ownedDataHist
! Owned pointer to underlying histogram
void ioStreamerPass2() override
Schema evolution: if histObsList wasn't filled from persistence (v1) then fill it here.
void initializeOwnedDataHist(std::unique_ptr< RooDataHist > &&dataHist)
double evaluate() const override
Return the current value: The value of the bin enclosing the current coordinates of the dependents,...
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...
Int_t getMaxVal(const RooArgSet &vars) const override
Only handle case of maximum in all variables.
double totVolume() const
Get total bin volume spanned by this hist function.
Int_t getBin() const
Compute bin number corresponding to current coordinates.
bool areIdentical(const RooDataHist &dh1, const RooDataHist &dh2)
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
std::vector< Int_t > getBins(RooFit::EvalContext &ctx) const
Compute bin numbers corresponding to all coordinates in evalData.
bool importWorkspaceHook(RooWorkspace &ws) override
Check if our datahist is already in the workspace.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Determine integration scenario.
RooDataHist * cloneAndOwnDataHist(const char *newname="")
Replaces underlying RooDataHist with a clone, which is now owned, and returns the clone.
Int_t _intOrder
Interpolation order.
static TClass * Class()
RooSetProxy _depList
List of observables mapped onto histogram observables.
RooArgSet _histObsList
List of observables defining dimensions of histogram.
bool forceAnalyticalInt(const RooAbsArg &dep) const override
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...
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Determine integration scenario.
Variable that can be changed from the outside.
Definition RooRealVar.h:37
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.
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:49
void Clear(Option_t *option="") override
Remove all objects from the array.
RooCmdArg Rename(const char *suffix)
RooCmdArg Embedded(bool flag=true)
std::string getColonSeparatedNameString(RooArgSet const &argSet, char delim=':')
static void output()