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 "RooHelpers.h"
35
36#include "TError.h"
37#include "TBuffer.h"
38
39#include <stdexcept>
40
41using namespace std;
42
44;
45
46
47
48////////////////////////////////////////////////////////////////////////////////
49/// Default constructor
50
52 _dataHist(0),
53 _intOrder(0),
54 _cdfBoundaries(kFALSE),
55 _totVolume(0),
56 _unitNorm(kFALSE)
57{
59}
60
61
62////////////////////////////////////////////////////////////////////////////////
63/// Constructor from a RooDataHist. The variable listed in 'vars' control the dimensionality of the
64/// function. Any additional dimensions present in 'dhist' will be projected out. RooDataHist dimensions
65/// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
66/// RooHistFunc neither owns or clone 'dhist' and the user must ensure the input histogram exists
67/// for the entire life span of this function.
68
69RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgSet& vars,
70 const RooDataHist& dhist, Int_t intOrder) :
71 RooAbsReal(name,title),
72 _depList("depList","List of dependents",this),
73 _dataHist((RooDataHist*)&dhist),
74 _codeReg(10),
75 _intOrder(intOrder),
76 _cdfBoundaries(kFALSE),
77 _totVolume(0),
78 _unitNorm(kFALSE)
79{
81 _depList.add(vars) ;
82
83 // Verify that vars and dhist.get() have identical contents
84 const RooArgSet* dvars = dhist.get() ;
85 if (vars.getSize()!=dvars->getSize()) {
86 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
87 << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
88 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
89 }
90
91 for (const auto arg : vars) {
92 if (!dvars->find(arg->GetName())) {
93 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
94 << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
95 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
96 }
97 }
98
100}
101
102
103
104////////////////////////////////////////////////////////////////////////////////
105/// Constructor from a RooDataHist. The variable listed in 'vars' control the dimensionality of the
106/// function. Any additional dimensions present in 'dhist' will be projected out. RooDataHist dimensions
107/// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
108/// RooHistFunc neither owns or clone 'dhist' and the user must ensure the input histogram exists
109/// for the entire life span of this function.
110
111RooHistFunc::RooHistFunc(const char *name, const char *title, const RooArgList& funcObs, const RooArgList& histObs,
112 const RooDataHist& dhist, Int_t intOrder) :
113 RooAbsReal(name,title),
114 _depList("depList","List of dependents",this),
115 _dataHist((RooDataHist*)&dhist),
116 _codeReg(10),
117 _intOrder(intOrder),
118 _cdfBoundaries(kFALSE),
119 _totVolume(0),
120 _unitNorm(kFALSE)
121{
122 _histObsList.addClone(histObs) ;
123 _depList.add(funcObs) ;
124
125 // Verify that vars and dhist.get() have identical contents
126 const RooArgSet* dvars = dhist.get() ;
127 if (histObs.getSize()!=dvars->getSize()) {
128 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
129 << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
130 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
131 }
132
133 for (const auto arg : histObs) {
134 if (!dvars->find(arg->GetName())) {
135 coutE(InputArguments) << "RooHistFunc::ctor(" << GetName()
136 << ") ERROR variable list and RooDataHist must contain the same variables." << endl ;
137 throw std::invalid_argument("RooHistFunc: ERROR variable list and RooDataHist must contain the same variables.");
138 }
139 }
140
142}
143
144
145
146////////////////////////////////////////////////////////////////////////////////
147/// Copy constructor
148
149RooHistFunc::RooHistFunc(const RooHistFunc& other, const char* name) :
150 RooAbsReal(other,name),
151 _depList("depList",this,other._depList),
152 _dataHist(other._dataHist),
153 _codeReg(other._codeReg),
154 _intOrder(other._intOrder),
155 _cdfBoundaries(other._cdfBoundaries),
156 _totVolume(other._totVolume),
157 _unitNorm(other._unitNorm)
158{
160
162}
163
164
165
166////////////////////////////////////////////////////////////////////////////////
167
169{
171}
172
173
174
175
176////////////////////////////////////////////////////////////////////////////////
177/// Return the current value: The value of the bin enclosing the current coordinates
178/// of the dependents, normalized by the histograms contents. Interpolation
179/// is applied if the RooHistFunc is configured to do that
180
182{
183 // Transfer values from
184 if (_depList.getSize()>0) {
185 for (auto i = 0u; i < _histObsList.size(); ++i) {
186 const auto harg = _histObsList[i];
187 const auto parg = _depList[i];
188
189 if (harg != parg) {
190 parg->syncCache() ;
191 harg->copyCache(parg,kTRUE) ;
192 if (!harg->inRange(0)) {
193 return 0 ;
194 }
195 }
196 }
197 }
198
200 return ret ;
201}
202
203
204void RooHistFunc::computeBatch(cudaStream_t*, double* output, size_t size, RooFit::Detail::DataMap const& dataMap) const {
205 std::vector<RooSpan<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 RooAbsCollection* common = _depList.selectCommon(vars) ;
242 if (common->getSize()==_depList.getSize()) {
243 delete common ;
244 return 1;
245 }
246 delete common ;
247 return 0 ;
248}
249
250////////////////////////////////////////////////////////////////////////////////
251
253{
254 R__ASSERT(code==1) ;
255
256 Double_t max(-1) ;
257 for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
258 _dataHist->get(i) ;
259 Double_t wgt = _dataHist->weight() ;
260 if (wgt>max) max=wgt ;
261 }
262
263 return max*1.05 ;
264}
265
266////////////////////////////////////////////////////////////////////////////////
267/// Return the total volume spanned by the observables of the RooDataHist
268
270{
271 // Return previously calculated value, if any
272 if (_totVolume>0) {
273 return _totVolume ;
274 }
275 _totVolume = 1. ;
276 for (const auto arg : _depList) {
277 RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
278 if (real) {
279 _totVolume *= (real->getMax()-real->getMin()) ;
280 } else {
281 RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
282 if (cat) {
283 _totVolume *= cat->numTypes() ;
284 }
285 }
286 }
287
288 return _totVolume ;
289}
290
291
292////////////////////////////////////////////////////////////////////////////////
293/// Determine integration scenario. If no interpolation is used,
294/// RooHistFunc can perform all integrals over its dependents
295/// analytically via partial or complete summation of the input
296/// histogram. If interpolation is used, only the integral
297/// over all RooHistPdf observables is implemented.
298
299Int_t RooHistFunc::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
300{
301 return RooHistPdf::getAnalyticalIntegral(allVars, analVars, rangeName, _histObsList, _depList, _intOrder);
302}
303
304
305////////////////////////////////////////////////////////////////////////////////
306/// Return integral identified by 'code'. The actual integration
307/// is deferred to RooDataHist::sum() which implements partial
308/// or complete summation over the histograms contents
309
310Double_t RooHistFunc::analyticalIntegral(Int_t code, const char* rangeName) const
311{
312 return RooHistPdf::analyticalIntegral(code, rangeName, _histObsList, _depList, *_dataHist, true);
313}
314
315
316////////////////////////////////////////////////////////////////////////////////
317/// Return sampling hint for making curves of (projections) of this function
318/// as the recursive division strategy of RooCurve cannot deal efficiently
319/// with the vertical lines that occur in a non-interpolated histogram
320
322{
323 // No hints are required when interpolation is used
324 if (_intOrder>1) {
325 return 0 ;
326 }
327
328
329 // Find histogram observable corresponding to pdf observable
330 RooAbsArg* hobs(0) ;
331 for (auto i = 0u; i < _histObsList.size(); ++i) {
332 const auto harg = _histObsList[i];
333 const auto parg = _depList[i];
334 if (string(parg->GetName())==obs.GetName()) {
335 hobs=harg ;
336 }
337 }
338 if (!hobs) {
339 return 0 ;
340 }
341
342 // Check that observable is in dataset, if not no hint is generated
343 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(hobs->GetName())) ;
344 if (!lvarg) {
345 return 0 ;
346 }
347
348 // Retrieve position of all bin boundaries
349 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
350 Double_t* boundaries = binning->array() ;
351
352 list<Double_t>* hint = new list<Double_t> ;
353
354 // Widen range slighty
355 xlo = xlo - 0.01*(xhi-xlo) ;
356 xhi = xhi + 0.01*(xhi-xlo) ;
357
358 Double_t delta = (xhi-xlo)*1e-8 ;
359
360 // Construct array with pairs of points positioned epsilon to the left and
361 // right of the bin boundaries
362 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
363 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
364 hint->push_back(boundaries[i]-delta) ;
365 hint->push_back(boundaries[i]+delta) ;
366 }
367 }
368
369 return hint ;
370}
371
372
373////////////////////////////////////////////////////////////////////////////////
374/// Return sampling hint for making curves of (projections) of this function
375/// as the recursive division strategy of RooCurve cannot deal efficiently
376/// with the vertical lines that occur in a non-interpolated histogram
377
378std::list<Double_t>* RooHistFunc::binBoundaries(RooAbsRealLValue& obs, Double_t xlo, Double_t xhi) const
379{
380 // No hints are required when interpolation is used
381 if (_intOrder>1) {
382 return 0 ;
383 }
384
385 // Find histogram observable corresponding to pdf observable
386 RooAbsArg* hobs(0) ;
387 for (auto i = 0u; i < _histObsList.size(); ++i) {
388 const auto harg = _histObsList[i];
389 const auto parg = _depList[i];
390 if (string(parg->GetName())==obs.GetName()) {
391 hobs=harg ;
392 }
393 }
394
395 // cout << "RooHistFunc::bb(" << GetName() << ") histObs = " << _histObsList << endl ;
396 // cout << "RooHistFunc::bb(" << GetName() << ") pdfObs = " << _depList << endl ;
397
398 RooAbsRealLValue* transform(0) ;
399 if (!hobs) {
400
401 // Considering alternate: input observable is histogram observable and pdf observable is transformation in terms of it
402 RooAbsArg* pobs(0) ;
403 for (auto i = 0u; i < _histObsList.size(); ++i) {
404 const auto harg = _histObsList[i];
405 const auto parg = _depList[i];
406 if (string(harg->GetName())==obs.GetName()) {
407 pobs=parg ;
408 hobs=harg ;
409 }
410 }
411
412 // Not found, or check that matching pdf observable is an l-value dependent on histogram observable fails
413 if (!hobs || !(pobs->dependsOn(obs) && dynamic_cast<RooAbsRealLValue*>(pobs))) {
414 cout << "RooHistFunc::binBoundaries(" << GetName() << ") obs = " << obs.GetName() << " hobs is not found, returning null" << endl ;
415 return 0 ;
416 }
417
418 // Now we are in business - we are in a situation where the pdf observable LV(x), mapping to a histogram observable x
419 // We can return bin boundaries by mapping the histogram boundaties through the inverse of the LV(x) transformation
420 transform = dynamic_cast<RooAbsRealLValue*>(pobs) ;
421 }
422
423
424 // cout << "hobs = " << hobs->GetName() << endl ;
425 // cout << "transform = " << (transform?transform->GetName():"<none>") << endl ;
426
427 // Check that observable is in dataset, if not no hint is generated
428 RooAbsArg* xtmp = _dataHist->get()->find(hobs->GetName()) ;
429 if (!xtmp) {
430 cout << "RooHistFunc::binBoundaries(" << GetName() << ") hobs = " << hobs->GetName() << " is not found in dataset?" << endl ;
431 _dataHist->get()->Print("v") ;
432 return 0 ;
433 }
434 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(hobs->GetName())) ;
435 if (!lvarg) {
436 cout << "RooHistFunc::binBoundaries(" << GetName() << ") hobs = " << hobs->GetName() << " but is not an LV, returning null" << endl ;
437 return 0 ;
438 }
439
440 // Retrieve position of all bin boundaries
441 const RooAbsBinning* binning = lvarg->getBinningPtr(0) ;
442 Double_t* boundaries = binning->array() ;
443
444 list<Double_t>* hint = new list<Double_t> ;
445
446 Double_t delta = (xhi-xlo)*1e-8 ;
447
448 // Construct array with pairs of points positioned epsilon to the left and
449 // right of the bin boundaries
450 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
451 if (boundaries[i]>xlo-delta && boundaries[i]<xhi+delta) {
452
453 Double_t boundary = boundaries[i] ;
454 if (transform) {
455 transform->setVal(boundary) ;
456 //cout << "transform bound " << boundary << " using " << transform->GetName() << " result " << obs.getVal() << endl ;
457 hint->push_back(obs.getVal()) ;
458 } else {
459 hint->push_back(boundary) ;
460 }
461 }
462 }
463
464 return hint ;
465}
466
467
468
469////////////////////////////////////////////////////////////////////////////////
470/// Check if our datahist is already in the workspace.
471/// In case of error, return true.
473{
474 // Check if dataset with given name already exists
475 RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName()) ;
476
477 if (wsdata) {
478 // If our data is exactly the same, we are done:
479 if (static_cast<RooDataHist*>(wsdata) == _dataHist)
480 return false;
481
482 // Yes it exists - now check if it is identical to our internal histogram
483 if (wsdata->InheritsFrom(RooDataHist::Class())) {
484
485 // Check if histograms are identical
486 if (areIdentical((RooDataHist&)*wsdata,*_dataHist)) {
487
488 // Exists and is of correct type, and identical -- adjust internal pointer to WS copy
489 _dataHist = (RooDataHist*) wsdata ;
490 } else {
491
492 // not identical, clone rename and import
493 TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
494 Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
495 if (flag) {
496 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
497 return kTRUE ;
498 }
499 _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
500 }
501
502 } else {
503
504 // Exists and is NOT of correct type: clone rename and import
505 TString uniqueName = Form("%s_%s",_dataHist->GetName(),GetName()) ;
506 Bool_t flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.Data()),RooFit::Embedded()) ;
507 if (flag) {
508 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << endl ;
509 return kTRUE ;
510 }
511 _dataHist = (RooDataHist*) ws.embeddedData(uniqueName.Data()) ;
512
513 }
514 return kFALSE ;
515 }
516
517 // We need to import our datahist into the workspace
518 ws.import(*_dataHist,RooFit::Embedded()) ;
519
520 // Redirect our internal pointer to the copy in the workspace
521 _dataHist = (RooDataHist*) ws.embeddedData(_dataHist->GetName()) ;
522 return kFALSE ;
523}
524
525
526////////////////////////////////////////////////////////////////////////////////
527
529{
530 if (fabs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return kFALSE ;
531 if (dh1.numEntries() != dh2.numEntries()) return kFALSE ;
532 for (int i=0 ; i < dh1.numEntries() ; i++) {
533 dh1.get(i) ;
534 dh2.get(i) ;
535 if (fabs(dh1.weight()-dh2.weight())>1e-8) return kFALSE ;
536 }
538 if (getColonSeparatedNameString(*dh1.get()) != getColonSeparatedNameString(*dh2.get())) return kFALSE ;
539 return kTRUE ;
540}
541
542
543
544////////////////////////////////////////////////////////////////////////////////
545/// Stream an object of class RooHistFunc.
546
548{
549 if (R__b.IsReading()) {
550 R__b.ReadClassBuffer(RooHistFunc::Class(),this);
551 // WVE - interim solution - fix proxies here
552 _proxyList.Clear() ;
554 } else {
555 R__b.WriteClassBuffer(RooHistFunc::Class(),this);
556 }
557}
558
559
560////////////////////////////////////////////////////////////////////////////////
561/// Schema evolution: if histObsList wasn't filled from persistence (v1)
562/// then fill it here. Can't be done in regular schema evolution in LinkDef
563/// as _depList content is not guaranteed to be initialized there
564
566{
567 if (_histObsList.getSize()==0) {
569 }
570}
571
572
573////////////////////////////////////////////////////////////////////////////////
574/// Compute bin number corresponding to current coordinates.
575/// \return If a bin is not in the current range of the observables, return -1.
577 if (!_depList.empty()) {
578 for (auto i = 0u; i < _histObsList.size(); ++i) {
579 const auto harg = _histObsList[i];
580 const auto parg = _depList[i];
581
582 if (harg != parg) {
583 parg->syncCache() ;
584 harg->copyCache(parg,kTRUE) ;
585 if (!harg->inRange(nullptr)) {
586 return -1;
587 }
588 }
589 }
590 }
591
592 return _dataHist->getIndex(_histObsList, true);
593}
594
595
596////////////////////////////////////////////////////////////////////////////////
597/// Compute bin numbers corresponding to all coordinates in `evalData`.
598/// \return Vector of bin numbers. If a bin is not in the current range of the observables, return -1.
599std::vector<Int_t> RooHistFunc::getBins(RooFit::Detail::DataMap const& dataMap) const {
600 std::vector<RooSpan<const double>> depData;
601 for (const auto dep : _depList) {
602 auto real = dynamic_cast<const RooAbsReal*>(dep);
603 if (real) {
604 depData.push_back(dataMap.at(real));
605 } else {
606 depData.emplace_back(nullptr, 0);
607 }
608 }
609
610 const auto batchSize = std::max_element(depData.begin(), depData.end(),
611 [](const RooSpan<const double>& a, const RooSpan<const double>& b){ return a.size() < b.size(); })->size();
612 std::vector<Int_t> results;
613
614 for (std::size_t evt = 0; evt < batchSize; ++evt) {
615 if (!_depList.empty()) {
616 for (auto i = 0u; i < _histObsList.size(); ++i) {
617 const auto harg = _histObsList[i];
618
619 if (evt < depData[i].size())
620 harg->setCachedValue(depData[i][evt], false);
621
622 if (!harg->inRange(nullptr)) {
623 results.push_back(-1);
624 continue;
625 }
626 }
627 }
628
629 results.push_back(_dataHist->getIndex(_histObsList, true));
630 }
631
632 return results;
633}
#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
const Bool_t kFALSE
Definition RtypesCore.h:101
const Bool_t kTRUE
Definition RtypesCore.h:100
#define ClassImp(name)
Definition Rtypes.h:364
#define R__ASSERT(e)
Definition TError.h:118
char name[80]
Definition TGX11.cxx:110
char * Form(const char *fmt,...)
RooAbsArg is the common abstract base class for objects that represent a value and a "shape" in RooFi...
Definition RooAbsArg.h:69
RooRefArray _proxyList
Definition RooAbsArg.h:655
Bool_t dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=0, Bool_t valueOnly=kFALSE) 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 &)
RooAbsBinning is the abstract base class for RooRealVar binning definitions.
virtual Double_t * array() const =0
virtual Int_t numBoundaries() const =0
Int_t numTypes(const char *=0) const
Return number of types defined (in range named rangeName if rangeName!=0)
RooAbsCollection is an abstract container object that can hold multiple RooAbsArg objects.
Int_t getSize() const
virtual RooAbsArg * addClone(const RooAbsArg &var, Bool_t silent=kFALSE)
Add a clone of the specified argument to list.
Storage_t::size_type size() const
virtual void Print(Option_t *options=0) const
This method must be overridden when a class wants to print itself.
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.
RooAbsData is the common abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:82
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_t getMax(const char *name=0) const
Get maximum of currently defined range.
virtual void setVal(Double_t value)=0
Set the current value of the object. Needs to be overridden by implementations.
virtual Double_t getMin(const char *name=0) const
Get miniminum of currently defined range.
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:64
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:94
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:35
RooCategory is an object to represent discrete states.
Definition RooCategory.h:27
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:45
Int_t getIndex(const RooAbsCollection &coord, Bool_t fast=false) const
Calculate bin number of the given coordinates.
double weight(std::size_t i) const
Return weight of i-th bin.
Double_t sumEntries() const override
Sum the weights of all bins.
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...
Int_t numEntries() const override
Return the number of bins.
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:84
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition DataMap.h:88
RooHistFunc implements a real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:30
virtual Int_t getMaxVal(const RooArgSet &vars) const
Only handle case of maximum in all variables.
Bool_t _cdfBoundaries
RooDataHist * _dataHist
double evaluate() const
Return the current value: The value of the bin enclosing the current coordinates of the dependents,...
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Determine integration scenario.
Bool_t importWorkspaceHook(RooWorkspace &ws)
Check if our datahist is already in the workspace.
std::vector< Int_t > getBins(RooFit::Detail::DataMap const &dataMap) const
Compute bin numbers corresponding to all coordinates in evalData.
virtual std::list< Double_t > * binBoundaries(RooAbsRealLValue &, Double_t, Double_t) const
Return sampling hint for making curves of (projections) of this function as the recursive division st...
RooHistFunc()
Default constructor.
Double_t totVolume() const
Get total bin volume spanned by this hist function.
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Return integral identified by 'code'.
virtual void ioStreamerPass2()
Schema evolution: if histObsList wasn't filled from persistence (v1) then fill it here.
Bool_t areIdentical(const RooDataHist &dh1, const RooDataHist &dh2)
Int_t getBin() const
Compute bin number corresponding to current coordinates.
virtual ~RooHistFunc()
virtual Double_t maxVal(Int_t code) const
Return maximum value for set of observables identified by code assigned in getMaxVal.
void computeBatch(cudaStream_t *, double *output, size_t size, RooFit::Detail::DataMap const &) const
Base function for computing multiple values of a RooAbsReal.
Int_t _intOrder
Auxiliary class keeping tracking of analytical integration code.
Double_t _totVolume
RooSetProxy _depList
RooArgSet _histObsList
virtual std::list< Double_t > * plotSamplingHint(RooAbsRealLValue &obs, Double_t xlo, Double_t xhi) const
Return sampling hint for making curves of (projections) of this function as the recursive division st...
static Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName, RooArgSet const &histObsList, RooSetProxy const &pdfObsList, Int_t intOrder)
static Double_t analyticalIntegral(Int_t code, const char *rangeName, RooArgSet const &histObsList, RooSetProxy const &pdfObsList, RooDataHist &dataHist, bool histFuncMode)
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:39
virtual Bool_t add(const RooAbsArg &var, Bool_t silent=kFALSE) override
Overloaded RooArgSet::add() method inserts 'var' into set and registers 'var' as server to owner with...
A simple container to hold a batch of data values.
Definition RooSpan.h:34
The RooWorkspace is a persistable container for RooFit projects.
Buffer base class used for serializing objects.
Definition TBuffer.h:43
virtual Int_t ReadClassBuffer(const TClass *cl, void *pointer, const TClass *onfile_class=0)=0
Bool_t IsReading() const
Definition TBuffer.h:86
virtual Int_t WriteClassBuffer(const TClass *cl, void *pointer)=0
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
virtual void Clear(Option_t *option="")
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:515
Basic string class.
Definition TString.h:136
const char * Data() const
Definition TString.h:369
RooCmdArg Embedded(Bool_t flag=kTRUE)
RooCmdArg Rename(const char *suffix)
std::string getColonSeparatedNameString(RooArgSet const &argSet)
Create a string with all sorted names of RooArgSet elements separated by colons.
static void output(int code)
Definition gifencode.c:226
void ws()
Definition ws.C:66