Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
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 "RooDataHist.h"
24
25#include <list>
26
27class RooRealVar;
28class RooAbsReal;
29
30class RooHistPdf : public RooAbsPdf {
31public:
33 RooHistPdf(const char *name, const char *title, const RooArgSet& vars, const RooDataHist& dhist, Int_t intOrder=0);
34 RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs, const RooArgList& histObs, const RooDataHist& dhist, Int_t intOrder=0);
35 RooHistPdf(const char *name, const char *title, const RooArgSet& vars,
36 std::unique_ptr<RooDataHist> dhist, int intOrder=0);
37 RooHistPdf(const char *name, const char *title, const RooArgList& pdfObs, const RooArgList& histObs,
38 std::unique_ptr<RooDataHist> dhist, int intOrder=0);
39 RooHistPdf(const RooHistPdf& other, const char* name=nullptr);
40 TObject* clone(const char* newname) const override { return new RooHistPdf(*this,newname); }
41 ~RooHistPdf() override ;
42
44 // Return RooDataHist that is represented
45 return *_dataHist ;
46 }
47 const RooDataHist& dataHist() const {
48 // Return RooDataHist that is represented
49 return *_dataHist ;
50 }
51
52 /// Replaces underlying RooDataHist with a clone, which is now owned, and returns the clone.
53 /// If the underlying RooDataHist is already owned, then that is returned instead of being cloned.
54 RooDataHist* cloneAndOwnDataHist(const char* newname="");
55
57 // Set histogram interpolation order
58 _intOrder = order ;
59 }
61 // Return histogram interpolation order
62 return _intOrder ;
63 }
64
65 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
66 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
67
68 bool forceAnalyticalInt(const RooAbsArg& dep) const override;
69
70 void setCdfBoundaries(bool flag) {
71 // Set use of special boundary conditions for c.d.f.s
72 _cdfBoundaries = flag ;
73 }
74 bool getCdfBoundaries() const {
75 // If true, special boundary conditions for c.d.f.s are used
76 return _cdfBoundaries ;
77 }
78
79 void setUnitNorm(bool flag) {
80 // Declare contents to have unit normalization
81 _unitNorm = flag ;
82 }
83 bool haveUnitNorm() const {
84 // Return true if contents is declared to be unit normalized
85 return _unitNorm ;
86 }
87
88 bool selfNormalized() const override { return _unitNorm ; }
89
90 Int_t getMaxVal(const RooArgSet& vars) const override ;
91 double maxVal(Int_t code) const override ;
92
93 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override ;
94 std::list<double>* binBoundaries(RooAbsRealLValue& /*obs*/, double /*xlo*/, double /*xhi*/) const override ;
95 bool isBinnedDistribution(const RooArgSet&) const override { return _intOrder==0 ; }
96
97 void computeBatch(double* output, size_t size, RooFit::Detail::DataMap const&) const override;
98
99 void translate(RooFit::Detail::CodeSquashContext &ctx) const override;
100 std::string
101 buildCallToAnalyticIntegral(int code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override;
102
103 protected:
104 bool areIdentical(const RooDataHist& dh1, const RooDataHist& dh2) ;
105
106 bool importWorkspaceHook(RooWorkspace& ws) override ;
107
108 double evaluate() const override;
109 double totalVolume() const ;
110 friend class RooAbsCachedPdf ;
111 double totVolume() const ;
112
113 RooArgSet _histObsList; ///< List of observables defining dimensions of histogram
114 RooSetProxy _pdfObsList; ///< List of observables mapped onto histogram observables
115 RooDataHist* _dataHist = nullptr; ///< Unowned pointer to underlying histogram
116 std::unique_ptr<RooDataHist> _ownedDataHist; ///<! Owned pointer to underlying histogram
117 mutable RooAICRegistry _codeReg ; ///<! Auxiliary class keeping tracking of analytical integration code
118 Int_t _intOrder = 0; ///< Interpolation order
119 bool _cdfBoundaries = false; ///< Use boundary conditions for CDFs.
120 mutable double _totVolume = 0.0; ///<! Total volume of space (product of ranges of observables)
121 bool _unitNorm = false; ///< Assume contents is unit normalized (for use as pdf cache)
122
123private:
124
125 friend class RooHistFunc;
126
127 static bool forceAnalyticalInt(RooArgSet const& pdfObsList, RooAbsArg const& dep);
128
129 static Int_t getAnalyticalIntegral(RooArgSet& allVars,
130 RooArgSet& analVars,
131 const char* rangeName,
132 RooArgSet const& histObsList,
133 RooArgSet const& pdfObsList,
134 Int_t intOrder) ;
135
136 static double analyticalIntegral(Int_t code,
137 const char* rangeName,
138 RooArgSet const& histObsList,
139 RooArgSet const& pdfObsList,
141 bool histFuncMode) ;
142
143 static std::list<double>* plotSamplingHint(RooDataHist const& dataHist,
144 RooArgSet const& pdfObsList,
145 RooArgSet const& histObsList,
146 int intOrder,
147 RooAbsRealLValue& obs,
148 double xlo,
149 double xhi);
150
151 static void rooHistTranslateImpl(RooAbsArg const *klass, RooFit::Detail::CodeSquashContext &ctx, int intOrder,
152 RooDataHist const *dataHist, const RooArgSet &obs, bool correctForBinSize);
153
154 static std::string rooHistIntegralTranslateImpl(int code, RooAbsArg const *klass, RooDataHist const *dataHist,
155 const RooArgSet &obs, bool histFuncMode);
156
157 ClassDefOverride(RooHistPdf,4) // Histogram based PDF
158};
159
160#endif
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
char name[80]
Definition TGX11.cxx:110
RooAICRegistry is a utility class for operator p.d.f classes that keeps track of analytical integrati...
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:79
Abstract base class for p.d.f.s that need or want to cache their evaluate() output in a RooHistPdf de...
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...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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
The RooDataHist is a container class to hold N-dimensional binned data.
Definition RooDataHist.h:39
A class to maintain the context for squashing of RooFit models into code.
RooHistFunc implements a real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:31
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)
Int_t getInterpolationOrder() const
Definition RooHistPdf.h:60
void setUnitNorm(bool flag)
Definition RooHistPdf.h:79
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.
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
bool isBinnedDistribution(const RooArgSet &) const override
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
Definition RooHistPdf.h:95
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 ...
RooAICRegistry _codeReg
! Auxiliary class keeping tracking of analytical integration code
Definition RooHistPdf.h:117
TObject * clone(const char *newname) const override
Definition RooHistPdf.h:40
void setCdfBoundaries(bool flag)
Definition RooHistPdf.h:70
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
bool haveUnitNorm() const
Definition RooHistPdf.h:83
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
double totalVolume() const
bool getCdfBoundaries() const
Definition RooHistPdf.h:74
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:88
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Determine integration scenario.
const RooDataHist & dataHist() const
Definition RooHistPdf.h:47
~RooHistPdf() override
Destructor.
double evaluate() const override
Return the current value: The value of the bin enclosing the current coordinates of the observables,...
void setInterpolationOrder(Int_t order)
Definition RooHistPdf.h:56
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
Persistable container for RooFit projects.
Mother of all ROOT objects.
Definition TObject.h:41
static void output()