Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBernstein.h
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2024, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11#ifndef RooFit_RooBernstein_h
12#define RooFit_RooBernstein_h
13
14#include <RooAbsPdf.h>
15#include <RooAbsRealLValue.h>
16#include <RooListProxy.h>
17#include <RooRealVar.h>
18#include <RooTemplateProxy.h>
19
20#include <string>
21
22class RooBernstein : public RooAbsPdf {
23public:
24 RooBernstein() = default;
25 RooBernstein(const char *name, const char *title, RooAbsRealLValue &_x, const RooArgList &_coefList);
26
27 RooBernstein(const RooBernstein &other, const char *name = nullptr);
28
29 TObject *clone(const char *newname) const override { return new RooBernstein(*this, newname); }
30
31 Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName = nullptr) const override;
32 double analyticalIntegral(Int_t code, const char *rangeName = nullptr) const override;
33 void selectNormalizationRange(const char *rangeName = nullptr, bool force = false) override;
34
35 RooAbsRealLValue const &x() const { return *_x; }
36 RooArgList const &coefList() const { return _coefList; }
37
38 // Implementation detail. Do not use.
39 void fillBuffer() const;
40 // Implementation detail. Do not use.
41 inline double xmin() const { return _buffer[_coefList.size()]; }
42 // Implementation detail. Do not use.
43 inline double xmax() const { return _buffer[_coefList.size() + 1]; }
44
45private:
46
49 std::string _refRangeName;
50 mutable std::vector<double> _buffer; ///<!
51
52 double evaluate() const override;
53 void doEval(RooFit::EvalContext &) const override;
54 inline bool canComputeBatchWithCuda() const override { return true; }
55
56 ClassDefOverride(RooBernstein, 2) // Bernstein polynomial PDF
57};
58
59#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:346
char name[80]
Definition TGX11.cxx:110
Storage_t::size_type size() const
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 ...
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:24
Bernstein basis polynomials are positive-definite in the range [0,1].
void selectNormalizationRange(const char *rangeName=nullptr, bool force=false) override
Force use of a given normalisation range.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
void fillBuffer() const
TObject * clone(const char *newname) const override
RooAbsRealLValue const & x() const
RooBernstein()=default
RooArgList const & coefList() const
bool canComputeBatchWithCuda() const override
std::string _refRangeName
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
void doEval(RooFit::EvalContext &) const override
Compute multiple values of Bernstein distribution.
double xmax() const
std::vector< double > _buffer
!
RooTemplateProxy< RooAbsRealLValue > _x
double xmin() const
RooListProxy _coefList
Mother of all ROOT objects.
Definition TObject.h:41