Logo ROOT   6.08/07
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 /**
11 \file RooChiSquarePdf.cxx
12 \class RooChiSquarePdf
13 \ingroup Roofit
14 
15 The PDF of the Chi Square distribution for n degrees of freedom.
16 Oddly, this is hard to find in ROOT (except via relation to GammaDist).
17 Here we also implement the analytic integral.
18 **/
19 
20 #include "RooFit.h"
21 
22 #include "Riostream.h"
23 #include "Riostream.h"
24 #include <math.h>
25 #include "TMath.h"
26 #include "RooChiSquarePdf.h"
27 #include "RooAbsReal.h"
28 #include "RooRealVar.h"
29 
30 #include "TError.h"
31 
32 using namespace std;
33 
35 ;
36 
37 
38 ////////////////////////////////////////////////////////////////////////////////
39 
41 {
42 }
43 
44 
45 ////////////////////////////////////////////////////////////////////////////////
46 
47 RooChiSquarePdf::RooChiSquarePdf(const char* name, const char* title,
48  RooAbsReal& x, RooAbsReal& ndof):
49  RooAbsPdf(name, title),
50  _x("x", "Dependent", this, x),
51  _ndof("ndof","ndof", this, ndof)
52 {
53 }
54 
55 
56 
57 ////////////////////////////////////////////////////////////////////////////////
58 
60  RooAbsPdf(other, name),
61  _x("x", this, other._x),
62  _ndof("ndof",this,other._ndof)
63 {
64 }
65 
66 
67 ////////////////////////////////////////////////////////////////////////////////
68 
70 {
71  if(_x <= 0) return 0;
72 
73  return pow(_x,(_ndof/2.)-1.) * exp(-_x/2.) / TMath::Gamma(_ndof/2.) / pow(2.,_ndof/2.);
74 
75 
76 }
77 
78 
79 ////////////////////////////////////////////////////////////////////////////////
80 /// No analytical calculation available (yet) of integrals over subranges
81 
82 Int_t RooChiSquarePdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
83 {
84  if (rangeName && strlen(rangeName)) {
85  return 0 ;
86  }
87 
88 
89  if (matchArgs(allVars, analVars, _x)) return 1;
90  return 0;
91 }
92 
93 
94 ////////////////////////////////////////////////////////////////////////////////
95 
96 Double_t RooChiSquarePdf::analyticalIntegral(Int_t code, const char* rangeName) const
97 {
98  R__ASSERT(code==1) ;
99  Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
100 
101  // TMath::Prob needs ndof to be an integer, or it returns 0.
102  // return TMath::Prob(xmin, _ndof) - TMath::Prob(xmax,_ndof);
103 
104  // cumulative is known based on lower incomplete gamma function, or regularized gamma function
105  // Wikipedia defines lower incomplete gamma function without the normalization 1/Gamma(ndof),
106  // but it is included in the ROOT implementation.
107  Double_t pmin = TMath::Gamma(_ndof/2,xmin/2);
108  Double_t pmax = TMath::Gamma(_ndof/2,xmax/2);
109 
110  // only use this if range is appropriate
111  return pmax-pmin;
112 }
113 
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:98
int Int_t
Definition: RtypesCore.h:41
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition: TMath.cxx:352
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:279
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