File indexing completed on 2024-11-15 09:56:48
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016 #ifndef ROO_WRAPPER_PDF
0017 #define ROO_WRAPPER_PDF
0018
0019 #include "RooAbsReal.h"
0020 #include "RooRealProxy.h"
0021 #include "RooAbsPdf.h"
0022 #include <list>
0023
0024 class RooWrapperPdf final : public RooAbsPdf {
0025 public:
0026
0027 RooWrapperPdf() { };
0028
0029
0030
0031
0032
0033
0034
0035
0036
0037
0038
0039 RooWrapperPdf(const char *name, const char *title, RooAbsReal& inputFunction, bool selfNormalized=false) :
0040 RooAbsPdf(name, title),
0041 _func("inputFunction", "Function to be converted into a PDF", this, inputFunction),
0042 _selfNormalized{selfNormalized} { }
0043
0044 RooWrapperPdf(const RooWrapperPdf& other, const char *name = nullptr) :
0045 RooAbsPdf(other, name),
0046 _func("inputFunction", this, other._func),
0047 _selfNormalized{other._selfNormalized} { }
0048
0049 TObject* clone(const char* newname) const override {
0050 return new RooWrapperPdf(*this, newname);
0051 }
0052
0053 bool selfNormalized() const override { return _selfNormalized; }
0054
0055
0056 bool forceAnalyticalInt(const RooAbsArg& dep) const override {
0057 return _func->forceAnalyticalInt(dep);
0058 }
0059 Int_t getAnalyticalIntegralWN(RooArgSet& allVars, RooArgSet& analVars, const RooArgSet* normSet,
0060 const char* rangeName=nullptr) const override {
0061 return _func->getAnalyticalIntegralWN(allVars, analVars, normSet, rangeName);
0062 }
0063 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& numVars,
0064 const char* rangeName=nullptr) const override {
0065 return _func->getAnalyticalIntegral(allVars, numVars, rangeName);
0066 }
0067 double analyticalIntegralWN(Int_t code, const RooArgSet* normSet, const char* rangeName) const override {
0068 return _func->analyticalIntegralWN(code, normSet, rangeName);
0069 }
0070 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override {
0071 return _func->analyticalIntegral(code, rangeName);
0072 }
0073
0074
0075
0076
0077
0078
0079
0080
0081
0082
0083
0084
0085
0086 Int_t getMaxVal(const RooArgSet& vars) const override {
0087 return _func.arg().getMaxVal(vars);
0088 }
0089 double maxVal(Int_t code) const override {
0090 return _func.arg().maxVal(code);
0091 }
0092 Int_t minTrialSamples(const RooArgSet& arGenObs) const override {
0093 return _func.arg().minTrialSamples(arGenObs);
0094 }
0095
0096
0097 bool isBinnedDistribution(const RooArgSet& obs) const override {
0098 return _func.arg().isBinnedDistribution(obs);
0099 }
0100 std::list<double>* binBoundaries(RooAbsRealLValue& obs, double xlo, double xhi) const override {
0101 return _func.arg().binBoundaries(obs, xlo, xhi);
0102 }
0103 std::list<double>* plotSamplingHint(RooAbsRealLValue& obs, double xlo, double xhi) const override {
0104 return _func.arg().plotSamplingHint(obs, xlo, xhi);
0105 }
0106
0107
0108
0109 private:
0110 RooRealProxy _func;
0111 bool _selfNormalized = false;
0112
0113 double evaluate() const override {
0114 return _func;
0115 }
0116
0117 ClassDefOverride(RooWrapperPdf,2);
0118 };
0119
0120 #endif