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 <cmath>
25using namespace std;
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) return 0;
59
60 return pow(_x,(_ndof/2.)-1.) * exp(-_x/2.) / TMath::Gamma(_ndof/2.) / pow(2.,_ndof/2.);
61}
62
63////////////////////////////////////////////////////////////////////////////////
64/// Compute multiple values of ChiSquare distribution.
65void RooChiSquarePdf::computeBatch(double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
66{
68 RooBatchCompute::compute(dataMap.config(this), RooBatchCompute::ChiSquare, output, nEvents, {dataMap.at(_x)}, extraArgs);
69}
70
71////////////////////////////////////////////////////////////////////////////////
72/// No analytical calculation available (yet) of integrals over subranges
73
74Int_t RooChiSquarePdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
75{
76 if (rangeName && strlen(rangeName)) {
77 return 0 ;
78 }
79
80 if (matchArgs(allVars, analVars, _x)) return 1;
81 return 0;
82}
83
84////////////////////////////////////////////////////////////////////////////////
85
86double RooChiSquarePdf::analyticalIntegral(Int_t code, const char* rangeName) const
87{
88 assert(1 == code); (void)code;
89 double xmin = _x.min(rangeName); double xmax = _x.max(rangeName);
90
91 // TMath::Prob needs ndof to be an integer, or it returns 0.
92 // return TMath::Prob(xmin, _ndof) - TMath::Prob(xmax,_ndof);
93
94 // cumulative is known based on lower incomplete gamma function, or regularized gamma function
95 // Wikipedia defines lower incomplete gamma function without the normalization 1/Gamma(ndof),
96 // but it is included in the ROOT implementation.
97 double pmin = TMath::Gamma(_ndof/2,xmin/2);
98 double pmax = TMath::Gamma(_ndof/2,xmax/2);
99
100 // only use this if range is appropriate
101 return pmax-pmin;
102}
#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.
void computeBatch(double *output, size_t nEvents, RooFit::Detail::DataMap const &) const override
Compute multiple values of ChiSquare distribution.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooRealProxy _ndof
RooRealProxy _x
RooBatchCompute::Config config(RooAbsArg const *arg) const
Definition DataMap.cxx:40
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
std::vector< double > ArgVector
void compute(Config cfg, Computer comp, RestrictArr output, size_t size, const VarVector &vars, ArgVector &extraArgs)
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition TMath.cxx:353
static void output()