Loading [MathJax]/extensions/tex2jax.js
ROOT  6.06/09
Reference Guide
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
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 //
12 // BEGIN_HTML
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 // END_HTML
17 //
18 
19 #include "RooFit.h"
20 
21 #include "Riostream.h"
22 #include "Riostream.h"
23 #include <math.h>
24 #include "TMath.h"
25 #include "RooChiSquarePdf.h"
26 #include "RooAbsReal.h"
27 #include "RooRealVar.h"
28 
29 #include "TError.h"
30 
31 using namespace std;
32 
34 ;
35 
36 
37 ////////////////////////////////////////////////////////////////////////////////
38 
40 {
41 }
42 
43 
44 ////////////////////////////////////////////////////////////////////////////////
45 
46 RooChiSquarePdf::RooChiSquarePdf(const char* name, const char* title,
47  RooAbsReal& x, RooAbsReal& ndof):
48  RooAbsPdf(name, title),
49  _x("x", "Dependent", this, x),
50  _ndof("ndof","ndof", this, ndof)
51 {
52 }
53 
54 
55 
56 ////////////////////////////////////////////////////////////////////////////////
57 
59  RooAbsPdf(other, name),
60  _x("x", this, other._x),
61  _ndof("ndof",this,other._ndof)
62 {
63 }
64 
65 
66 ////////////////////////////////////////////////////////////////////////////////
67 
69 {
70  if(_x <= 0) return 0;
71 
72  return pow(_x,(_ndof/2.)-1.) * exp(-_x/2.) / TMath::Gamma(_ndof/2.) / pow(2.,_ndof/2.);
73 
74 
75 }
76 
77 
78 ////////////////////////////////////////////////////////////////////////////////
79 /// No analytical calculation available (yet) of integrals over subranges
80 
81 Int_t RooChiSquarePdf::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* rangeName) const
82 {
83  if (rangeName && strlen(rangeName)) {
84  return 0 ;
85  }
86 
87 
88  if (matchArgs(allVars, analVars, _x)) return 1;
89  return 0;
90 }
91 
92 
93 ////////////////////////////////////////////////////////////////////////////////
94 
95 Double_t RooChiSquarePdf::analyticalIntegral(Int_t code, const char* rangeName) const
96 {
97  R__ASSERT(code==1) ;
98  Double_t xmin = _x.min(rangeName); Double_t xmax = _x.max(rangeName);
99 
100  // TMath::Prob needs ndof to be an integer, or it returns 0.
101  // return TMath::Prob(xmin, _ndof) - TMath::Prob(xmax,_ndof);
102 
103  // cumulative is known based on lower incomplete gamma function, or regularized gamma function
104  // Wikipedia defines lower incomplete gamma function without the normalization 1/Gamma(ndof),
105  // but it is included in the ROOT implementation.
106  Double_t pmin = TMath::Gamma(_ndof/2,xmin/2);
107  Double_t pmax = TMath::Gamma(_ndof/2,xmax/2);
108 
109  // only use this if range is appropriate
110  return pmax-pmin;
111 }
112 
float xmin
Definition: THbookFile.cxx:93
#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 x[n]
Definition: legend1.C:17
RooRealProxy _ndof
double pow(double, double)
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
No analytical calculation available (yet) of integrals over subranges.
ClassImp(RooChiSquarePdf)
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
float xmax
Definition: THbookFile.cxx:93
Double_t evaluate() const
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
#define name(a, b)
Definition: linkTestLib0.cpp:5
RooAbsPdf is the abstract interface for all probability density functions The class provides hybrid a...
Definition: RooAbsPdf.h:41
Double_t min(const char *rname=0) const
Definition: RooRealProxy.h:56
RooRealProxy _x
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
double exp(double)
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57