Logo ROOT  
Reference Guide
 
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Modules Pages
Loading...
Searching...
No Matches
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/** \class RooNonCentralChiSquare
11 \ingroup Roofit
12
13The PDF of the Non-Central Chi Square distribution for n degrees of freedom.
14It is the asymptotic distribution of the profile likelihood ratio test q_mu
15when a different mu' is true. It is Wald's generalization of Wilks' Theorem.
16
17See:
18
19 Asymptotic formulae for likelihood-based tests of new physics
20
21 By Glen Cowan, Kyle Cranmer, Eilam Gross, Ofer Vitells
22 http://arXiv.org/abs/arXiv:1007.1727
23
24 [Wikipedia](http://en.wikipedia.org/wiki/Noncentral_chi-square_distribution#Approximation)
25
26It requires MathMore to evaluate for non-integer degrees of freedom, k.
27
28When the Mathmore library is available we can use the MathMore libraries implemented using GSL.
29It makes use of the modified Bessel function of the first kind (for k > 2). For k < 2 it uses
30the hypergeometric function 0F1.
31When is not available we use explicit summation of normal chi-squared distributions
32The usage of the sum can be forced by calling SetForceSum(true);
33
34This implementation could be improved. BOOST has a nice implementation:
35
36http://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
37
38http://wesnoth.repositoryhosting.com/trac/wesnoth_wesnoth/browser/trunk/include/boost/math/distributions/non_central_chi_squared.hpp?rev=6
39**/
40
41#include "Riostream.h"
42
44#include "RooAbsReal.h"
45#include "RooAbsCategory.h"
46#include <cmath>
47#include "TMath.h"
48//#include "RooNumber.h"
49#include "Math/DistFunc.h"
50
51
52#include "RooMsgService.h"
53
54#include "TError.h"
55
56using std::endl;
57
58
59////////////////////////////////////////////////////////////////////////////////
60
62 RooAbsReal &_lambda)
63 : RooAbsPdf(name, title),
64 x("x", "x", this, _x),
65 k("k", "k", this, _k),
66 lambda("lambda", "lambda", this, _lambda),
67 fErrorTol(1E-3),
68 fMaxIters(10),
69 fForceSum(false),
70 fHasIssuedConvWarning(false),
71 fHasIssuedSumWarning(false)
72{
73 ccoutD(InputArguments) << "RooNonCentralChiSquare::ctor(" << GetName() <<
74 "MathMore Available, will use Bessel function expressions unless SetForceSum(true) "<< std::endl ;
75}
76
77////////////////////////////////////////////////////////////////////////////////
78
81 x("x", this, other.x),
82 k("k", this, other.k),
83 lambda("lambda", this, other.lambda),
84 fErrorTol(other.fErrorTol),
85 fMaxIters(other.fMaxIters),
86 fForceSum(other.fForceSum),
87 fHasIssuedConvWarning(false),
88 fHasIssuedSumWarning(false)
89{
90 ccoutD(InputArguments) << "RooNonCentralChiSquare::ctor(" << GetName() <<
91 "MathMore Available, will use Bessel function expressions unless SetForceSum(true) "<< std::endl ;
92}
93
94////////////////////////////////////////////////////////////////////////////////
95
99
100////////////////////////////////////////////////////////////////////////////////
101
103{
104 // ENTER EXPRESSION IN TERMS OF VARIABLE ARGUMENTS HERE
105
106
107 // chi^2(0,k) gives inf and causes various problems
108 // truncate
109 double xmin = x.min();
110 double xmax = x.max();
111 double _x = x;
112 if(_x<=0){
113 // options for dealing with this
114 // return 0; // gives a funny dip
115 // _x = 1./RooNumber::infinity(); // too tall
116 _x = xmin + 1e-3*(xmax-xmin); // very small fraction of range
117 }
118
119 // special case (form below doesn't work when lambda==0)
120 if(lambda==0){
121 return ROOT::Math::chisquared_pdf(_x,k);
122 }
123
124 // three forms
125 // FIRST FORM
126 // \sum_i=0^\infty exp(-lambda/2) (\lamda/2)^i chi2(x,k+2i) / i!
127 // could truncate sum
128
129 if ( fForceSum ){
131 coutI(InputArguments) << "RooNonCentralChiSquare sum being forced" << std::endl ;
133 }
134 double sum = 0;
135 double ithTerm = 0;
136 double errorTol = fErrorTol;
137 int MaxIters = fMaxIters;
138 int iDominant = (int) std::floor(lambda/2);
139 // std::cout <<"iDominant: " << iDominant << std::endl;
140
141 // do 0th term last
142 // if(iDominant==0) iDominant = 1;
143 for(int i = iDominant; ; ++i){
144 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)*ROOT::Math::chisquared_pdf(_x,k+2*i)/TMath::Gamma(i+1);
145 sum+=ithTerm;
146 // std::cout <<"progress: " << i << " " << ithTerm/sum << std::endl;
147 if(ithTerm/sum < errorTol)
148 break;
149
150 if( i>iDominant+MaxIters){
153 coutW(Eval) << "RooNonCentralChiSquare did not converge: for x=" << x <<" k="<<k
154 << ", lambda="<<lambda << " fractional error = " << ithTerm/sum
155 << "\n either adjust tolerance with SetErrorTolerance(tol) or max_iter with SetMaxIter(max_it)"
156 << std::endl;
157 }
158 break;
159 }
160 }
161
162 for(int i = iDominant - 1; i >= 0; --i){
163 // std::cout <<"Progress: " << i << " " << ithTerm/sum << std::endl;
164 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)*ROOT::Math::chisquared_pdf(_x,k+2*i)/TMath::Gamma(i+1);
165 sum+=ithTerm;
166 }
167
168
169 return sum;
170 }
171
172 // SECOND FORM (use MathMore function based on Bessel function (if k>2) or
173 // or regularized confluent hypergeometric limit function.
175}
176
177////////////////////////////////////////////////////////////////////////////////
178
180{
181 if (matchArgs(allVars,analVars,x)) return 1 ;
182 return 0 ;
183}
184
185////////////////////////////////////////////////////////////////////////////////
186
188{
189 R__ASSERT(code==1 );
190 // std::cout << "evaluating analytic integral" << std::endl;
191 double xmin = x.min(rangeName);
192 double xmax = x.max(rangeName);
193
194 // if xmin~0 and xmax big, then can return 1. b/c evaluate is normalized.
195
196 // special case (form below doesn't work when lambda==0)
197 if(lambda==0){
199 }
200
201 // three forms
202 // FIRST FORM
203 // \sum_i=0^\infty exp(-lambda/2) (\lamda/2)^i chi2(x,k+2i) / i!
204 // could truncate sum
205
206
207 double sum = 0;
208 double ithTerm = 0;
209 double errorTol = fErrorTol; // for normalization allow slightly larger error
210 int MaxIters = fMaxIters; // for normalization use more terms
211
212 int iDominant = (int) std::floor(lambda/2);
213 // std::cout <<"iDominant: " << iDominant << std::endl;
214 // iDominant=0;
215 for(int i = iDominant; ; ++i){
216 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)
219 sum+=ithTerm;
220 // std::cout <<"progress: " << i << " " << ithTerm << " " << sum << std::endl;
221 if(ithTerm/sum < errorTol)
222 break;
223
224 if( i>iDominant+MaxIters){
227 coutW(Eval) << "RooNonCentralChiSquare Normalization did not converge: for k="<<k
228 << ", lambda="<<lambda << " fractional error = " << ithTerm/sum
229 << "\n either adjust tolerance with SetErrorTolerance(tol) or max_iter with SetMaxIter(max_it)"
230 << std::endl;
231 }
232 break;
233 }
234 }
235
236 for(int i = iDominant - 1; i >= 0; --i){
237 ithTerm =exp(-lambda/2.)*pow(lambda/2.,i)
240 sum+=ithTerm;
241 }
242 return sum;
243}
#define e(i)
Definition RSha256.hxx:103
#define coutI(a)
#define coutW(a)
#define ccoutD(a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
#define R__ASSERT(e)
Checks condition e and reports a fatal error if it's false.
Definition TError.h:125
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:59
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:24
The PDF of the Non-Central 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.
Int_t getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
double chisquared_pdf(double x, double r, double x0=0)
Probability density function of the distribution with degrees of freedom.
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 chisquared_cdf(double x, double r, double x0=0)
Cumulative distribution function of the distribution with degrees of freedom (lower tail).
Double_t x[n]
Definition legend1.C:17
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition TMath.cxx:353
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345