Logo ROOT   6.14/05
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 
13 The PDF of the Chi Square distribution for n degrees of freedom.
14 Oddly, this is hard to find in ROOT (except via relation to GammaDist).
15 Here we also implement the analytic integral.
16 **/
17 
18 #include "RooFit.h"
19 
20 #include "Riostream.h"
21 #include "Riostream.h"
22 #include <math.h>
23 #include "TMath.h"
24 #include "RooChiSquarePdf.h"
25 #include "RooAbsReal.h"
26 #include "RooRealVar.h"
27 
28 #include "TError.h"
29 
30 using namespace std;
31 
33 
34 ////////////////////////////////////////////////////////////////////////////////
35 
37 {
38 }
39 
40 ////////////////////////////////////////////////////////////////////////////////
41 
42 RooChiSquarePdf::RooChiSquarePdf(const char* name, const char* title,
43  RooAbsReal& x, RooAbsReal& ndof):
44  RooAbsPdf(name, title),
45  _x("x", "Dependent", this, x),
46  _ndof("ndof","ndof", this, ndof)
47 {
48 }
49 
50 ////////////////////////////////////////////////////////////////////////////////
51 
53  RooAbsPdf(other, name),
54  _x("x", this, other._x),
55  _ndof("ndof",this,other._ndof)
56 {
57 }
58 
59 ////////////////////////////////////////////////////////////////////////////////
60 
62 {
63  if(_x <= 0) return 0;
64 
65  return pow(_x,(_ndof/2.)-1.) * exp(-_x/2.) / TMath::Gamma(_ndof/2.) / pow(2.,_ndof/2.);
66 }
67 
68 ////////////////////////////////////////////////////////////////////////////////
69 /// No analytical calculation available (yet) of integrals over subranges
70 
71 Int_t RooChiSquarePdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
72 {
73  if (rangeName && strlen(rangeName)) {
74  return 0 ;
75  }
76 
77  if (matchArgs(allVars, analVars, _x)) return 1;
78  return 0;
79 }
80 
81 ////////////////////////////////////////////////////////////////////////////////
82 
83 Double_t RooChiSquarePdf::analyticalIntegral(Int_t code, const char* rangeName) const
84 {
85  R__ASSERT(code==1) ;
86  Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
87 
88  // TMath::Prob needs ndof to be an integer, or it returns 0.
89  // return TMath::Prob(xmin, _ndof) - TMath::Prob(xmax,_ndof);
90 
91  // cumulative is known based on lower incomplete gamma function, or regularized gamma function
92  // Wikipedia defines lower incomplete gamma function without the normalization 1/Gamma(ndof),
93  // but it is included in the ROOT implementation.
94  Double_t pmin = TMath::Gamma(_ndof/2,xmin/2);
95  Double_t pmax = TMath::Gamma(_ndof/2,xmax/2);
96 
97  // only use this if range is appropriate
98  return pmax-pmin;
99 }
float xmin
Definition: THbookFile.cxx:93
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
Double_t evaluate() const
#define R__ASSERT(e)
Definition: TError.h:96
int Int_t
Definition: RtypesCore.h:41
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition: TMath.cxx:348
STL namespace.
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
No analytical calculation available (yet) of integrals over subranges.
Double_t x[n]
Definition: legend1.C:17
RooRealProxy _ndof
double pow(double, double)
float xmax
Definition: THbookFile.cxx:93
#define ClassImp(name)
Definition: Rtypes.h:359
The PDF of the Chi Square distribution for n degrees of freedom.
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
double Double_t
Definition: RtypesCore.h:55
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition: RooAbsReal.h:53
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
RooRealProxy _x
double exp(double)
char name[80]
Definition: TGX11.cxx:109