Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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
22RooHistPdf implements a 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((RooDataHist*)&dhist),
59 _codeReg(10),
60 _intOrder(intOrder),
61 _cdfBoundaries(false)
62{
64 _pdfObsList.add(vars) ;
65
66 // Verify that vars and dhist.get() have identical contents
67 const RooArgSet* dvars = dhist.get() ;
68 if (vars.size()!=dvars->size()) {
69 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
70 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
71 assert(0) ;
72 }
73 for (const auto arg : vars) {
74 if (!dvars->find(arg->GetName())) {
75 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
76 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
77 assert(0) ;
78 }
79 }
80
81
82 // Adjust ranges of _histObsList to those of _dataHist
83 for (const auto hobs : _histObsList) {
84 // Guaranteed to succeed, since checked above in ctor
85 RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
86 RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
87 if (dhreal){
88 ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
89 }
90 }
91
92}
93
94
95
96
97////////////////////////////////////////////////////////////////////////////////
98/// Constructor from a RooDataHist. The first list of observables are the p.d.f.
99/// observables, which may any RooAbsReal (function or variable). The second list
100/// are the corresponding observables in the RooDataHist which must be of type
101/// RooRealVar or RooCategory This constructor thus allows to apply a coordinate transformation
102/// on the histogram data to be applied.
103
104RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs,
105 const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder) :
106 RooAbsPdf(name,title),
107 _pdfObsList("pdfObs","List of p.d.f. observables",this),
108 _dataHist((RooDataHist*)&dhist),
109 _codeReg(10),
110 _intOrder(intOrder),
111 _cdfBoundaries(false)
112{
113 _histObsList.addClone(histObs) ;
114 _pdfObsList.add(pdfObs) ;
115
116 // Verify that vars and dhist.get() have identical contents
117 const RooArgSet* dvars = dhist.get() ;
118 if (histObs.size()!=dvars->size()) {
119 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
120 << ") ERROR histogram variable list and RooDataHist must contain the same variables." << std::endl ;
121 throw(std::string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
122 }
123
124 for (const auto arg : histObs) {
125 if (!dvars->find(arg->GetName())) {
126 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
127 << ") ERROR variable list and RooDataHist must contain the same variables." << std::endl ;
128 throw(std::string("RooHistPdf::ctor() ERROR: histogram variable list and RooDataHist must contain the same variables")) ;
129 }
130 if (!arg->isFundamental()) {
131 coutE(InputArguments) << "RooHistPdf::ctor(" << GetName()
132 << ") ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory." << std::endl ;
133 throw(std::string("RooHistPdf::ctor() ERROR all elements of histogram observables set must be of type RooRealVar or RooCategory.")) ;
134 }
135 }
136
137
138 // Adjust ranges of _histObsList to those of _dataHist
139 for (const auto hobs : _histObsList) {
140 // Guaranteed to succeed, since checked above in ctor
141 RooAbsArg* dhobs = dhist.get()->find(hobs->GetName()) ;
142 RooRealVar* dhreal = dynamic_cast<RooRealVar*>(dhobs) ;
143 if (dhreal){
144 ((RooRealVar*)hobs)->setRange(dhreal->getMin(),dhreal->getMax()) ;
145 }
146 }
147}
148
149
150RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgSet& vars,
151 std::unique_ptr<RooDataHist> dhist, int intOrder)
152 : RooHistPdf{name, title, vars, *dhist, intOrder}
153{
154 _ownedDataHist = std::move(dhist);
155}
156RooHistPdf::RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs, const RooArgList& histObs,
157 std::unique_ptr<RooDataHist> dhist, int intOrder)
158 : RooHistPdf{name, title, pdfObs, histObs, *dhist, intOrder}
159{
160 _ownedDataHist = std::move(dhist);
161}
162
163
164////////////////////////////////////////////////////////////////////////////////
165/// Copy constructor
166
167RooHistPdf::RooHistPdf(const RooHistPdf& other, const char* name) :
168 RooAbsPdf(other,name),
169 _pdfObsList("pdfObs",this,other._pdfObsList),
170 _dataHist(other._dataHist),
171 _codeReg(other._codeReg),
172 _intOrder(other._intOrder),
173 _cdfBoundaries(other._cdfBoundaries),
174 _totVolume(other._totVolume),
175 _unitNorm(other._unitNorm)
176{
178
179}
180
181
182
183
184////////////////////////////////////////////////////////////////////////////////
185/// Destructor
186
188{
189
190}
191
193 if (_ownedDataHist) return _ownedDataHist.get();
194 _ownedDataHist.reset(static_cast<RooDataHist*>(_dataHist->Clone(newname)));
196 return _dataHist;
197}
198
199void RooHistPdf::computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const {
200
201 // For interpolation and histograms of higher dimension, use base function
202 if(_pdfObsList.size() > 1) {
203 RooAbsReal::computeBatch(output, nEvents, dataMap);
204 return;
205 }
206
207 auto xVals = dataMap.at(_pdfObsList[0]);
209}
210
211
212////////////////////////////////////////////////////////////////////////////////
213/// Return the current value: The value of the bin enclosing the current coordinates
214/// of the observables, normalized by the histograms contents. Interpolation
215/// is applied if the RooHistPdf is configured to do that.
216
218{
219 // Transfer values from
220 for (unsigned int i=0; i < _pdfObsList.size(); ++i) {
221 RooAbsArg* harg = _histObsList[i];
222 RooAbsArg* parg = _pdfObsList[i];
223
224 if (harg != parg) {
225 parg->syncCache() ;
226 harg->copyCache(parg,true) ;
227 if (!harg->inRange(nullptr)) {
228 return 0 ;
229 }
230 }
231 }
232
234
235 return std::max(ret, 0.0);
236}
237
239 RooDataHist const *dataHist, const RooArgSet &obs, bool correctForBinSize)
240{
241 if (intOrder != 0) {
242 ooccoutE(klass, InputArguments) << "RooHistPdf::weight(" << klass->GetName()
243 << ") ERROR: Code Squashing currently only supports non-interpolation cases."
244 << std::endl;
245 return;
246 }
247
248 std::string const &idxName = dataHist->calculateTreeIndexForCodeSquash(klass, ctx, obs);
249 std::string const &weightName = dataHist->declWeightArrayForCodeSquash(klass, ctx, correctForBinSize);
250 ctx.addResult(klass, weightName + "[" + idxName + "]");
251}
252
254{
256}
257
258////////////////////////////////////////////////////////////////////////////////
259/// Return the total volume spanned by the observables of the RooHistPdf
260
262{
263 // Return previously calculated value, if any
264 if (_totVolume>0) {
265 return _totVolume ;
266 }
267 _totVolume = 1. ;
268
269 for (const auto arg : _histObsList) {
270 RooRealVar* real = dynamic_cast<RooRealVar*>(arg) ;
271 if (real) {
272 _totVolume *= (real->getMax()-real->getMin()) ;
273 } else {
274 RooCategory* cat = dynamic_cast<RooCategory*>(arg) ;
275 if (cat) {
276 _totVolume *= cat->numTypes() ;
277 }
278 }
279 }
280
281 return _totVolume ;
282}
283
284namespace {
285
286bool fullRange(const RooAbsArg& x, const RooAbsArg& y ,const char* range)
287{
288 const RooAbsRealLValue *_x = dynamic_cast<const RooAbsRealLValue*>(&x);
289 const RooAbsRealLValue *_y = dynamic_cast<const RooAbsRealLValue*>(&y);
290 if (!_x || !_y) return false;
291 if (!range || !strlen(range) || !_x->hasRange(range) ||
292 _x->getBinningPtr(range)->isParameterized()) {
293 // parameterized ranges may be full range now, but that might change,
294 // so return false
295 if (range && strlen(range) && _x->getBinningPtr(range)->isParameterized())
296 return false;
297 return (_x->getMin() == _y->getMin() && _x->getMax() == _y->getMax());
298 }
299 return (_x->getMin(range) == _y->getMin() && _x->getMax(range) == _y->getMax());
300}
301
302bool okayForAnalytical(RooAbsArg const& obs, RooArgSet const& allVars)
303{
304 auto lobs = dynamic_cast<RooAbsRealLValue const*>(&obs);
305 if(lobs == nullptr) return false;
306
307 bool isOkayForAnalyticalInt = false;
308
309 for(RooAbsArg *var : allVars) {
310 if(obs.dependsOn(*var)) {
311 if(!lobs->isJacobianOK(*var)) return false;
312 isOkayForAnalyticalInt = true;
313 }
314 }
315
316 return isOkayForAnalyticalInt;
317}
318
319} // namespace
320
321
323 RooArgSet& analVars,
324 const char* rangeName,
325 RooArgSet const& histObsList,
326 RooArgSet const& pdfObsList,
327 Int_t intOrder)
328{
329 // First make list of pdf observables to histogram observables
330 // and select only those for which the integral is over the full range
331
332 Int_t code = 0;
333 Int_t frcode = 0;
334 for (unsigned int n=0; n < pdfObsList.size() && n < histObsList.size(); ++n) {
335 const auto pa = pdfObsList[n];
336 const auto ha = histObsList[n];
337
338 if (okayForAnalytical(*pa, allVars)) {
339 code |= 2 << n;
340 analVars.add(*pa);
341 if (fullRange(*pa, *ha, rangeName)) {
342 frcode |= 2 << n;
343 }
344 }
345 }
346
347 if (code == frcode) {
348 // integrate over full range of all observables - use bit 0 to indicate
349 // full range integration over all observables
350 code |= 1;
351 }
352
353 // Disable partial analytical integrals if interpolation is used, and we
354 // integrate over sub-ranges, but leave them enabled when we integrate over
355 // the full range of one or several variables
356 if (intOrder > 1 && !(code & 1)) {
357 analVars.removeAll();
358 return 0;
359 }
360 return (code >= 2) ? code : 0;
361}
362
363
365 const char* rangeName,
366 RooArgSet const& histObsList,
367 RooArgSet const& pdfObsList,
368 RooDataHist& dataHist,
369 bool histFuncMode) {
370 // Simplest scenario, full-range integration over all dependents
371 if (((2 << histObsList.size()) - 1) == code) {
372 return dataHist.sum(histFuncMode);
373 }
374
375 // Partial integration scenario, retrieve set of variables, calculate partial
376 // sum, figure out integration ranges (if needed)
377 RooArgSet intSet;
378 std::map<const RooAbsArg*, std::pair<double, double> > ranges;
379 for (unsigned int n=0; n < pdfObsList.size() && n < histObsList.size(); ++n) {
380 const auto pa = pdfObsList[n];
381 const auto ha = histObsList[n];
382
383 if (code & (2 << n)) {
384 intSet.add(*ha);
385 }
386 if (!(code & 1)) {
387 ranges[ha] = RooHelpers::getRangeOrBinningInterval(pa, rangeName);
388 }
389 // WVE must sync hist slice list values to pdf slice list
390 // Transfer values from
391 if (ha != pa) {
392 pa->syncCache();
393 ha->copyCache(pa,true);
394 }
395 }
396
397 double ret = (code & 1) ? dataHist.sum(intSet,histObsList,true,!histFuncMode) :
398 dataHist.sum(intSet,histObsList,true,!histFuncMode, ranges);
399
400 return ret ;
401}
402
403std::string RooHistPdf::rooHistIntegralTranslateImpl(int code, RooAbsArg const *klass, RooDataHist const *dataHist,
404 const RooArgSet &obs, bool histFuncMode)
405{
406 if (((2 << obs.size()) - 1) != code) {
407 oocoutE(klass, InputArguments)
408 << "RooHistPdf::integral(" << klass->GetName()
409 << ") ERROR: AD currently only supports integrating over all histogram observables." << std::endl;
410 return "";
411 }
412 return std::to_string(dataHist->sum(histFuncMode));
413}
414
415std::string RooHistPdf::buildCallToAnalyticIntegral(int code, const char * /*rangeName */,
416 RooFit::Detail::CodeSquashContext & /* ctx */) const
417{
418 return rooHistIntegralTranslateImpl(code, this, _dataHist, _pdfObsList, false);
419}
420
421////////////////////////////////////////////////////////////////////////////////
422/// Determine integration scenario. If no interpolation is used,
423/// RooHistPdf can perform all integrals over its dependents
424/// analytically via partial or complete summation of the input
425/// histogram. If interpolation is used on the integral over
426/// all histogram observables is supported
427
428Int_t RooHistPdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
429{
430 return getAnalyticalIntegral(allVars, analVars, rangeName, _histObsList, _pdfObsList, _intOrder);
431}
432
433
434////////////////////////////////////////////////////////////////////////////////
435/// Return integral identified by 'code'. The actual integration
436/// is deferred to RooDataHist::sum() which implements partial
437/// or complete summation over the histograms contents.
438
439double RooHistPdf::analyticalIntegral(Int_t code, const char* rangeName) const
440{
441 return analyticalIntegral(code, rangeName, _histObsList, _pdfObsList, *_dataHist, false);
442}
443
444
445bool RooHistPdf::forceAnalyticalInt(RooArgSet const& pdfObsList, const RooAbsArg& dep)
446{
447 bool isOkayForAnalyticalInt = false;
448
449 for (RooAbsArg * obs : pdfObsList) {
450 if(obs->dependsOn(dep)) {
451 // If the observable doesn't depend linearly on the integration
452 // variable we will not do analytical integration.
453 auto lvalue = dynamic_cast<RooAbsRealLValue const*>(obs);
454 if(!(lvalue && lvalue->isJacobianOK(dep))) return false;
455 isOkayForAnalyticalInt = true;
456 }
457 }
458
459 return isOkayForAnalyticalInt;
460}
461
462
464{
465 return forceAnalyticalInt(_pdfObsList, dep);
466}
467
468
469////////////////////////////////////////////////////////////////////////////////
470/// Return sampling hint for making curves of (projections) of this function
471/// as the recursive division strategy of RooCurve cannot deal efficiently
472/// with the vertical lines that occur in a non-interpolated histogram
473
474std::list<double>* RooHistPdf::plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const
475{
477}
478
479
480std::list<double>* RooHistPdf::plotSamplingHint(RooDataHist const& dataHist,
481 RooArgSet const& pdfObsList,
482 RooArgSet const& histObsList,
483 int intOrder,
484 RooAbsRealLValue& obs,
485 double xlo,
486 double xhi)
487{
488 // No hints are required when interpolation is used
489 if (intOrder>0) {
490 return nullptr;
491 }
492
493 // Check that observable is in dataset, if not no hint is generated
494 RooAbsArg* dhObs = nullptr;
495 for (unsigned int i=0; i < pdfObsList.size(); ++i) {
496 RooAbsArg* histObs = histObsList[i];
497 RooAbsArg* pdfObs = pdfObsList[i];
498 if (std::string(obs.GetName())==pdfObs->GetName()) {
499 dhObs = dataHist.get()->find(histObs->GetName()) ;
500 break;
501 }
502 }
503
504 if (!dhObs) {
505 return nullptr;
506 }
507 RooAbsLValue* lval = dynamic_cast<RooAbsLValue*>(dhObs) ;
508 if (!lval) {
509 return nullptr;
510 }
511
512 // Retrieve position of all bin boundaries
513
514 const RooAbsBinning* binning = lval->getBinningPtr(nullptr);
515 std::span<const double> boundaries{binning->array(), static_cast<std::size_t>(binning->numBoundaries())};
516
517 // Use the helper function from RooCurve to make sure to get sampling hints
518 // that work with the RooFitPlotting.
519 return RooCurve::plotSamplingHintForBinBoundaries(boundaries, xlo, xhi);
520}
521
522
523////////////////////////////////////////////////////////////////////////////////
524/// Return sampling hint for making curves of (projections) of this function
525/// as the recursive division strategy of RooCurve cannot deal efficiently
526/// with the vertical lines that occur in a non-interpolated histogram
527
528std::list<double>* RooHistPdf::binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const
529{
530 // No hints are required when interpolation is used
531 if (_intOrder>0) {
532 return nullptr;
533 }
534
535 // Check that observable is in dataset, if not no hint is generated
536 RooAbsLValue* lvarg = dynamic_cast<RooAbsLValue*>(_dataHist->get()->find(obs.GetName())) ;
537 if (!lvarg) {
538 return nullptr ;
539 }
540
541 // Retrieve position of all bin boundaries
542 const RooAbsBinning* binning = lvarg->getBinningPtr(nullptr);
543 double* boundaries = binning->array() ;
544
545 auto hint = new std::list<double> ;
546
547 // Construct array with pairs of points positioned epsilon to the left and
548 // right of the bin boundaries
549 for (Int_t i=0 ; i<binning->numBoundaries() ; i++) {
550 if (boundaries[i]>=xlo && boundaries[i]<=xhi) {
551 hint->push_back(boundaries[i]) ;
552 }
553 }
554
555 return hint ;
556}
557
558
559
560
561////////////////////////////////////////////////////////////////////////////////
562/// Only handle case of maximum in all variables
563
565{
566 std::unique_ptr<RooAbsCollection> common{_pdfObsList.selectCommon(vars)};
567 if (common->size()==_pdfObsList.size()) {
568 return 1;
569 }
570 return 0 ;
571}
572
573
574////////////////////////////////////////////////////////////////////////////////
575
576double RooHistPdf::maxVal(Int_t code) const
577{
578 R__ASSERT(code==1) ;
579
580 double max(-1) ;
581 for (Int_t i=0 ; i<_dataHist->numEntries() ; i++) {
582 _dataHist->get(i) ;
583 double wgt = _dataHist->weight() ;
584 if (wgt>max) max=wgt ;
585 }
586
587 return max*1.05 ;
588}
589
590
591
592
593////////////////////////////////////////////////////////////////////////////////
594
596{
597 if (std::abs(dh1.sumEntries()-dh2.sumEntries())>1e-8) return false ;
598 if (dh1.numEntries() != dh2.numEntries()) return false ;
599 for (int i=0 ; i < dh1.numEntries() ; i++) {
600 dh1.get(i) ;
601 dh2.get(i) ;
602 if (std::abs(dh1.weight()-dh2.weight())>1e-8) return false ;
603 }
604 return true ;
605}
606
607
608
609////////////////////////////////////////////////////////////////////////////////
610/// Check if our datahist is already in the workspace
611
613{
614 for(auto const& data : ws.allData()) {
615 // If your dataset is already in this workspace nothing needs to be done
616 if (data == _dataHist) {
617 return false ;
618 }
619 }
620
621 // Check if dataset with given name already exists
622 if (RooAbsData* wsdata = ws.embeddedData(_dataHist->GetName())) {
623
624 // Yes it exists - now check if it is identical to our internal histogram
625 if (wsdata->InheritsFrom(RooDataHist::Class())) {
626
627 // Check if histograms are identical
628 if (areIdentical((RooDataHist&)*wsdata,*_dataHist)) {
629
630 // Exists and is of correct type, and identical -- adjust internal pointer to WS copy
631 _dataHist = (RooDataHist*) wsdata ;
632 } else {
633
634 // not identical, clone rename and import
635 auto uniqueName = std::string(_dataHist->GetName()) + "_" + GetName();
636 bool flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.c_str()),RooFit::Embedded()) ;
637 if (flag) {
638 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << std::endl ;
639 return true ;
640 }
641 _dataHist = (RooDataHist*) ws.embeddedData(uniqueName) ;
642 }
643
644 } else {
645
646 // Exists and is NOT of correct type: clone rename and import
647 auto uniqueName = std::string(_dataHist->GetName()) + "_" + GetName();
648 bool flag = ws.import(*_dataHist,RooFit::Rename(uniqueName.c_str()),RooFit::Embedded()) ;
649 if (flag) {
650 coutE(ObjectHandling) << " RooHistPdf::importWorkspaceHook(" << GetName() << ") unable to import clone of underlying RooDataHist with unique name " << uniqueName << ", abort" << std::endl ;
651 return true ;
652 }
653 _dataHist = static_cast<RooDataHist*>(ws.embeddedData(uniqueName));
654
655 }
656 return false ;
657 }
658
659 // We need to import our datahist into the workspace
661
662 // Redirect our internal pointer to the copy in the workspace
664 return false ;
665}
666
667
668////////////////////////////////////////////////////////////////////////////////
669/// Stream an object of class RooHistPdf.
670
672{
673 if (R__b.IsReading()) {
675 // WVE - interim solution - fix proxies here
676 //_proxyList.Clear() ;
677 //registerProxy(_pdfObsList) ;
678 } else {
680 }
681}
#define e(i)
Definition RSha256.hxx:103
#define oocoutE(o, a)
#define coutE(a)
#define ooccoutE(o, a)
#define ClassImp(name)
Definition Rtypes.h:377
#define R__ASSERT(e)
Definition TError.h:118
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:79
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:378
RooAbsBinning is the 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.
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.
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
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 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 computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) 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: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...
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:897
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
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:55
std::string calculateTreeIndexForCodeSquash(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
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...
std::string declWeightArrayForCodeSquash(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, bool correctForBinSize) const
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.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
Definition DataMap.cxx:22
RooHistPdf implements a propability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
static void rooHistTranslateImpl(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, int intOrder, RooDataHist const *dataHist, const RooArgSet &obs, bool correctForBinSize)
RooArgSet _histObsList
List of observables defining dimensions of histogram.
Definition RooHistPdf.h:113
Int_t _intOrder
Interpolation order.
Definition RooHistPdf.h:118
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:115
bool _cdfBoundaries
Use boundary conditions for CDFs.
Definition RooHistPdf.h:119
static std::string rooHistIntegralTranslateImpl(int code, RooAbsArg const *klass, RooDataHist const *dataHist, const RooArgSet &obs, bool histFuncMode)
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.
bool importWorkspaceHook(RooWorkspace &ws) override
Check if our datahist is already in the workspace.
static TClass * Class()
std::string buildCallToAnalyticIntegral(int code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
RooSetProxy _pdfObsList
List of observables mapped onto histogram observables.
Definition RooHistPdf.h:114
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
void computeBatch(double *output, size_t size, RooFit::Detail::DataMap const &) const override
Base function for computing multiple values of a RooAbsReal.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Return integral identified by 'code'.
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 ...
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:43
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:120
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:116
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Determine integration scenario.
~RooHistPdf() override
Destructor.
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:121
RooRealVar represents a 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()