STK++ 0.9.13
STK_Funct_gammaRatio.cpp
Go to the documentation of this file.
1/*--------------------------------------------------------------------*/
2/* Copyright (C) 2004-2016 Serge Iovleff, Université Lille 1, Inria
3
4 This program is free software; you can redistribute it and/or modify
5 it under the terms of the GNU Lesser General Public License as
6 published by the Free Software Foundation; either version 2 of the
7 License, or (at your option) any later version.
8
9 This program is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU Lesser General Public License for more details.
13
14 You should have received a copy of the GNU Lesser General Public
15 License along with this program; if not, write to the
16 Free Software Foundation, Inc.,
17 59 Temple Place,
18 Suite 330,
19 Boston, MA 02111-1307
20 USA
21
22 Contact : S..._Dot_I..._At_stkpp_Dot_org (see copyright for ...)
23*/
24
25/*
26 * Project: Analysis
27 * Purpose: Implementation of functions around the gamma function ratio
28 * Author: Serge Iovleff, S..._Dot_I..._At_stkpp_Dot_org (see copyright for ...)
29 **/
30
36#ifndef IS_RTKPP_LIBRARY
37
38
39#include "../include/STK_Funct_gammaRatio.h"
40#include "../include/STK_Funct_raw.h"
41
42namespace STK
43{
44
45namespace Funct
46{
62static Real apois(Real const& a, Real const& b)
63{
64 // special value
65 if (b==0.) return( 0. );
66 // stirling approximation and deviance
67 return( sqrt(a)*exp(-lgammaStirlingError(a)-dev0(a, b))* (Const::_1_SQRT2PI_));
68}
69
77class Seriedl: public ISerie<Seriedl>
78{
79 private:
80 Real const& x_;
81 mutable int n_;
82 mutable Real xn_factn_;// x^n/n!
83 mutable Real aplusn_; // a+n
84
85 public:
87 inline Seriedl( Real const& a, Real const& x)
88 : x_(x), n_(1), xn_factn_(x), aplusn_(a+1.)
89 {}
91 inline Real firstImpl() const { return x_/aplusn_;}
93 inline Real nextImpl() const
94 {
95 // n and x^n/ n!
96 xn_factn_ *= (x_/(++n_));
97 // x^n/((a+n)n!)
98 return xn_factn_/(++aplusn_);
99 }
100};
101
119static Real gammaRatio_dl( Real const& a, Real const& x, bool lower_tail)
120{
121 Seriedl f(a, x);
122 if (lower_tail)
123 return (1.- a*sumAlternateSerie(f))*poisson_pdf_raw(a, x)*exp(x);
124 else
125 return apois(a, x) * sumAlternateSerie(f) * exp(x)
126 - expm1(a*log(x)-lgamma(a+1.));
127}
128
146static Real gammaRatio_cf( Real const& a, Real const& x, bool lower_tail)
147{
148 // initialize b_n
149 Real bn = x + 1. -a; // b_1
150 // initialize numerator
151 Real Nold = 0, Ncur = 1.; // = a_1 = 1
152 // initialize denominator
153 Real Dold = 1, Dcur = bn; // = b_1
154
155 // normalize if necessary by 2^32
156 if (std::abs(Dcur) > 4294967296.0)
157 {
158 Ncur /= Dcur;
159 Dold /= Dcur;
160 Dcur = 1.;
161 }
162
163 // auxiliary variables
164 Real aminusn = a;
165 // result
166 Real cf = Ncur/Dcur;
167 // iterations
168 for (int n=1; ; n++)
169 {
170 // compute a_n
171 Real an = n*(--aminusn);
172 // check trivial case
173 // (note a_n =0 => b_n != 0 for x>0 and thus cf!=N_n/D_n exists)
174 if (an == 0.) break;
175 // update b_n
176 bn += 2.;
177 // compute numerator and denominator
178 Real anew = bn * Ncur + an * Nold;
179 Nold = Ncur;
180 Ncur = anew;
181 anew = bn * Dcur + an * Dold;
182 Dold = Dcur;
183 Dcur = anew;
184 // normalize if necessary with 2^32
185 if (std::abs(Dcur) > 4294967296.0)
186 {
187 Ncur /= Dcur;
188 Nold /= Dcur;
189 Dold /= Dcur;
190 Dcur = 1.;
191 }
192 // normalize if necessary with 2^32
193 if (std::abs(Ncur) > 4294967296.0)
194 {
195 Ncur /= 4294967296.0;
196 Nold /= 4294967296.0;
197 Dold /= 4294967296.0;
198 Dcur /= 4294967296.0;
199 }
200 // test D_n not too small
201 if (std::abs(Dcur) != 0)
202 {
203 Real cfold = cf;
204 cf = Ncur/Dcur;
205 // check cv: relative; absolute for small cf
206 if (std::abs(cf - cfold) < std::abs(cf)*Arithmetic<Real>::epsilon())
207 break;
208 }
209 }
210 // cv
211 return (lower_tail) ? (1 - apois(a, x)*cf)
212 : (apois(a, x)*cf);
213}
214
229static Real gammaRatio_sr( Real const& a, Real const& x, bool lower_tail)
230{
231 Real y = a+1., term = x / (a+1.), sum = term;
232 /* sum = 1+\sum_{n=1}^\infty x^n / (a+1)*(a+2)*...*(a+n)) */
233 do
234 {
235 y++;
236 sum += (term *= x / y);
237 }
238 while (term > sum * Arithmetic<Real>::epsilon());
239 sum+=1;
240
241 if (lower_tail)
242 return poisson_pdf_raw(a, x)*sum;
243 else
244 return 1.-poisson_pdf_raw(a, x)*sum;
245}
246
267static Real gammaRatio_ae( Real const& a, Real const& x, bool lower_tail)
268{
269 Real term = 1, sum = 0, b=a-1;
270 // sum = 1+\sum_{n=1}^\infty (a-1)*...*(a - n) / x^n
271 while ((b > 1.) && (term > sum * Arithmetic<Real>::epsilon()))
272 {
273 sum += (term *= b / x);
274 b--;
275 }
276 sum+=1.;
277 // achieved cv
278 if (b>1.)
279 return lower_tail ? (1.-poisson_pdf_raw(a-1, x)*sum)
280 : ( poisson_pdf_raw(a-1, x)*sum);
281 // case (b>x) and (b <=1)
282 if (b>x)
283 return lower_tail ? (gammaRatio_sr(b, x, lower_tail)-poisson_pdf_raw(a-1, x)*sum)
284 : (gammaRatio_sr(b, x, lower_tail)+poisson_pdf_raw(a-1, x)*sum);
285 // case (b <x) and (b <=1)
286 return lower_tail ? gammaRatio_cf(b, x, lower_tail)-poisson_pdf_raw(a-1, x)*sum
287 : gammaRatio_cf(b, x, lower_tail)+poisson_pdf_raw(a-1, x)*sum;
288}
289
303static inline Real poisson_ae( Real const& a1, Real const& apd, bool lower_tail = true)
304{
305 // odd coefficients
306 static const Real coefs_i[8] =
307 {
308 0., /* placeholder used for 1-indexing */
309 2/3.,
310 -4/135.,
311 8/2835.,
312 16/8505.,
313 -8992/12629925.,
314 -334144/492567075.,
315 698752/1477701225.
316 };
317 // stirling coefficients
318 static const Real coefs_s[8] =
319 {
320 0., /* placeholder used for 1-indexing */
321 1/12.,
322 1/288.,
323 -139/51840.,
324 -571/2488320.,
325 163879/209018880.,
326 5246819/75246796800.,
327 -534703531/902961561600.
328 };
329 // compute D = -x(log(1+d/x) - d/x) = xlog(x/(x+d))+(x+d)-x
330 Real D = dev0(a1, apd);
331 // compute sqrt(2D)
332 Real sqrt2D = sqrt (2. * D);
333 // compute D/x
334 Real D_x = D/a1;
335
336 // treat negative difference
337 if (apd - a1 < 0) sqrt2D = -sqrt2D;
338
339 // variables for the numerator
340 Real num = 0;
341 Real sum1, num1_term, sum2, num2_term;
342 num1_term = sum1 = sqrt(a1);
343 num2_term = sum2 = sqrt2D;
344
345 // variables for the denominator
346 Real den = a1;
347 Real den_term = 1;
348 // computation
349 for (int i = 1; i < 8; i++)
350 {
351 // compute the numerator
352 num += num1_term * coefs_i[i];
353 num += num2_term * coefs_s[i];
354 // first sum
355 sum1 *= (D_x / i) ; // sqrt(x) * D^n/(n! x^n)
356 num1_term = num1_term/a1 + sum1;
357 // second sum
358 sum2 *= (D_x / (i + 0.5)); // sqrt(2D) * D^n/((1+0.5)*...(n+0.5)*x^n)
359 num2_term = num2_term/a1 + sum2;
360 // compute the denominator
361 den += den_term * coefs_s[i];
362 den_term /= a1;
363 }
364 // result (P or Q ?)
365 return lower_tail ? normal_cdf_raw( sqrt2D)
366 - (num / den) * normal_pdf_raw(sqrt2D)
367 : normal_cdf_raw(-sqrt2D)
368 + (num / den) * normal_pdf_raw(sqrt2D);
369}
370
371/* Compute the incomplete gamma function ratio P(a,x)
372 * \f[ P(a, x) = \frac{1}{\Gamma(a)}
373 * \int_0^x e^{-t} t^{a-1} dt
374 * \f]
375 * @param a parameter of the gamma ratio function
376 * @param x value to evaluate the gamma ratio function
377 * @param lower_tail @c true if we want the lower tail, @c false otherwise
378 **/
379Real gammaRatio(Real const& a, Real const& x, bool lower_tail)
380{
381 // Check if a and x are available
382 if (isNA(a)||isNA(x)) return a;
383 // Negative parameter not allowed
384 if (a<=0)
386 return gammaRatio_raw(a,x,lower_tail);
387}
388 // trivial case
389Real gammaRatio_raw(Real const& a, Real const& x, bool lower_tail)
390{
391 // trivial case
392 if (x<=0) return lower_tail ? 0 : 1.;
393 if (Arithmetic<Real>::isInfinite(x)) return lower_tail ? 1. : 0.;
394 // small values of x
395 if (x <= 1.) return gammaRatio_dl(a, x, lower_tail);
396 // large a compared to x
397 if (0.75*a> x) return gammaRatio_sr(a, x, lower_tail);
398 // large x compared to a
399 if (0.75*x > a)
400 { // if a>1 we can use the a.e. before computing
401 if (a>1) return gammaRatio_ae(a, x, lower_tail);
402 else return gammaRatio_cf(a, x, lower_tail);
403 }
404 // (a<100) and (x ~ a)
405 if (a<100.)
406 {
407 // a greater than x
408 if (a>x) return gammaRatio_sr(a, x, lower_tail);
409 else
410 { // if a>1 we can use the a.e. before computing
411 if (a>1) return gammaRatio_ae(a, x, lower_tail);
412 else return gammaRatio_cf(a, x, lower_tail);}
413 }
414 // (x ~ a) and (a>=100)
415 return poisson_ae(a-1, x, lower_tail);
416}
417} // namespace Funct
418
419} // namespace STK
420
421#endif /* IS_RTKPP_LIB */
const double a1
#define STKDOMAIN_ERROR_2ARG(Where, Arg1, Arg2, Error)
Definition STK_Macros.h:147
This Serie computes.
Real firstImpl() const
compute the first term of the serie
Real nextImpl() const
compute the next terms
Seriedl(Real const &a, Real const &x)
constructor
Interface base class for Series.
Definition STK_ISerie.h:67
The MultidimRegression class allows to regress a multidimensional output variable among a multivariat...
Real dev0(Real const &a, Real const &b)
compute the partial deviance .
Real expm1(Real const &x)
compute the function .
Real gammaRatio(Real const &a, Real const &x, bool lower_tail)
Compute the incomplete gamma functions ratio.
Real normal_pdf_raw(Real const &x)
compute the probability distribution function of the normal density.
Real sumAlternateSerie(const ISerie< Serie > &f, int const &n=0)
Sum an alternating series using the Chebichev polynomials.
Definition STK_Algo.h:65
Real poisson_pdf_raw(Real const &x, Real const &lambda)
Compute the Poisson density.
Real lgammaStirlingError(Real const &z)
Compute the error when we compute using the Stirling's formula.
Real normal_cdf_raw(Real const &x)
Compute the cumulative distribution function of the normal density.
Real lgamma(Real const &)
This function computes the function .
Real gammaRatioP(Real const &a, Real const &x)
Compute the incomplete gamma function ratio P(a,x).
bool isNA(Type const &x)
utility method allowing to know if a value is a NA (Not Available) value
Arrays::SumOp< Lhs, Rhs >::result_type sum(Lhs const &lhs, Rhs const &rhs)
convenience function for summing two arrays
double Real
STK fundamental type of Real values.
Real gammaRatio_raw(Real const &a, Real const &x, bool lower_tail)
The namespace STK is the main domain space of the Statistical ToolKit project.
Arithmetic properties of STK fundamental types.