Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooPoisson.h
Go to the documentation of this file.
1 /*****************************************************************************
2 * Project: RooFit *
3 * *
4 * Simple Poisson PDF
5 * author: Kyle Cranmer <cranmer@cern.ch>
6 * *
7 *****************************************************************************/
8
9#ifndef ROOPOISSON
10#define ROOPOISSON
11
12#include "RooAbsPdf.h"
13#include "RooRealProxy.h"
14#include "RooCategoryProxy.h"
15#include "RooAbsReal.h"
16#include "RooAbsCategory.h"
17#include "RooTrace.h"
18
19class RooPoisson : public RooAbsPdf {
20public:
21 RooPoisson() : _noRounding{false} {}
22 // Original constructor without RooAbsReal::Ref for backwards compatibility.
23 inline RooPoisson(const char *name, const char *title, RooAbsReal& _x, RooAbsReal& _mean, bool noRounding=false)
24 : RooPoisson{name, title, RooAbsReal::Ref{_x}, RooAbsReal::Ref{_mean}, noRounding} {}
25 RooPoisson(const char *name, const char *title, RooAbsReal::Ref _x, RooAbsReal::Ref _mean, bool noRounding=false);
26 RooPoisson(const RooPoisson& other, const char* name=nullptr) ;
27 TObject* clone(const char* newname) const override { return new RooPoisson(*this,newname); }
28
29 Int_t getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName=nullptr) const override;
30 double analyticalIntegral(Int_t code, const char* rangeName=nullptr) const override;
31
32 Int_t getGenerator(const RooArgSet& directVars, RooArgSet &generateVars, bool staticInitOK=true) const override;
33 void generateEvent(Int_t code) override;
34
35 /// Switch off/on rounding of `x` to the nearest integer.
36 void setNoRounding(bool flag = true) {_noRounding = flag;}
37 bool getNoRounding() const { return _noRounding; }
38
39 /// Switch on or off protection against negative means.
40 void protectNegativeMean(bool flag = true) {_protectNegative = flag;}
41
42 /// Get the x variable.
43 RooAbsReal const& getX() const { return x.arg(); }
44
45 /// Get the mean parameter.
46 RooAbsReal const& getMean() const { return mean.arg(); }
47
48 void translate(RooFit::Detail::CodeSquashContext &ctx) const override;
49 std::string
50 buildCallToAnalyticIntegral(int code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override;
51
52protected:
53
57 bool _protectNegative{true};
58
59 double evaluate() const override;
60 void doEval(RooFit::EvalContext &) const override;
61 inline bool canComputeBatchWithCuda() const override { return true; }
62
63 ClassDefOverride(RooPoisson,3) // A Poisson PDF
64};
65
66#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
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:68
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:55
A class to maintain the context for squashing of RooFit models into code.
Poisson pdf.
Definition RooPoisson.h:19
RooRealProxy x
Definition RooPoisson.h:54
bool _noRounding
Definition RooPoisson.h:56
bool canComputeBatchWithCuda() const override
Definition RooPoisson.h:61
RooAbsReal const & getX() const
Get the x variable.
Definition RooPoisson.h:43
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
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 ...
RooPoisson(const char *name, const char *title, RooAbsReal &_x, RooAbsReal &_mean, bool noRounding=false)
Definition RooPoisson.h:23
bool _protectNegative
Definition RooPoisson.h:57
void setNoRounding(bool flag=true)
Switch off/on rounding of x to the nearest integer.
Definition RooPoisson.h:36
void generateEvent(Int_t code) override
Implement internal generator using TRandom::Poisson.
void doEval(RooFit::EvalContext &) const override
Compute multiple values of the Poisson distribution.
void protectNegativeMean(bool flag=true)
Switch on or off protection against negative means.
Definition RooPoisson.h:40
Int_t getGenerator(const RooArgSet &directVars, RooArgSet &generateVars, bool staticInitOK=true) const override
Advertise internal generator in x.
TObject * clone(const char *newname) const override
Definition RooPoisson.h:27
double evaluate() const override
Implementation in terms of the TMath::Poisson() function.
RooRealProxy mean
Definition RooPoisson.h:55
bool getNoRounding() const
Definition RooPoisson.h:37
std::string buildCallToAnalyticIntegral(int code, const char *rangeName, RooFit::Detail::CodeSquashContext &ctx) const override
This function defines the analytical integral translation for the class.
RooAbsReal const & getMean() const
Get the mean parameter.
Definition RooPoisson.h:46
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
const T & arg() const
Return reference to object held in proxy.
Mother of all ROOT objects.
Definition TObject.h:41