TMB Documentation  v1.9.11
TMB/inst/include/tiny_ad/gamma/gamma.cpp
1 /*
2  * Mathlib : A C Library of Special Functions
3  * Copyright (C) 1998 Ross Ihaka
4  * Copyright (C) 2000-2013 The R Core Team
5  * Copyright (C) 2002-2004 The R Foundation
6  *
7  * This program is free software; you can redistribute it and/or modify
8  * it under the terms of the GNU General Public License as published by
9  * the Free Software Foundation; either version 2 of the License, or
10  * (at your option) any later version.
11  *
12  * This program is distributed in the hope that it will be useful,
13  * but WITHOUT ANY WARRANTY; without even the implied warranty of
14  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15  * GNU General Public License for more details.
16  *
17  * You should have received a copy of the GNU General Public License
18  * along with this program; if not, a copy is available at
19  * https://www.R-project.org/Licenses/
20  *
21  * SYNOPSIS
22  *
23  * #include <Rmath.h>
24  * double gammafn(double x);
25  *
26  * DESCRIPTION
27  *
28  * This function computes the value of the gamma function.
29  *
30  * NOTES
31  *
32  * This function is a translation into C of a Fortran subroutine
33  * by W. Fullerton of Los Alamos Scientific Laboratory.
34  * (e.g. http://www.netlib.org/slatec/fnlib/gamma.f)
35  *
36  * The accuracy of this routine compares (very) favourably
37  * with those of the Sun Microsystems portable mathematical
38  * library.
39  *
40  * MM specialized the case of n! for n < 50 - for even better precision
41  */
42 
43 template<class Float>
44 Float gammafn(Float x)
45 {
46  const static double gamcs[42] = {
47  +.8571195590989331421920062399942e-2,
48  +.4415381324841006757191315771652e-2,
49  +.5685043681599363378632664588789e-1,
50  -.4219835396418560501012500186624e-2,
51  +.1326808181212460220584006796352e-2,
52  -.1893024529798880432523947023886e-3,
53  +.3606925327441245256578082217225e-4,
54  -.6056761904460864218485548290365e-5,
55  +.1055829546302283344731823509093e-5,
56  -.1811967365542384048291855891166e-6,
57  +.3117724964715322277790254593169e-7,
58  -.5354219639019687140874081024347e-8,
59  +.9193275519859588946887786825940e-9,
60  -.1577941280288339761767423273953e-9,
61  +.2707980622934954543266540433089e-10,
62  -.4646818653825730144081661058933e-11,
63  +.7973350192007419656460767175359e-12,
64  -.1368078209830916025799499172309e-12,
65  +.2347319486563800657233471771688e-13,
66  -.4027432614949066932766570534699e-14,
67  +.6910051747372100912138336975257e-15,
68  -.1185584500221992907052387126192e-15,
69  +.2034148542496373955201026051932e-16,
70  -.3490054341717405849274012949108e-17,
71  +.5987993856485305567135051066026e-18,
72  -.1027378057872228074490069778431e-18,
73  +.1762702816060529824942759660748e-19,
74  -.3024320653735306260958772112042e-20,
75  +.5188914660218397839717833550506e-21,
76  -.8902770842456576692449251601066e-22,
77  +.1527474068493342602274596891306e-22,
78  -.2620731256187362900257328332799e-23,
79  +.4496464047830538670331046570666e-24,
80  -.7714712731336877911703901525333e-25,
81  +.1323635453126044036486572714666e-25,
82  -.2270999412942928816702313813333e-26,
83  +.3896418998003991449320816639999e-27,
84  -.6685198115125953327792127999999e-28,
85  +.1146998663140024384347613866666e-28,
86  -.1967938586345134677295103999999e-29,
87  +.3376448816585338090334890666666e-30,
88  -.5793070335782135784625493333333e-31
89  };
90 
91  int i, n;
92  Float y;
93  Float sinpiy, value;
94 
95 #ifdef NOMORE_FOR_THREADS
96  static int ngam = 0;
97  static double xmin = 0, xmax = 0., xsml = 0., dxrel = 0.;
98 
99  /* Initialize machine dependent constants, the first time gamma() is called.
100  FIXME for threads ! */
101  if (ngam == 0) {
102  ngam = chebyshev_init(gamcs, 42, DBL_EPSILON/20);/*was .1*d1mach(3)*/
103  gammalims(&xmin, &xmax);/*-> ./gammalims.c */
104  xsml = exp(fmax2(log(DBL_MIN), -log(DBL_MAX)) + 0.01);
105  /* = exp(.01)*DBL_MIN = 2.247e-308 for IEEE */
106  dxrel = sqrt(DBL_EPSILON);/*was sqrt(d1mach(4)) */
107  }
108 #else
109 /* For IEEE double precision DBL_EPSILON = 2^-52 = 2.220446049250313e-16 :
110  * (xmin, xmax) are non-trivial, see ./gammalims.c
111  * xsml = exp(.01)*DBL_MIN
112  * dxrel = sqrt(DBL_EPSILON) = 2 ^ -26
113 */
114 # define ngam 22
115 # define xmin -170.5674972726612
116 # define xmax 171.61447887182298
117 # define xsml 2.2474362225598545e-308
118 # define dxrel 1.490116119384765696e-8
119 #endif
120 
121  if(ISNAN(x)) return x;
122 
123  /* If the argument is exactly zero or a negative integer
124  * then return NaN. */
125  if (x == 0 || (x < 0 && x == round(x))) {
126  ML_ERROR(ME_DOMAIN, "gammafn");
127  return ML_NAN;
128  }
129 
130  y = fabs(x);
131 
132  if (y <= 10) {
133 
134  /* Compute gamma(x) for -10 <= x <= 10
135  * Reduce the interval and find gamma(1 + y) for 0 <= y < 1
136  * first of all. */
137 
138  n = (int) trunc(x);
139  if(x < 0) --n;
140  y = x - n;/* n = floor(x) ==> y in [ 0, 1 ) */
141  --n;
142  value = chebyshev_eval(y * 2 - 1, gamcs, ngam) + .9375;
143  if (n == 0)
144  return value;/* x = 1.dddd = 1+y */
145 
146  if (n < 0) {
147  /* compute gamma(x) for -10 <= x < 1 */
148 
149  /* exact 0 or "-n" checked already above */
150 
151  /* The answer is less than half precision */
152  /* because x too near a negative integer. */
153  if (x < -0.5 && fabs(x - (int)trunc(x - 0.5) / x) < dxrel) {
154  ML_ERROR(ME_PRECISION, "gammafn");
155  }
156 
157  /* The argument is so close to 0 that the result would overflow. */
158  if (y < xsml) {
159  ML_ERROR(ME_RANGE, "gammafn");
160  if(x > 0) return ML_POSINF;
161  else return ML_NEGINF;
162  }
163 
164  n = -n;
165 
166  for (i = 0; i < n; i++) {
167  value /= (x + i);
168  }
169  return value;
170  }
171  else {
172  /* gamma(x) for 2 <= x <= 10 */
173 
174  for (i = 1; i <= n; i++) {
175  value *= (y + i);
176  }
177  return value;
178  }
179  }
180  else {
181  /* gamma(x) for y = |x| > 10. */
182 
183  if (x > xmax) { /* Overflow */
184  ML_ERROR(ME_RANGE, "gammafn");
185  return ML_POSINF;
186  }
187 
188  if (x < xmin) { /* Underflow */
189  ML_ERROR(ME_UNDERFLOW, "gammafn");
190  return 0.;
191  }
192 
193  if(y <= 50 && y == (int)trunc(y)) { /* compute (n - 1)! */
194  value = 1.;
195  for (i = 2; i < y; i++) value *= i;
196  }
197  else { /* normal case */
198  value = exp((y - 0.5) * log(y) - y + M_LN_SQRT_2PI +
199  ((2*y == (int)trunc(2*y))? stirlerr(y) : lgammacor(y)));
200  }
201  if (x > 0)
202  return value;
203 
204  if (fabs((x - (int)trunc(x - 0.5))/x) < dxrel){
205 
206  /* The answer is less than half precision because */
207  /* the argument is too near a negative integer. */
208 
209  ML_ERROR(ME_PRECISION, "gammafn");
210  }
211 
212  sinpiy = sinpi(y);
213  if (sinpiy == 0) { /* Negative integer arg - overflow */
214  ML_ERROR(ME_RANGE, "gammafn");
215  return ML_POSINF;
216  }
217 
218  return -M_PI / (y * sinpiy * value);
219  }
220 }
License: GPL v2