TMB Documentation  v1.9.11
lgamma.cpp
1 /*
2  * Mathlib : A C Library of Special Functions
3  * Copyright (C) 1998 Ross Ihaka
4  * Copyright (C) 2000-2012 The R Core Team
5  *
6  * This program is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2 of the License, or
9  * (at your option) any later version.
10  *
11  * This program is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this program; if not, a copy is available at
18  * https://www.R-project.org/Licenses/
19  *
20  * SYNOPSIS
21  *
22  * #include <Rmath.h>
23  * double lgammafn_sign(double x, int *sgn);
24  * double lgammafn(double x);
25  *
26  * DESCRIPTION
27  *
28  * The function lgammafn computes log|gamma(x)|. The function
29  * lgammafn_sign in addition assigns the sign of the gamma function
30  * to the address in the second argument if this is not NULL.
31  *
32  * NOTES
33  *
34  * This routine is a translation into C of a Fortran subroutine
35  * by W. Fullerton of Los Alamos Scientific Laboratory.
36  *
37  * The accuracy of this routine compares (very) favourably
38  * with those of the Sun Microsystems portable mathematical
39  * library.
40  */
41 
42 template<class Float>
43 Float lgammafn_sign(Float x, int *sgn)
44 {
45  Float ans, y, sinpiy;
46 
47 #ifdef NOMORE_FOR_THREADS
48  static double xmax = 0.;
49  static double dxrel = 0.;
50 
51  if (xmax == 0) {/* initialize machine dependent constants _ONCE_ */
52  xmax = d1mach(2)/log(d1mach(2));/* = 2.533 e305 for IEEE double */
53  dxrel = sqrt (d1mach(4));/* sqrt(Eps) ~ 1.49 e-8 for IEEE double */
54  }
55 #else
56 /* For IEEE double precision DBL_EPSILON = 2^-52 = 2.220446049250313e-16 :
57  xmax = DBL_MAX / log(DBL_MAX) = 2^1024 / (1024 * log(2)) = 2^1014 / log(2)
58  dxrel = sqrt(DBL_EPSILON) = 2^-26 = 5^26 * 1e-26 (is *exact* below !)
59  */
60 #define xmax 2.5327372760800758e+305
61 #define dxrel 1.490116119384765625e-8
62 #endif
63 
64  if (sgn != NULL) *sgn = 1;
65 
66 #ifdef IEEE_754
67  if(ISNAN(x)) return x;
68 #endif
69 
70  if (sgn != NULL && x < 0 && fmod(floor(-x), 2.) == 0)
71  *sgn = -1;
72 
73  if (x <= 0 && x == trunc(x)) { /* Negative integer argument */
74  ML_ERROR(ME_RANGE, "lgamma");
75  return ML_POSINF;/* +Inf, since lgamma(x) = log|gamma(x)| */
76  }
77 
78  y = fabs(x);
79 
80  if (y < 1e-306) return -log(y); // denormalized range, R change
81  if (y <= 10) return log(fabs(gammafn(x)));
82  /*
83  ELSE y = |x| > 10 ---------------------- */
84 
85  if (y > xmax) {
86  ML_ERROR(ME_RANGE, "lgamma");
87  return ML_POSINF;
88  }
89 
90  if (x > 0) { /* i.e. y = x > 10 */
91 #ifdef IEEE_754
92  if(x > 1e17)
93  return(x*(log(x) - 1.));
94  else if(x > 4934720.)
95  return(M_LN_SQRT_2PI + (x - 0.5) * log(x) - x);
96  else
97 #endif
98  return M_LN_SQRT_2PI + (x - 0.5) * log(x) - x + lgammacor(x);
99  }
100  /* else: x < -10; y = -x */
101  sinpiy = fabs(sinpi(y));
102 
103  if (sinpiy == 0) { /* Negative integer argument ===
104  Now UNNECESSARY: caught above */
105  MATHLIB_WARNING(" ** should NEVER happen! *** [lgamma.c: Neg.int, y=%g]\n",y);
106  ML_ERR_return_NAN;
107  }
108 
109  ans = M_LN_SQRT_PId2 + (x - 0.5) * log(y) - x - log(sinpiy) - lgammacor(y);
110 
111  if(fabs((x - trunc(x - 0.5)) * ans / x) < dxrel) {
112 
113  /* The answer is less than half precision because
114  * the argument is too near a negative integer. */
115 
116  ML_ERROR(ME_PRECISION, "lgamma");
117  }
118 
119  return ans;
120 }
121 
122 template<class Float>
123 Float lgammafn(Float x)
124 {
125  return lgammafn_sign(x, NULL);
126 }
License: GPL v2