Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooHistPdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2005, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18\file RooHistPdf.cxx
19\class RooHistPdf
20\ingroup Roofitcore
21
22A propability density function sampled from a
23multidimensional histogram. The histogram distribution is explicitly
24normalized by RooHistPdf and can have an arbitrary number of real or
25discrete dimensions.
26**/
27
28#include "Riostream.h"
29
30#include "RooCategory.h"
31#include "RooCurve.h"
32#include "RooDataHist.h"
33#include "RooFitImplHelpers.h"
34#include "RooGlobalFunc.h"
35#include "RooHistPdf.h"
36#include "RooMsgService.h"
37#include "RooRealVar.h"
38#include "RooUniformBinning.h"
39#include "RooWorkspace.h"
40
41#include "TError.h"
42#include "TBuffer.h"
43
44
46
47
48////////////////////////////////////////////////////////////////////////////////
49/// Constructor from a RooDataHist. RooDataHist dimensions
50/// can be either real or discrete. See RooDataHist::RooDataHist for details on the binning.
51/// RooHistPdf neither owns or clone 'dhist' and the user must ensure the input histogram exists
52/// for the entire life span of this PDF.
53
54RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgSet& vars,
55 const RooDataHist& dhist, Int_t intOrder) :
56 RooAbsPdf(name,title),
57 _pdfObsList("pdfObs","List of p.d.f. observables",this),
58 _dataHist(const_cast<RooDataHist*>(&dhist)),
59 _codeReg(10),
60 _intOrder(intOrder)
61{
63 _pdfObsList.add(vars) ;
64
65 // Verify that vars and dhist.get() have identical contents
66 const RooArgSet* dvars = dhist.get() ;
67 if (vars.size()!=dvars->size()) {
68 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
69 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
70 assert(0) ;
71 }
72 for (const auto arg : vars) {
73 if (!dvars->find(arg->GetName())) {
74 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
75 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
76 assert(0) ;
77 }
78 }
79
80
81 // Adjust ranges of _histObsList to those of _dataHist
82 for (const auto hobs : _histObsList) {
83 // Guaranteed to succeed, since checked above in constructor
84 RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
85 RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
86 if (dhreal){
87 (static_cast<RooRealVar*>(hobs))->setRange(dhreal->getMin(),dhreal->getMax()) ;
88 }
89 }
90
91}
92
93
94
95
96////////////////////////////////////////////////////////////////////////////////
97/// Constructor from a RooDataHist. The first list of observables are the p.d.f.
98/// observables, which may any RooAbsReal (function or variable). The second list
99/// are the corresponding observables in the RooDataHist which must be of type
100/// RooRealVar or RooCategory This constructor thus allows to apply a coordinate transformation
101/// on the histogram data to be applied.
102
103RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs,
104 const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder) :
105 RooAbsPdf(name,title),
106 _pdfObsList("pdfObs","List of p.d.f. observables",this),
107 _dataHist(const_cast<RooDataHist*>(&dhist)),
108 _codeReg(10),
109 _intOrder(intOrder)
110{
111 _histObsList.addClone(histObs) ;
112 _pdfObsList.add(pdfObs) ;
113
114 // Verify that vars and dhist.get() have identical contents
115 const RooArgSet* dvars = dhist.get() ;
116 if (histObs.size()!=dvars->size()) {
117 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
118 << ") ERROR histogram variable list and RooDataHist must contain the same variables." << std::endl ;
119 throw(std::string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
120 }
121
122 for (const auto arg : histObs) {
123 if (!dvars->find(arg->GetName())) {
124 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
125 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
126 throw(std::string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
127 }
128 if (!arg->isFundamental()) {
129 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
130 << ") ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory." << std::endl ;
131 throw(std::string("RooHistPdf::ctor() ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory.")) ;
132 }
133 }
134
135
136 // Adjust ranges of _histObsList to those of _dataHist
137 for (const auto hobs : _histObsList) {
138 // Guaranteed to succeed, since checked above in constructor
139 RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
140 RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
141 if (dhreal){
142 (static_cast<RooRealVar*>(hobs))->setRange(dhreal->getMin(),dhreal->getMax()) ;
143 }
144 }
145}
146
147RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgSet &vars, std::unique_ptr<RooDataHist> dhist,
148 int intOrder)
149 : RooHistPdf{name, title, vars, *dhist, intOrder}
150{
151 initializeOwnedDataHist(std::move(dhist));
152}
153RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList &pdfObs, const RooArgList &histObs,
154 std::unique_ptr<RooDataHist> dhist, int intOrder)
155 : RooHistPdf{name, title, pdfObs, histObs, *dhist, intOrder}
156{
157 initializeOwnedDataHist(std::move(dhist));
158}
159
160
161////////////////////////////////////////////////////////////////////////////////
162/// Copy constructor
163
164RooHistPdf::RooHistPdf(const RooHistPdf& other, const char* name) :
165 RooAbsPdf(other,name),
166 _pdfObsList("pdfObs",this,other._pdfObsList),
167 _dataHist(other._dataHist),
168 _codeReg(other._codeReg),
169 _intOrder(other._intOrder),
170 _cdfBoundaries(other._cdfBoundaries),
171 _totVolume(other._totVolume),
172 _unitNorm(other._unitNorm)
173{
175}
176
178 if (_ownedDataHist) return _ownedDataHist.get();
179 _ownedDataHist.reset(static_cast<RooDataHist*>(_dataHist->Clone(newname)));
181 return _dataHist;
182}
183
185{
186 std::span<double> output = ctx.output();
187
188 // For interpolation and histograms of higher dimension, use base function
189 if (_pdfObsList.size() > 1) {
191 return;
192 }
193
194 auto xVals = ctx.at(_pdfObsList[0]);
195 _dataHist->weights(output.data(), xVals, _intOrder, true, _cdfBoundaries);
196}
197
198
199////////////////////////////////////////////////////////////////////////////////
200/// Return the current value: The value of the bin enclosing the current coordinates
201/// of the observables, normalized by the histograms contents. Interpolation
202/// is applied if the RooHistPdf is configured to do that.
203
205{
206 // Transfer values from
207 for (unsigned int i=0; i < _pdfObsList.size(); ++i) {
208 RooAbsArg* harg = _histObsList[i];
209 RooAbsArg* parg = _pdfObsList[i];
210
211 if (harg != parg) {
212 parg->syncCache() ;
213 harg->copyCache(parg,true) ;
214 if (!harg->inRange(nullptr)) {
215 return 0 ;
216 }
217 }
218 }
219
221
222 return std::max(ret, 0.0);
223}
224
225////////////////////////////////////////////////////////////////////////////////
226/// Return the total volume spanned by the observables of the RooHistPdf
227
229{
230 // Return previously calculated value, if any
231 if (_totVolume>0) {
232 return _totVolume ;
233 }
234 _totVolume = 1. ;
235
236 for (const auto arg : _histObsList) {
237 RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
238 if (real) {
239 _totVolume *= (real->getMax()-real->getMin()) ;
240 } else {
241 RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
242 if (cat) {
243 _totVolume *= cat->numTypes() ;
244 }
245 }
246 }
247
248 return _totVolume ;
249}
250
251namespace {
252
253bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
254{
255 const RooAbsRealLValue *_x = dynamic_cast<const RooAbsRealLValue*>(&x);
256 const RooAbsRealLValue *_y = dynamic_cast<const RooAbsRealLValue*>(&y);
257 if (!_x || !_y) return false;
258 if (!range || !strlen(range) || !_x->hasRange(range) ||
259 _x->getBinningPtr(range)->isParameterized()) {
260 // parameterized ranges may be full range now, but that might change,
261 // so return false
262 if (range && strlen(range) && _x->getBinningPtr(range)->isParameterized())
263 return false;
264 return (_x->getMin() == _y->getMin() && _x->getMax() == _y->getMax());
265 }
266 return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
267}
268
269bool okayForAnalytical(RooAbsArg const& obs, RooArgSet const& allVars)
270{
271 auto lobs = dynamic_cast<RooAbsRealLValue const*>(&obs);
272 if(lobs == nullptr) return false;
273
274 bool isOkayForAnalyticalInt = false;
275
276 for(RooAbsArg *var : allVars) {
277 if(obs.dependsOn(*var)) {
278 if(!lobs->isJacobianOK(*var)) return false;
279 isOkayForAnalyticalInt = true;
280 }
281 }
282
283 return isOkayForAnalyticalInt;
284}
285
286} // namespace
287
288
290 RooArgSet& analVars,
291 const char* rangeName,
292 RooArgSet const& histObsList,
293 RooArgSet const& pdfObsList,
294 Int_t intOrder)
295{
296 // First make list of pdf observables to histogram observables
297 // and select only those for which the integral is over the full range
298
299 Int_t code = 0;
300 Int_t frcode = 0;
301 for (unsigned int n=0; n < pdfObsList.size() && n < histObsList.size(); ++n) {
302 const auto pa = pdfObsList[n];
303 const auto ha = histObsList[n];
304
305 if (okayForAnalytical(*pa, allVars)) {
306 code |= 2 << n;
307 analVars.add(*pa);
308 if (fullRange(*pa, *ha, rangeName)) {
309 frcode |= 2 << n;
310 }
311 }
312 }
313
314 if (code == frcode) {
315 // integrate over full range of all observables - use bit 0 to indicate
316 // full range integration over all observables
317 code |= 1;
318 }
319
320 // Disable partial analytical integrals if interpolation is used, and we
321 // integrate over sub-ranges, but leave them enabled when we integrate over
322 // the full range of one or several variables
323 if (intOrder > 1 && !(code & 1)) {
324 analVars.removeAll();
325 return 0;
326 }
327 return (code >= 2) ? code : 0;
328}
329
330
332 const char* rangeName,
333 RooArgSet const& histObsList,
334 RooArgSet const& pdfObsList,
335 RooDataHist& dataHist,
336 bool histFuncMode) {
337 // Simplest scenario, full-range integration over all dependents
338 if (((2 << histObsList.size()) - 1) == code) {
339 return dataHist.sum(histFuncMode);
340 }
341
342 // Partial integration scenario, retrieve set of variables, calculate partial
343 // sum, figure out integration ranges (if needed)
344 RooArgSet intSet;
345 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
346 for (unsigned int n=0; n < pdfObsList.size() && n < histObsList.size(); ++n) {
347 const auto pa = pdfObsList[n];
348 const auto ha = histObsList[n];
349
350 if (code & (2 << n)) {
351 intSet.add(*ha);
352 }
353 if (!(code & 1)) {
354 ranges[ha] = RooHelpers::getRangeOrBinningInterval(pa, rangeName);
355 }
356 // WVE must sync hist slice list values to pdf slice list
357 // Transfer values from
358 if (ha != pa) {
359 pa->syncCache();
360 ha->copyCache(pa,true);
361 }
362 }
363
364 double ret = (code & 1) ? dataHist.sum(intSet,histObsList,true,!histFuncMode) :
365 dataHist.sum(intSet,histObsList,true,!histFuncMode, ranges);
366
367 return ret ;
368}
369
370////////////////////////////////////////////////////////////////////////////////
371/// Determine integration scenario. If no interpolation is used,
372/// RooHistPdf can perform all integrals over its dependents
373/// analytically via partial or complete summation of the input
374/// histogram. If interpolation is used on the integral over
375/// all histogram observables is supported
376
377Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
378{
379 return getAnalyticalIntegral(allVars, analVars, rangeName, _histObsList, _pdfObsList, _intOrder);
380}
381
382
383////////////////////////////////////////////////////////////////////////////////
384/// Return integral identified by 'code'. The actual integration
385/// is deferred to RooDataHist::sum() which implements partial
386/// or complete summation over the histograms contents.
387
388double RooHistPdf::analyticalIntegral(Int_t code, const char* rangeName) const
389{
390 return analyticalIntegral(code, rangeName, _histObsList, _pdfObsList, *_dataHist, false);
391}
392
393
394bool RooHistPdf::forceAnalyticalInt(RooArgSet const& pdfObsList, const RooAbsArg& dep)
395{
396 bool isOkayForAnalyticalInt = false;
397
398 for (RooAbsArg * obs : pdfObsList) {
399 if(obs->dependsOn(dep)) {
400 // If the observable doesn't depend linearly on the integration
401 // variable we will not do analytical integration.
402 auto lvalue = dynamic_cast<RooAbsRealLValue const*>(obs);
403 if(!(lvalue && lvalue->isJacobianOK(dep))) return false;
404 isOkayForAnalyticalInt = true;
405 }
406 }
407
408 return isOkayForAnalyticalInt;
409}
410
411
413{
414 return forceAnalyticalInt(_pdfObsList, dep);
415}
416
417
418////////////////////////////////////////////////////////////////////////////////
419/// Return sampling hint for making curves of (projections) of this function
420/// as the recursive division strategy of RooCurve cannot deal efficiently
421/// with the vertical lines that occur in a non-interpolated histogram
422
423std::list<double>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
424{
426}
427
428
429std::list<double>* RooHistPdf::plotSamplingHint(RooDataHist const& dataHist,
430 RooArgSet const& pdfObsList,
431 RooArgSet const& histObsList,
432 int intOrder,
433 RooAbsRealLValue& obs,
434 double xlo,
435 double xhi)
436{
437 // No hints are required when interpolation is used
438 if (intOrder>0) {
439 return nullptr;
440 }
441
442 // Check that observable is in dataset, if not no hint is generated
443 RooAbsArg* dhObs = nullptr;
444 for (unsigned int i=0; i < pdfObsList.size(); ++i) {
445 RooAbsArg* histObs = histObsList[i];
446 RooAbsArg* pdfObs = pdfObsList[i];
447 if (std::string(obs.GetName())==pdfObs->GetName()) {
448 dhObs = dataHist.get()->find(histObs->GetName()) ;
449 break;
450 }
451 }
452
453 if (!dhObs) {
454 return nullptr;
455 }
456 RooAbsLValue* lval = dynamic_cast<RooAbsLValue*>(dhObs) ;
457 if (!lval) {
458 return nullptr;
459 }
460
461 // Retrieve position of all bin boundaries
462
463 const RooAbsBinning* binning = lval->getBinningPtr(nullptr);
464 std::span<const double> boundaries{binning->array(), static_cast<std::size_t>(binning->numBoundaries())};
465
466 // Use the helper function from RooCurve to make sure to get sampling hints
467 // that work with the RooFitPlotting.
468 return RooCurve::plotSamplingHintForBinBoundaries(boundaries, xlo, xhi);
469}
470
471
472////////////////////////////////////////////////////////////////////////////////
473/// Return sampling hint for making curves of (projections) of this function
474/// as the recursive division strategy of RooCurve cannot deal efficiently
475/// with the vertical lines that occur in a non-interpolated histogram
476
477std::list<double>* RooHistPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
478{
479 // No hints are required when interpolation is used
480 if (_intOrder>0) {
481 return nullptr;
482 }
483
484 // Check that observable is in dataset, if not no hint is generated
485 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
486 if (!lvarg) {
487 return nullptr ;
488 }
489
490 // Retrieve position of all bin boundaries
491 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
492 double* boundaries = binning->array() ;
493
494 auto hint = new std::list<double> ;
495
496 // Construct array with pairs of points positioned epsilon to the left and
497 // right of the bin boundaries
498 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
499 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
500 hint->push_back(boundaries[i]) ;
501 }
502 }
503
504 return hint ;
505}
506
507
508
509
510////////////////////////////////////////////////////////////////////////////////
511/// Only handle case of maximum in all variables
512
514{
515 std::unique_ptr<RooAbsCollection> common{_pdfObsList.selectCommon(vars)};
516 if (common->size()==_pdfObsList.size()) {
517 return 1;
518 }
519 return 0 ;
520}
521
522
523////////////////////////////////////////////////////////////////////////////////
524
525double RooHistPdf::maxVal(Int_t code) const
526{
527 R__ASSERT(code==1) ;
528
529 double max(-1) ;
530 for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
531 _dataHist->get(i) ;
532 double wgt = _dataHist->weight() ;
533 if (wgt>max) max=wgt ;
534 }
535
536 return max*1.05 ;
537}
538
539
540
541
542////////////////////////////////////////////////////////////////////////////////
543
545{
546 if (std::abs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return false ;
547 if (dh1.numEntries() != dh2.numEntries()) return false ;
548 for (int i=0 ; i < dh1.numEntries() ; i++) {
549 dh1.get(i) ;
550 dh2.get(i) ;
551 if (std::abs(dh1.weight()-dh2.weight())>1e-8) return false ;
552 }
553 return true ;
554}
555
556
557
558////////////////////////////////////////////////////////////////////////////////
559/// Check if our datahist is already in the workspace
560
562{
563 for(auto const& data : ws.allData()) {
564 // If your dataset is already in this workspace nothing needs to be done
565 if (data == _dataHist) {
566 return false ;
567 }
568 }
569
570 // Check if dataset with given name already exists
571 if (RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName())) {
572
573 // Yes it exists - now check if it is identical to our internal histogram
574 if (wsdata->InheritsFrom(RooDataHist::Class())) {
575
576 // Check if histograms are identical
577 if (areIdentical(static_cast<RooDataHist&>(*wsdata),*_dataHist)) {
578
579 // Exists and is of correct type, and identical -- adjust internal pointer to WS copy
580 _dataHist = static_cast<RooDataHist*>(wsdata) ;
581 } else {
582
583 // not identical, clone rename and import
584 auto uniqueName = std::string(_dataHist->GetName()) + "_" + GetName();
585 bool flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.c_str()),RooFit::Embedded()) ;
586 if (flag) {
587 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << std::endl ;
588 return true ;
589 }
590 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(uniqueName)) ;
591 }
592
593 } else {
594
595 // Exists and is NOT of correct type: clone rename and import
596 auto uniqueName = std::string(_dataHist->GetName()) + "_" + GetName();
597 bool flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.c_str()),RooFit::Embedded()) ;
598 if (flag) {
599 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << std::endl ;
600 return true ;
601 }
602 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(uniqueName));
603
604 }
605 return false ;
606 }
607
608 // We need to import our datahist into the workspace
610
611 // Redirect our internal pointer to the copy in the workspace
612 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(_dataHist->GetName())) ;
613 return false ;
614}
615
616
617////////////////////////////////////////////////////////////////////////////////
618/// Stream an object of class RooHistPdf.
619
621{
622 if (R__b.IsReading()) {
624 // WVE - interim solution - fix proxies here
625 //_proxyList.Clear() ;
626 //registerProxy(_pdfObsList) ;
627 } else {
629 }
630}
#define e(i)
Definition RSha256.hxx:103
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:382
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
Option_t Option_t TPoint TPoint const char GetTextMagnitude GetFillStyle GetLineColor GetLineWidth GetMarkerStyle GetTextAlign GetTextColor GetTextSize void data
char name[80]
Definition TGX11.cxx:110
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
virtual void copyCache(const RooAbsArg *source, bool valueOnly=false, bool setValDirty=true)=0
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.
virtual void syncCache(const RooArgSet *nset=nullptr)=0
friend void RooRefArray::Streamer(TBuffer &)
virtual bool inRange(const char *) const
Definition RooAbsArg.h:319
Abstract base class for RooRealVar binning definitions.
virtual bool isParameterized() const
Interface function.
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)
virtual void removeAll()
Remove all arguments from our set, deleting them if we own them.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
Storage_t::size_type size() const
virtual RooAbsArg * addClone(const RooAbsArg &var, bool silent=false)
Add a clone of the specified argument to list.
RooAbsArg * find(const char *name) const
Find object with given name in list.
Abstract base class for binned and unbinned datasets.
Definition RooAbsData.h:57
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
Abstract base class for objects that are lvalues, i.e.
virtual const RooAbsBinning * getBinningPtr(const char *rangeName) const =0
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
virtual double getMax(const char *name=nullptr) const
Get maximum of currently defined range.
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
const RooAbsBinning * getBinningPtr(const char *rangeName) const override
bool hasRange(const char *name) const override
Check if variable has a binning with given name.
virtual void doEval(RooFit::EvalContext &) const
Base function for computing multiple values of a RooAbsReal.
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
RooArgSet * selectCommon(const RooAbsCollection &refColl) const
Use RooAbsCollection::selecCommon(), but return as RooArgSet.
Definition RooArgSet.h:154
Object to represent discrete states.
Definition RooCategory.h:28
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
static std::list< double > * plotSamplingHintForBinBoundaries(std::span< const double > boundaries, double xlo, double xhi)
Returns sampling hints for a histogram with given boundaries.
Definition RooCurve.cxx:893
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
double sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
void weights(double *output, std::span< double const > xVals, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A vectorized version of RooDataHist::weight() for one dimensional histograms with up to one dimension...
static TClass * Class()
TObject * Clone(const char *newname="") const override
Make a clone of an object using the Streamer facility.
Definition RooDataHist.h:61
double weight(std::size_t i) const
Return weight of i-th bin.
double weightFast(const RooArgSet &bin, int intOrder, bool correctForBinSize, bool cdfBoundaries)
A faster version of RooDataHist::weight that assumes the passed arguments are aligned with the histog...
const RooArgSet * get() const override
Get bin centre of current bin.
Definition RooDataHist.h:82
double sumEntries() const override
Sum the weights of all bins.
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
A propability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
RooArgSet _histObsList
List of observables defining dimensions of histogram.
Definition RooHistPdf.h:110
Int_t _intOrder
Interpolation order.
Definition RooHistPdf.h:115
bool forceAnalyticalInt(const RooAbsArg &dep) const override
bool areIdentical(const RooDataHist &dh1, const RooDataHist &dh2)
RooDataHist * _dataHist
Unowned pointer to underlying histogram.
Definition RooHistPdf.h:112
bool _cdfBoundaries
Use boundary conditions for CDFs.
Definition RooHistPdf.h:116
std::list< double > * binBoundaries(RooAbsRealLValue &, double, double) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
double totVolume() const
Return the total volume spanned by the observables of the RooHistPdf.
void initializeOwnedDataHist(std::unique_ptr< RooDataHist > &&dataHist)
Definition RooHistPdf.h:148
bool importWorkspaceHook(RooWorkspace &ws) override
Check if our datahist is already in the workspace.
static TClass * Class()
RooSetProxy _pdfObsList
List of observables mapped onto histogram observables.
Definition RooHistPdf.h:111
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Return integral identified by 'code'.
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Return sampling hint for making curves of (projections) of this function as the recursive division st...
RooDataHist & dataHist()
Definition RooHistPdf.h:42
Int_t getMaxVal(const RooArgSet &vars) const override
Only handle case of maximum in all variables.
double _totVolume
! Total volume of space (product of ranges of observables)
Definition RooHistPdf.h:117
RooDataHist * cloneAndOwnDataHist(const char *newname="")
Replaces underlying RooDataHist with a clone, which is now owned, and returns the clone.
std::unique_ptr< RooDataHist > _ownedDataHist
! Owned pointer to underlying histogram
Definition RooHistPdf.h:113
void doEval(RooFit::EvalContext &) const override
Base function for computing multiple values of a RooAbsReal.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Determine integration scenario.
double evaluate() const override
Return the current value: The value of the bin enclosing the current coordinates of the observables,...
bool _unitNorm
Assume contents is unit normalized (for use as pdf cache)
Definition RooHistPdf.h:118
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setRange(const char *name, double min, double max)
Set a fit or plotting range.
Persistable container for RooFit projects.
RooAbsData * embeddedData(RooStringView name) const
Retrieve dataset (binned or unbinned) with given name. A null pointer is returned if not found.
std::list< RooAbsData * > allData() const
Return list of all dataset in the workspace.
bool import(const RooAbsArg &arg, const RooCmdArg &arg1={}, const RooCmdArg &arg2={}, const RooCmdArg &arg3={}, const RooCmdArg &arg4={}, const RooCmdArg &arg5={}, const RooCmdArg &arg6={}, const RooCmdArg &arg7={}, const RooCmdArg &arg8={}, const RooCmdArg &arg9={})
Import a RooAbsArg object, e.g.
Buffer base class used for serializing objects.
Definition TBuffer.h:43
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
RooCmdArg Rename(const char *suffix)
RooCmdArg Embedded(bool flag=true)
Double_t y[n]
Definition legend1.C:17
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
std::pair< double, double > getRangeOrBinningInterval(RooAbsArg const *arg, const char *rangeName)
static void output()