TMB Documentation  v1.9.11
lgammacor.cpp
1 /*
2  * Mathlib : A C Library of Special Functions
3  * Copyright (C) 1998 Ross Ihaka
4  * Copyright (C) 2000-2001 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 lgammacor(double x);
24  *
25  * DESCRIPTION
26  *
27  * Compute the log gamma correction factor for x >= 10 so that
28  *
29  * log(gamma(x)) = .5*log(2*pi) + (x-.5)*log(x) -x + lgammacor(x)
30  *
31  * [ lgammacor(x) is called Del(x) in other contexts (e.g. dcdflib)]
32  *
33  * NOTES
34  *
35  * This routine is a translation into C of a Fortran subroutine
36  * written by W. Fullerton of Los Alamos Scientific Laboratory.
37  *
38  * SEE ALSO
39  *
40  * Loader(1999)'s stirlerr() {in ./stirlerr.c} is *very* similar in spirit,
41  * is faster and cleaner, but is only defined "fast" for half integers.
42  */
43 
44 template<class Float>
45 Float attribute_hidden lgammacor(Float x)
46 {
47  const static double algmcs[15] = {
48  +.1666389480451863247205729650822e+0,
49  -.1384948176067563840732986059135e-4,
50  +.9810825646924729426157171547487e-8,
51  -.1809129475572494194263306266719e-10,
52  +.6221098041892605227126015543416e-13,
53  -.3399615005417721944303330599666e-15,
54  +.2683181998482698748957538846666e-17,
55  -.2868042435334643284144622399999e-19,
56  +.3962837061046434803679306666666e-21,
57  -.6831888753985766870111999999999e-23,
58  +.1429227355942498147573333333333e-24,
59  -.3547598158101070547199999999999e-26,
60  +.1025680058010470912000000000000e-27,
61  -.3401102254316748799999999999999e-29,
62  +.1276642195630062933333333333333e-30
63  };
64 
65  Float tmp;
66 
67 /* For IEEE double precision DBL_EPSILON = 2^-52 = 2.220446049250313e-16 :
68  * xbig = 2 ^ 26.5
69  * xmax = DBL_MAX / 48 = 2^1020 / 3 */
70 #define nalgm 5
71 #define xbig 94906265.62425156
72 #define xmax 3.745194030963158e306
73 
74  if (x < 10)
75  ML_ERR_return_NAN
76  else if (x >= xmax) {
77  ML_ERROR(ME_UNDERFLOW, "lgammacor");
78  /* allow to underflow below */
79  }
80  else if (x < xbig) {
81  tmp = 10 / x;
82  return chebyshev_eval(tmp * tmp * 2 - 1, algmcs, nalgm) / x;
83  }
84  return 1 / (x * 12);
85 }
License: GPL v2