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