Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooWrapperPdf.h
Go to the documentation of this file.
1// Author: Stephan Hageboeck, CERN
2/*****************************************************************************
3 * Project: RooFit *
4 * Package: RooFitCore *
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-2018, 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_WRAPPER_PDF
17#define ROO_WRAPPER_PDF
18
19#include "RooAbsReal.h"
20#include "RooRealProxy.h"
21#include "RooAbsPdf.h"
22#include <list>
23
24class RooWrapperPdf final : public RooAbsPdf {
25public:
26
28 /// Construct a new RooWrapperPdf.
29 /// \param[in] name A name to identify this object.
30 /// \param[in] title Title (for e.g. plotting)
31 /// \param[in] inputFunction Any RooAbsReal that should be converted into a PDF. Although it's possible
32 /// \param[in] selfNormalized The return value the RooAbsPdf::selfNormalized() function for the wrapped PDF object.
33 /// If it is `true`, then no automatic normalization will be
34 /// performed when evaluating the function. In this case, the
35 /// effect RooWrapperPdf is not to change the evaluated values,
36 /// but only to wrap the function in something that is of type
37 /// RooAbsPdf, which can be useful if some interface reqiures it.
38 /// to pass a PDF, it only makes sense for non-PDF functions.
39 RooWrapperPdf(const char *name, const char *title, RooAbsReal& inputFunction, bool selfNormalized=false) :
40 RooAbsPdf(name, title),
41 _func("inputFunction", "Function to be converted into a PDF", this, inputFunction),
43
44 RooWrapperPdf(const RooWrapperPdf& other, const char *name = nullptr) :
45 RooAbsPdf(other, name),
46 _func("inputFunction", this, other._func),
48
49 TObject* clone(const char* newname) const override {
50 return new RooWrapperPdf(*this, newname);
51 }
52
53 bool selfNormalized() const override { return _selfNormalized; }
54
55 // Analytical Integration handling
56 bool forceAnalyticalInt(const RooAbsArg& dep) const override {
57 return _func->forceAnalyticalInt(dep);
58 }
59 Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet,
60 const char* rangeName=nullptr) const override {
61 return _func->getAnalyticalIntegralWN(allVars, analVars, normSet, rangeName);
62 }
64 const char* rangeName=nullptr) const override {
65 return _func->getAnalyticalIntegral(allVars, numVars, rangeName);
66 }
67 double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const override {
68 return _func->analyticalIntegralWN(code, normSet, rangeName);
69 }
70 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override {
71 return _func->analyticalIntegral(code, rangeName);
72 }
73
74
75 // Internal toy generation. Since our _func is not a PDF (if it is, it doesn't make sense to use this wrapper),
76 // we cannot do anything.
77 /// Get specialised generator. Since the underlying function is not a PDF, this will always return zero.
78// Int_t getGenerator(const RooArgSet& /*directVars*/, RooArgSet& /*generateVars*/,
79// bool /*staticInitOK = true*/) const override { return 0; }
80// void initGenerator(Int_t /*code*/) override { }
81// void generateEvent(Int_t /*code*/) override { }
82// bool isDirectGenSafe(const RooAbsArg& /*arg*/) const override { return false; }
83
84
85 // Hints for optimized brute-force sampling
86 Int_t getMaxVal(const RooArgSet& vars) const override {
87 return _func.arg().getMaxVal(vars);
88 }
89 double maxVal(Int_t code) const override {
90 return _func.arg().maxVal(code);
91 }
92 Int_t minTrialSamples(const RooArgSet& arGenObs) const override {
93 return _func.arg().minTrialSamples(arGenObs);
94 }
95
96 // Plotting and binning hints
97 bool isBinnedDistribution(const RooArgSet& obs) const override {
98 return _func.arg().isBinnedDistribution(obs);
99 }
100 std::list<double>* binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const override {
101 return _func.arg().binBoundaries(obs, xlo, xhi);
102 }
103 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override {
104 return _func.arg().plotSamplingHint(obs, xlo, xhi);
105 }
106
107
108
109private:
111 bool _selfNormalized = false;
112
113 double evaluate() const override {
114 return _func;
115 }
116
118};
119
120#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
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
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 ...
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
The RooWrapperPdf is a class that can be used to convert a function into a PDF.
RooWrapperPdf(const char *name, const char *title, RooAbsReal &inputFunction, bool selfNormalized=false)
Construct a new RooWrapperPdf.
RooWrapperPdf(const RooWrapperPdf &other, const char *name=nullptr)
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
Int_t getMaxVal(const RooArgSet &vars) const override
Get specialised generator. Since the underlying function is not a PDF, this will always return zero.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &numVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
Int_t getAnalyticalIntegralWN(RooArgSet &allVars, RooArgSet &analVars, const RooArgSet *normSet, const char *rangeName=nullptr) const override
Variant of getAnalyticalIntegral that is also passed the normalization set that should be applied to ...
RooRealProxy _func
bool isBinnedDistribution(const RooArgSet &obs) const override
Tests if the distribution is binned. Unless overridden by derived classes, this always returns false.
std::list< double > * plotSamplingHint(RooAbsRealLValue &obs, double xlo, double xhi) const override
Interface for returning an optional hint for initial sampling points when constructing a curve projec...
bool selfNormalized() const override
Shows if a PDF is self-normalized, which means that no attempt is made to add a normalization term.
TObject * clone(const char *newname) const override
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
bool forceAnalyticalInt(const RooAbsArg &dep) const override
double analyticalIntegralWN(Int_t code, const RooArgSet *normSet, const char *rangeName) const override
Analytical integral with normalization (see RooAbsReal::analyticalIntegralWN() for further informatio...
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
std::list< double > * binBoundaries(RooAbsRealLValue &obs, double xlo, double xhi) const override
Retrieve bin boundaries if this distribution is binned in obs.
Int_t minTrialSamples(const RooArgSet &arGenObs) const override
Mother of all ROOT objects.
Definition TObject.h:41