Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooStudentT.cxx
Go to the documentation of this file.
1/**
2 * \class RooStudentT
3 * \ingroup Roofit
4 *
5 * Location-scale Student's t-distribution
6 * \see https://en.wikipedia.org/wiki/Student's_t-distribution#Location-scale_t-distribution
7 */
8
9#include "RooStudentT.h"
10#include "RooHelpers.h"
11
12#include "TMath.h"
13
14RooStudentT::RooStudentT(const char *name, const char *title, RooAbsReal::Ref x, RooAbsReal::Ref mean,
16 : RooAbsPdf(name, title),
17 _x("x", "Observable", this, x),
18 _mean("mean", "Mean", this, mean),
19 _sigma("sigma", "Width", this, sigma),
20 _ndf("ndf", " Degrees of freedom", this, ndf)
21{
22 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal &>(sigma)}, 0);
23 RooHelpers::checkRangeOfParameters(this, {&static_cast<RooAbsReal &>(ndf)}, 1);
24}
25
28 _x("x", this, other._x),
29 _mean("mean", this, other._mean),
30 _sigma("sigma", this, other._sigma),
31 _ndf("ndf", this, other._ndf)
32{
33}
34
36{
37 const double t = (_x - _mean) / _sigma;
38 return TMath::Student(t, _ndf);
39}
40
41Int_t RooStudentT::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
42{
43 return matchArgs(allVars, analVars, _x) ? 1 : 0;
44}
45
46double RooStudentT::analyticalIntegral(Int_t code, const char *rangeName) const
47{
48 R__ASSERT(code == 1);
49
50 double tmin = (_x.min(rangeName) - _mean) / _sigma;
51 double tmax = (_x.max(rangeName) - _mean) / _sigma;
52 return (TMath::StudentI(tmax, _ndf) - TMath::StudentI(tmin, _ndf)) * _sigma;
53}
54
55/// Advertise that we know the maximum of self for given (mean,ndf,sigma).
57{
58 RooArgSet dummy;
59
60 if (matchArgs(vars, dummy, _x)) {
61 return 1;
62 }
63 return 0;
64}
65
66double RooStudentT::maxVal(Int_t code) const
67{
68 R__ASSERT(code == 1);
69
70 return TMath::Student(0, _ndf);
71}
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
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
bool matchArgs(const RooArgSet &allDeps, RooArgSet &analDeps, const RooArgProxy &a, const Proxies &... proxies) const
Definition RooAbsReal.h:428
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
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
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
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
void checkRangeOfParameters(const RooAbsReal *callingClass, std::initializer_list< const RooAbsReal * > pars, double min=-std::numeric_limits< double >::max(), double max=std::numeric_limits< double >::max(), bool limitsInAllowedRange=false, std::string const &extraMessage="")
Check if the parameters have a range, and warn if the range extends below / above the set limits.
Double_t Student(Double_t T, Double_t ndf)
Computes density function for Student's t- distribution.
Definition TMath.cxx:2625
Double_t StudentI(Double_t T, Double_t ndf)
Calculates the cumulative distribution function of Student's t-distribution second parameter stands f...
Definition TMath.cxx:2648