Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooBernstein.h
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * File: $Id$
5 * Authors: *
6 * Kyle Cranmer
7 * *
8 * *
9 * Redistribution and use in source and binary forms, *
10 * with or without modification, are permitted according to the terms *
11 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
12 *****************************************************************************/
13#ifndef ROO_BERNSTEIN
14#define ROO_BERNSTEIN
15
16#include "RooAbsPdf.h"
17#include "RooTemplateProxy.h"
18#include "RooRealVar.h"
19#include "RooListProxy.h"
20#include "RooAbsRealLValue.h"
21#include <string>
22
23class RooRealVar;
24class RooArgList;
25
26class RooBernstein : public RooAbsPdf {
27public:
28
30 RooBernstein(const char *name, const char *title,
32
33 RooBernstein(const RooBernstein &other, const char *name = nullptr);
34
35 TObject* clone(const char* newname) const override { return new RooBernstein(*this, newname); }
36
37 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override ;
38 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override ;
39 void selectNormalizationRange(const char* rangeName=nullptr, bool force=false) override ;
40
41private:
42
45 std::string _refRangeName ;
46
47 double evaluate() const override;
48 void computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const&) const override;
49 inline bool canComputeBatchWithCuda() const override { return true; }
50
51 ClassDefOverride(RooBernstein,2) // Bernstein polynomial PDF
52};
53
54#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:341
char name[80]
Definition TGX11.cxx:110
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...
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
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 computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of Bernstein distribution.
TObject * clone(const char *newname) const override
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.
RooTemplateProxy< RooAbsRealLValue > _x
RooListProxy _coefList
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:37
Mother of all ROOT objects.
Definition TObject.h:41
static void output()