ROOT  6.06/09
Reference Guide
RooNonCentralChiSquare.cxx
Go to the documentation of this file.
1 /*****************************************************************************
2  * Project: RooFit *
3  * Package: RooFitModels *
4  * @(#)root/roofit:$Id: RooNonCentralChiSquare *
5  * Authors: *
6  * Kyle Cranmer
7  * *
8  *****************************************************************************/
9 
10 //////////////////////////////////////////////////////////////////////////////
11 //
12 // BEGIN_HTML
13 // The PDF of the Non-Central Chi Square distribution for n degrees of freedom.
14 // It is the asymptotic distribution of the profile likeihood ratio test q_mu
15 // when a different mu' is true. It is Wald's generalization of Wilks' Theorem.
16 //
17 // See:
18 // Asymptotic formulae for likelihood-based tests of new physics
19 // By Glen Cowan, Kyle Cranmer, Eilam Gross, Ofer Vitells
20 // http://arXiv.org/abs/arXiv:1007.1727
21 //
22 // Wikipedia:
23 // http://en.wikipedia.org/wiki/Noncentral_chi-square_distribution#Approximation
24 //
25 // It requires MathMore to evaluate for non-integer degrees of freedom, k.
26 //
27 // When the Mathmore library is available we can use the MathMore libraries impelmented using GSL.
28 // It makes use of the modified Bessel function of the first kind (for k > 2). For k < 2 it uses
29 // the hypergeometric function 0F1.
30 // When is not available we use explicit summation of normal chi-squared distributions
31 // The usage of the sum can be forced by calling SetForceSum(true);
32 //
33 // This implementation could be improved. BOOST has a nice implementation:
34 // http://live.boost.org/doc/libs/1_42_0/libs/math/doc/sf_and_dist/html/math_toolkit/dist/dist_ref/dists/nc_chi_squared_dist.html
35 // http://wesnoth.repositoryhosting.com/trac/wesnoth_wesnoth/browser/trunk/include/boost/math/distributions/non_central_chi_squared.hpp?rev=6
36 // END_HTML
37 //
38 
39 #include "Riostream.h"
40 
41 #include "RooNonCentralChiSquare.h"
42 #include "RooAbsReal.h"
43 #include "RooAbsCategory.h"
44 #include <math.h>
45 #include "TMath.h"
46 //#include "RooNumber.h"
47 #include "Math/DistFunc.h"
48 
49 
50 #include "RooMsgService.h"
51 
52 #include "TError.h"
53 
54 using namespace std;
55 
57 
58 RooNonCentralChiSquare::RooNonCentralChiSquare(const char *name, const char *title,
59  RooAbsReal& _x,
60  RooAbsReal& _k,
61  RooAbsReal& _lambda) :
62  RooAbsPdf(name,title),
63  x("x","x",this,_x),
64  k("k","k",this,_k),
65  lambda("lambda","lambda",this,_lambda),
66  fErrorTol(1E-3),
67  fMaxIters(10),
68  fHasIssuedConvWarning(false),
69  fHasIssuedSumWarning(false)
70 {
71 #ifdef R__HAS_MATHMORE
72  ccoutD(InputArguments) << "RooNonCentralChiSquare::ctor(" << GetName() <<
73  "MathMore Available, will use Bessel function expressions unless SetForceSum(true) "<< endl ;
74  fForceSum = false;
75 #else
76  fForceSum = true;
77 #endif
78 }
79 
81  RooAbsPdf(other,name),
82  x("x",this,other.x),
83  k("k",this,other.k),
84  lambda("lambda",this,other.lambda),
85  fErrorTol(other.fErrorTol),
86  fMaxIters(other.fMaxIters),
87  fHasIssuedConvWarning(false),
88  fHasIssuedSumWarning(false)
89 {
90 #ifdef R__HAS_MATHMORE
91  ccoutD(InputArguments) << "RooNonCentralChiSquare::ctor(" << GetName() <<
92  "MathMore Available, will use Bessel function expressions unless SetForceSum(true) "<< endl ;
93  fForceSum = other.fForceSum;
94 #else
95  fForceSum = true;
96 #endif
97 }
98 
99 
101  fForceSum = flag;
102 #ifndef R__HAS_MATHMORE
103  if (!fForceSum) {
104  ccoutD(InputArguments) << "RooNonCentralChiSquare::SetForceSum" << GetName() <<
105  "MathMore is not available- ForceSum must be on "<< endl ;
106  fForceSum = true;
107  }
108 #endif
109 
110 }
111 
112 
114 {
115  // ENTER EXPRESSION IN TERMS OF VARIABLE ARGUMENTS HERE
116 
117 
118  // chi^2(0,k) gives inf and causes various problems
119  // truncate
120  Double_t xmin = x.min();
121  Double_t xmax = x.max();
122  double _x = x;
123  if(_x<=0){
124  // options for dealing with this
125  // return 0; // gives a funny dip
126  // _x = 1./RooNumber::infinity(); // too tall
127  _x = xmin + 1e-3*(xmax-xmin); // very small fraction of range
128  }
129 
130  // special case (form below doesn't work when lambda==0)
131  if(lambda==0){
132  return ROOT::Math::chisquared_pdf(_x,k);
133  }
134 
135  // three forms
136  // FIRST FORM
137  // \sum_i=0^\infty exp(-lambda/2) (\lamda/2)^i chi2(x,k+2i) / i!
138  // could truncate sum
139 
140  if ( fForceSum ){
142  coutI(InputArguments) << "RooNonCentralChiSquare sum being forced" << endl ;
144  }
145  double sum = 0;
146  double ithTerm = 0;
147  double errorTol = fErrorTol;
148  int MaxIters = fMaxIters;
149  int iDominant = (int) TMath::Floor(lambda/2);
150  // cout <<"iDominant: " << iDominant << endl;
151 
152  // do 0th term last
153  // if(iDominant==0) iDominant = 1;
154  for(int i = iDominant; ; ++i){
155  ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)*ROOT::Math::chisquared_pdf(_x,k+2*i)/TMath::Gamma(i+1);
156  sum+=ithTerm;
157  // cout <<"progress: " << i << " " << ithTerm/sum << endl;
158  if(ithTerm/sum < errorTol)
159  break;
160 
161  if( i>iDominant+MaxIters){
164  coutW(Eval) << "RooNonCentralChiSquare did not converge: for x=" << x <<" k="<<k
165  << ", lambda="<<lambda << " fractional error = " << ithTerm/sum
166  << "\n either adjust tolerance with SetErrorTolerance(tol) or max_iter with SetMaxIter(max_it)"
167  << endl;
168  }
169  break;
170  }
171  }
172 
173  for(int i = iDominant - 1; i >= 0; --i){
174  // cout <<"Progress: " << i << " " << ithTerm/sum << endl;
175  ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)*ROOT::Math::chisquared_pdf(_x,k+2*i)/TMath::Gamma(i+1);
176  sum+=ithTerm;
177  }
178 
179 
180  return sum;
181  }
182 
183  // SECOND FORM (use MathMore function based on Bessel function (if k>2) or
184  // or regularized confluent hypergeometric limit function.
185 #ifdef R__HAS_MATHMORE
187 #else
188  coutF(Eval) << "RooNonCentralChisquare: ForceSum must be set" << endl;
189  return 0;
190 #endif
191 
192 }
193 
194 ////////////////////////////////////////////////////////////////////////////////
195 
196 Int_t RooNonCentralChiSquare::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char* /*rangeName*/) const
197 {
198  if (matchArgs(allVars,analVars,x)) return 1 ;
199  return 0 ;
200 }
201 
202 
203 
204 ////////////////////////////////////////////////////////////////////////////////
205 
206 Double_t RooNonCentralChiSquare::analyticalIntegral(Int_t code, const char* rangeName) const
207 {
208  R__ASSERT(code==1 );
209  // cout << "evaluating analytic integral" << endl;
210  Double_t xmin = x.min(rangeName);
211  Double_t xmax = x.max(rangeName);
212 
213  // if xmin~0 and xmax big, then can return 1. b/c evaluate is normalized.
214 
215  // special case (form below doesn't work when lambda==0)
216  if(lambda==0){
218  }
219 
220  // three forms
221  // FIRST FORM
222  // \sum_i=0^\infty exp(-lambda/2) (\lamda/2)^i chi2(x,k+2i) / i!
223  // could truncate sum
224 
225 
226  double sum = 0;
227  double ithTerm = 0;
228  double errorTol = fErrorTol; // for nomralization allow slightly larger error
229  int MaxIters = fMaxIters; // for normalization use more terms
230 
231  int iDominant = (int) TMath::Floor(lambda/2);
232  // cout <<"iDominant: " << iDominant << endl;
233  // iDominant=0;
234  for(int i = iDominant; ; ++i){
235  ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)
236  *( ROOT::Math::chisquared_cdf(xmax,k+2*i)/TMath::Gamma(i+1)
237  - ROOT::Math::chisquared_cdf(xmin,k+2*i)/TMath::Gamma(i+1) );
238  sum+=ithTerm;
239  // cout <<"progress: " << i << " " << ithTerm << " " << sum << endl;
240  if(ithTerm/sum < errorTol)
241  break;
242 
243  if( i>iDominant+MaxIters){
246  coutW(Eval) << "RooNonCentralChiSquare Normalization did not converge: for k="<<k
247  << ", lambda="<<lambda << " fractional error = " << ithTerm/sum
248  << "\n either adjust tolerance with SetErrorTolerance(tol) or max_iter with SetMaxIter(max_it)"
249  << endl;
250  }
251  break;
252  }
253  }
254 
255  for(int i = iDominant - 1; i >= 0; --i){
256  ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)
257  *( ROOT::Math::chisquared_cdf(xmax,k+2*i)/TMath::Gamma(i+1)
258  -ROOT::Math::chisquared_cdf(xmin,k+2*i)/TMath::Gamma(i+1) );
259  sum+=ithTerm;
260  }
261  return sum;
262 }
263 
264 
265 
266 
267 
float xmin
Definition: THbookFile.cxx:93
Double_t Floor(Double_t x)
Definition: TMath.h:473
double chisquared_pdf(double x, double r, double x0=0)
Probability density function of the distribution with degrees of freedom.
#define coutI(a)
Definition: RooMsgService.h:32
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=0) const
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported...
#define R__ASSERT(e)
Definition: TError.h:98
int Int_t
Definition: RtypesCore.h:41
bool Bool_t
Definition: RtypesCore.h:59
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition: TMath.cxx:352
STL namespace.
#define coutW(a)
Definition: RooMsgService.h:34
ClassImp(TIterator) Bool_t TIterator return false
Compare two iterator objects.
Definition: TIterator.cxx:20
Double_t x[n]
Definition: legend1.C:17
double pow(double, double)
Double_t analyticalIntegral(Int_t code, const char *rangeName=0) const
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral. ...
Double_t E()
Definition: TMath.h:54
ClassImp(RooNonCentralChiSquare) RooNonCentralChiSquare
virtual const char * GetName() const
Returns name of object.
Definition: TNamed.h:51
float xmax
Definition: THbookFile.cxx:93
#define coutF(a)
Definition: RooMsgService.h:36
double noncentral_chisquared_pdf(double x, double r, double lambda)
Probability density function of the non central distribution with degrees of freedom and the noon-c...
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
Bool_t matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
double chisquared_cdf(double x, double r, double x0=0)
Cumulative distribution function of the distribution with degrees of freedom (lower tail)...
double exp(double)
#define ccoutD(a)
Definition: RooMsgService.h:38
Double_t max(const char *rname=0) const
Definition: RooRealProxy.h:57