Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooStudentT.h
Go to the documentation of this file.
1#ifndef RooFit_RooFit_RooStudentT_h
2#define RooFit_RooFit_RooStudentT_h
3
4#include "RooAbsPdf.h"
5#include "RooRealProxy.h"
6
7class RooStudentT : public RooAbsPdf {
8public:
12 RooStudentT(const RooStudentT &other, const char *name = nullptr);
13 TObject *clone(const char *newname = nullptr) const override { return new RooStudentT(*this, newname); }
14
15 Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName = nullptr) const override;
16 double analyticalIntegral(Int_t code, const char *rangeName = nullptr) const override;
17
18 /// Get the _x variable.
19 RooAbsReal const &x() const { return _x.arg(); }
20
21 /// Get the _mean parameter.
22 RooAbsReal const &mean() const { return _mean.arg(); }
23
24 /// Get the standard deviation parameter.
25 RooAbsReal const &sigma() const { return _sigma.arg(); }
26
27 /// Get the degrees of freedom parameter.
28 RooAbsReal const &ndf() const { return _ndf.arg(); }
29
30 Int_t getMaxVal(const RooArgSet &vars) const override;
31 double maxVal(Int_t code) const override;
32
33protected:
34 double evaluate() const override;
35
36private:
37 RooRealProxy _x; ///< variable
39 RooRealProxy _sigma; ///< standard deviation
40 RooRealProxy _ndf; ///< degrees of freedom
41
42 ClassDefOverride(RooStudentT, 1) // Location-scale Student's t-distribution PDF
43};
44
45#endif
#define ClassDefOverride(name, id)
Definition Rtypes.h:348
char name[80]
Definition TGX11.cxx:110
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
A RooAbsReal::Ref can be constructed from a RooAbsReal& or a double that will be implicitly converted...
Definition RooAbsReal.h:72
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
Location-scale Student's t-distribution.
Definition RooStudentT.h:7
RooAbsReal const & mean() const
Get the _mean parameter.
Definition RooStudentT.h:22
RooAbsReal const & x() const
Get the _x variable.
Definition RooStudentT.h:19
RooRealProxy _ndf
degrees of freedom
Definition RooStudentT.h:40
RooRealProxy _sigma
standard deviation
Definition RooStudentT.h:39
RooRealProxy _x
variable
Definition RooStudentT.h:37
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getMaxVal(const RooArgSet &vars) const override
Advertise that we know the maximum of self for given (mean,ndf,sigma).
RooRealProxy _mean
mean
Definition RooStudentT.h:38
RooAbsReal const & sigma() const
Get the standard deviation parameter.
Definition RooStudentT.h:25
TObject * clone(const char *newname=nullptr) const override
Definition RooStudentT.h:13
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
double maxVal(Int_t code) const override
Return maximum value for set of observables identified by code assigned in getMaxVal.
RooAbsReal const & ndf() const
Get the degrees of freedom parameter.
Definition RooStudentT.h:28
const T & arg() const
Return reference to object held in proxy.
Mother of all ROOT objects.
Definition TObject.h:42