Logo ROOT  
Reference Guide
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 "RooAbsReal.h"
20#include "RooRealVar.h"
21#include "RooBatchCompute.h"
22
23#include "TMath.h"
24
25#include <cmath>
26using namespace std;
27
29
30////////////////////////////////////////////////////////////////////////////////
31
33{
34}
35
36////////////////////////////////////////////////////////////////////////////////
37
38RooChiSquarePdf::RooChiSquarePdf(const char* name, const char* title,
39 RooAbsReal& x, RooAbsReal& ndof):
40 RooAbsPdf(name, title),
41 _x("x", "Dependent", this, x),
42 _ndof("ndof","ndof", this, ndof)
43{
44}
45
46////////////////////////////////////////////////////////////////////////////////
47
49 RooAbsPdf(other, name),
50 _x("x", this, other._x),
51 _ndof("ndof",this,other._ndof)
52{
53}
54
55////////////////////////////////////////////////////////////////////////////////
56
58{
59 if(_x <= 0) return 0;
60
61 return pow(_x,(_ndof/2.)-1.) * exp(-_x/2.) / TMath::Gamma(_ndof/2.) / pow(2.,_ndof/2.);
62}
63
64////////////////////////////////////////////////////////////////////////////////
65/// Compute multiple values of ChiSquare distribution.
66void RooChiSquarePdf::computeBatch(cudaStream_t* stream, double* output, size_t nEvents, RooFit::Detail::DataMap const& dataMap) const
67{
69 dispatch->compute(stream, RooBatchCompute::ChiSquare, output, nEvents, {dataMap.at(_x)}, {_ndof});
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:375
char name[80]
Definition: TGX11.cxx:110
float xmin
Definition: THbookFile.cxx:95
float xmax
Definition: THbookFile.cxx:95
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:64
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:57
virtual void compute(cudaStream_t *, Computer, RestrictArr, size_t, const VarVector &, const ArgVector &={})=0
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.
void computeBatch(cudaStream_t *, 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=0) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const override
No analytical calculation available (yet) of integrals over subranges.
RooRealProxy _ndof
RooRealProxy _x
auto & at(RooAbsArg const *arg, RooAbsArg const *=nullptr)
Definition: DataMap.h:88
double min(const char *rname=0) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
double max(const char *rname=0) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
RVec< PromoteTypes< T0, T1 > > pow(const T0 &x, const RVec< T1 > &v)
Definition: RVec.hxx:1753
RVec< PromoteType< T > > exp(const RVec< T > &v)
Definition: RVec.hxx:1744
Double_t x[n]
Definition: legend1.C:17
void(off) SmallVectorTemplateBase< T
R__EXTERN RooBatchComputeInterface * dispatchCUDA
R__EXTERN RooBatchComputeInterface * dispatchCPU
This dispatch pointer points to an implementation of the compute library, provided one has been loade...
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition: TMath.cxx:353
static void output()