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