Logo ROOT  
Reference Guide
RooHistPdf.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitCore *
4 * File: $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#ifndef ROO_HIST_PDF
17#define ROO_HIST_PDF
18
19#include "RooAbsPdf.h"
20#include "RooRealProxy.h"
21#include "RooSetProxy.h"
22#include "RooAICRegistry.h"
23#include <list>
24
25class RooRealVar;
26class RooAbsReal;
27class RooDataHist ;
28
29class RooHistPdf : public RooAbsPdf {
30public:
31 RooHistPdf() ;
32 RooHistPdf(const char *name, const char *title, const RooArgSet& vars, const RooDataHist& dhist, Int_t intOrder=0);
33 RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs, const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder=0);
34 RooHistPdf(const RooHistPdf& other, const char* name=0);
35 TObject* clone(const char* newname) const override { return new RooHistPdf(*this,newname); }
36 ~RooHistPdf() override ;
37
39 // Return RooDataHist that is represented
40 return *_dataHist ;
41 }
42 const RooDataHist& dataHist() const {
43 // Return RooDataHist that is represented
44 return *_dataHist ;
45 }
46
48 // Set histogram interpolation order
49 _intOrder = order ;
50 }
52 // Return histogram interpolation order
53 return _intOrder ;
54 }
55
57 RooArgSet& analVars,
58 const char* rangeName,
59 RooArgSet const& histObsList,
60 RooSetProxy const& pdfObsList,
61 Int_t intOrder) ;
62
63 static double analyticalIntegral(Int_t code,
64 const char* rangeName,
65 RooArgSet const& histObsList,
66 RooSetProxy const& pdfObsList,
68 bool histFuncMode) ;
69
70 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=0) const override ;
71 double analyticalIntegral(Int_t code, const char* rangeName=0) const override ;
72
73 void setCdfBoundaries(bool flag) {
74 // Set use of special boundary conditions for c.d.f.s
75 _cdfBoundaries = flag ;
76 }
77 bool getCdfBoundaries() const {
78 // If true, special boundary conditions for c.d.f.s are used
79 return _cdfBoundaries ;
80 }
81
82 void setUnitNorm(bool flag) {
83 // Declare contents to have unit normalization
84 _unitNorm = flag ;
85 }
86 bool haveUnitNorm() const {
87 // Return true if contents is declared to be unit normalized
88 return _unitNorm ;
89 }
90
91 bool selfNormalized() const override { return _unitNorm ; }
92
93 Int_t getMaxVal(const RooArgSet& vars) const override ;
94 double maxVal(Int_t code) const override ;
95
96 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override ;
97 std::list<double>* binBoundaries(RooAbsRealLValue& /*obs*/, double /*xlo*/, double /*xhi*/) const override ;
98 bool isBinnedDistribution(const RooArgSet&) const override { return _intOrder==0 ; }
99
100
101protected:
102
103 bool areIdentical(const RooDataHist& dh1, const RooDataHist& dh2) ;
104
105 bool importWorkspaceHook(RooWorkspace& ws) override ;
106
107 double evaluate() const override;
108 double totalVolume() const ;
109 friend class RooAbsCachedPdf ;
110 double totVolume() const ;
111
112 RooArgSet _histObsList ; ///< List of observables defining dimensions of histogram
113 RooSetProxy _pdfObsList ; ///< List of observables mapped onto histogram observables
114 RooDataHist* _dataHist ; ///< Unowned pointer to underlying histogram
115 mutable RooAICRegistry _codeReg ; ///<! Auxiliary class keeping tracking of analytical integration code
116 Int_t _intOrder ; ///< Interpolation order
117 bool _cdfBoundaries ; ///< Use boundary conditions for CDFs.
118 mutable double _totVolume ; ///<! Total volume of space (product of ranges of observables)
119 bool _unitNorm ; ///< Assume contents is unit normalized (for use as pdf cache)
120
121 ClassDefOverride(RooHistPdf,4) // Histogram based PDF
122};
123
124#endif
#define ClassDefOverride(name, id)
Definition: Rtypes.h:339
char name[80]
Definition: TGX11.cxx:110
RooAICRegistry is a utility class for operator p.d.f classes that keeps track of analytical integrati...
RooAbsCachedPdf is the abstract base class for p.d.f.s that need or want to cache their evaluate() ou...
RooAbsRealLValue is the common abstract base class for objects that represent a real value that may a...
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:64
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:57
The RooDataHist is a container class to hold N-dimensional binned data.
Definition: RooDataHist.h:45
RooHistPdf implements a probablity density function sampled from a multidimensional histogram.
Definition: RooHistPdf.h:29
Int_t getInterpolationOrder() const
Definition: RooHistPdf.h:51
void setUnitNorm(bool flag)
Definition: RooHistPdf.h:82
RooArgSet _histObsList
List of observables defining dimensions of histogram.
Definition: RooHistPdf.h:112
Int_t _intOrder
Interpolation order.
Definition: RooHistPdf.h:116
RooHistPdf()
Default constructor coverity[UNINIT_CTOR].
Definition: RooHistPdf.cxx:53
bool areIdentical(const RooDataHist &dh1, const RooDataHist &dh2)
Definition: RooHistPdf.cxx:506
RooDataHist * _dataHist
Unowned pointer to underlying histogram.
Definition: RooHistPdf.h:114
bool _cdfBoundaries
Use boundary conditions for CDFs.
Definition: RooHistPdf.h:117
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...
Definition: RooHistPdf.cxx:437
double totVolume() const
Return the total volume spanned by the observables of the RooHistPdf.
Definition: RooHistPdf.cxx:228
bool importWorkspaceHook(RooWorkspace &ws) override
Check if our datahist is already in the workspace.
Definition: RooHistPdf.cxx:523
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...
Definition: RooHistPdf.cxx:379
RooSetProxy _pdfObsList
List of observables mapped onto histogram observables.
Definition: RooHistPdf.h:113
static double analyticalIntegral(Int_t code, const char *rangeName, RooArgSet const &histObsList, RooSetProxy const &pdfObsList, RooDataHist &dataHist, bool histFuncMode)
Definition: RooHistPdf.cxx:311
bool isBinnedDistribution(const RooArgSet &) const override
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition: RooHistPdf.h:98
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
Definition: RooHistPdf.cxx:487
RooAICRegistry _codeReg
! Auxiliary class keeping tracking of analytical integration code
Definition: RooHistPdf.h:115
TObject * clone(const char *newname) const override
Definition: RooHistPdf.h:35
void setCdfBoundaries(bool flag)
Definition: RooHistPdf.h:73
RooDataHist & dataHist()
Definition: RooHistPdf.h:38
bool haveUnitNorm() const
Definition: RooHistPdf.h:86
Int_t getMaxVal(const RooArgSet &vars) const override
Only handle case of maximum in all variables.
Definition: RooHistPdf.cxx:473
double _totVolume
! Total volume of space (product of ranges of observables)
Definition: RooHistPdf.h:118
double totalVolume() const
static Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName, RooArgSet const &histObsList, RooSetProxy const &pdfObsList, Int_t intOrder)
Definition: RooHistPdf.cxx:270
bool getCdfBoundaries() const
Definition: RooHistPdf.h:77
bool selfNormalized() const override
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
Definition: RooHistPdf.h:91
const RooDataHist & dataHist() const
Definition: RooHistPdf.h:42
~RooHistPdf() override
Destructor.
Definition: RooHistPdf.cxx:189
double evaluate() const override
Return the current value: The value of the bin enclosing the current coordinates of the observables,...
Definition: RooHistPdf.cxx:203
void setInterpolationOrder(Int_t order)
Definition: RooHistPdf.h:47
bool _unitNorm
Assume contents is unit normalized (for use as pdf cache)
Definition: RooHistPdf.h:119
RooRealVar represents a variable that can be changed from the outside.
Definition: RooRealVar.h:40
The RooWorkspace is a persistable container for RooFit projects.
Definition: RooWorkspace.h:43
Mother of all ROOT objects.
Definition: TObject.h:37
void ws()
Definition: ws.C:66