Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooChiSquarePdf.cxx
Go to the documentation of this file.
1/*****************************************************************************
2 * Project: RooFit *
3 * Package: RooFitModels *
4 * @(#)root/roofit:$Id$
5 * Authors: *
6 * Kyle Cranmer
7 * *
8 *****************************************************************************/
9
10/** \class RooChiSquarePdf
11 \ingroup Roofit
12
13The PDF of the Chi Square distribution for n degrees of freedom.
14Oddly, this is hard to find in ROOT (except via relation to GammaDist).
15Here we also implement the analytic integral.
16**/
17
18#include "RooChiSquarePdf.h"
19#include "RooRealVar.h"
20#include "RooBatchCompute.h"
21
22#include "TMath.h"
23
24#include <array>
25#include <cmath>
26
28
29////////////////////////////////////////////////////////////////////////////////
30
32{
33}
34
35////////////////////////////////////////////////////////////////////////////////
36
37RooChiSquarePdf::RooChiSquarePdf(const char* name, const char* title,
38 RooAbsReal& x, RooAbsReal& ndof):
39 RooAbsPdf(name, title),
40 _x("x", "Dependent", this, x),
41 _ndof("ndof","ndof", this, ndof)
42{
43}
44
45////////////////////////////////////////////////////////////////////////////////
46
48 RooAbsPdf(other, name),
49 _x("x", this, other._x),
50 _ndof("ndof",this,other._ndof)
51{
52}
53
54////////////////////////////////////////////////////////////////////////////////
55
57{
58 if (_x <= 0)
59 return 0;
60
61 return pow(_x, (_ndof / 2.) - 1.) * std::exp(-_x / 2.) / TMath::Gamma(_ndof / 2.) / std::pow(2., _ndof / 2.);
62}
63
64////////////////////////////////////////////////////////////////////////////////
65/// Compute multiple values of ChiSquare distribution.
67{
68 std::array<double, 1> extraArgs{_ndof};
69 RooBatchCompute::compute(ctx.config(this), RooBatchCompute::ChiSquare, ctx.output(), {ctx.at(_x)}, extraArgs);
70}
71
72////////////////////////////////////////////////////////////////////////////////
73/// No analytical calculation available (yet) of integrals over subranges
74
75Int_t RooChiSquarePdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
76{
77 if (rangeName && strlen(rangeName)) {
78 return 0 ;
79 }
80
81 if (matchArgs(allVars, analVars, _x)) return 1;
82 return 0;
83}
84
85////////////////////////////////////////////////////////////////////////////////
86
87double RooChiSquarePdf::analyticalIntegral(Int_t code, const char* rangeName) const
88{
89 assert(1 == code); (void)code;
90 double xmin = _x.min(rangeName); double xmax = _x.max(rangeName);
91
92 // TMath::Prob needs ndof to be an integer, or it returns 0.
93 // return TMath::Prob(xmin, _ndof) - TMath::Prob(xmax,_ndof);
94
95 // cumulative is known based on lower incomplete gamma function, or regularized gamma function
96 // Wikipedia defines lower incomplete gamma function without the normalization 1/Gamma(ndof),
97 // but it is included in the ROOT implementation.
98 double pmin = TMath::Gamma(_ndof/2,xmin/2);
99 double pmax = TMath::Gamma(_ndof/2,xmax/2);
100
101 // only use this if range is appropriate
102 return pmax-pmin;
103}
#define ClassImp(name)
Definition Rtypes.h:377
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:55
The PDF of the Chi Square distribution for n degrees of freedom.
double evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
No analytical calculation available (yet) of integrals over subranges.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy _ndof
RooRealProxy _x
void doEval(RooFit::EvalContext &) const override
Compute multiple values of ChiSquare distribution.
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
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.
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition TMath.cxx:353