Gamma Function

The Gamma function is defined by the following integral,

\Gamma(x) = \int_0^t dt  t^{x-1} \exp(-t)

Further information on the Gamma function can be found in Abramowitz & Stegun, Chapter 6. The functions described in this section are declared in the header file `gsl_sf_gamma.h'.

Function: double gsl_sf_gamma (double x)
Function: int gsl_sf_gamma_e (double x, gsl_sf_result * result)
These routines compute the Gamma function \Gamma(x), subject to x not being a negative integer. The function is computed using the real Lanczos method. The maximum value of x such that \Gamma(x) is not considered an overflow is given by the macro GSL_SF_GAMMA_XMAX and is 171.0.

Function: double gsl_sf_lngamma (double x)
Function: int gsl_sf_lngamma_e (double x, gsl_sf_result * result)
These routines compute the logarithm of the Gamma function, \log(\Gamma(x)), subject to x not a being negative integer. For x<0 the real part of \log(\Gamma(x)) is returned, which is equivalent to \log(|\Gamma(x)|). The function is computed using the real Lanczos method.

Function: int gsl_sf_lngamma_sgn_e (double x, gsl_sf_result * result_lg, double * sgn)
This routine computes the sign of the gamma function and the logarithm its magnitude, subject to x not being a negative integer. The function is computed using the real Lanczos method. The value of the gamma function can be reconstructed using the relation \Gamma(x) = sgn * \exp(resultlg).

Function: double gsl_sf_gammastar (double x)
Function: int gsl_sf_gammastar_e (double x, gsl_sf_result * result)
These routines compute the regulated Gamma Function \Gamma^*(x) for x > 0. The regulated gamma function is given by,

\Gamma^*(x) = \Gamma(x)/(\sqrt{2\pi} x^{(x-1/2)} \exp(-x))
            = (1 + (1/12x) + ...)  for x \to \infty

and is a useful suggestion of Temme.

Function: double gsl_sf_gammainv (double x)
Function: int gsl_sf_gammainv_e (double x, gsl_sf_result * result)
These routines compute the reciprocal of the gamma function, 1/\Gamma(x) using the real Lanczos method.

Function: int gsl_sf_lngamma_complex_e (double zr, double zi, gsl_sf_result * lnr, gsl_sf_result * arg)
This routine computes \log(\Gamma(z)) for complex z=z_r+i z_i and z not a negative integer, using the complex Lanczos method. The returned parameters are lnr = \log|\Gamma(z)| and arg = \arg(\Gamma(z)) in (-\pi,\pi]. Note that the phase part (arg) is not well-determined when |z| is very large, due to inevitable roundoff in restricting to (-\pi,\pi]. This will result in a GSL_ELOSS error when it occurs. The absolute value part (lnr), however, never suffers from loss of precision.

Function: double gsl_sf_taylorcoeff (int n, double x)
Function: int gsl_sf_taylorcoeff_e (int n, double x, gsl_sf_result * result)
These routines compute the Taylor coefficient x^n / n! for x >= 0, n >= 0.

Function: double gsl_sf_fact (unsigned int n)
Function: int gsl_sf_fact_e (unsigned int n, gsl_sf_result * result)
These routines compute the factorial n!. The factorial is related to the Gamma function by n! = \Gamma(n+1).

Function: double gsl_sf_doublefact (unsigned int n)
Function: int gsl_sf_doublefact_e (unsigned int n, gsl_sf_result * result)
These routines compute the double factorial n!! = n(n-2)(n-4) \dots.

Function: double gsl_sf_lnfact (unsigned int n)
Function: int gsl_sf_lnfact_e (unsigned int n, gsl_sf_result * result)
These routines compute the logarithm of the factorial of n, \log(n!). The algorithm is faster than computing \ln(\Gamma(n+1)) via gsl_sf_lngamma for n < 170, but defers for larger n.

Function: double gsl_sf_lndoublefact (unsigned int n)
Function: int gsl_sf_lndoublefact_e (unsigned int n, gsl_sf_result * result)
These routines compute the logarithm of the double factorial of n, \log(n!!).

Function: double gsl_sf_choose (unsigned int n, unsigned int m)
Function: int gsl_sf_choose_e (unsigned int n, unsigned int m, gsl_sf_result * result)
These routines compute the combinatorial factor n choose m = n!/(m!(n-m)!)

Function: double gsl_sf_lnchoose (unsigned int n, unsigned int m)
Function: int gsl_sf_lnchoose_e (unsigned int n, unsigned int m, gsl_sf_result * result)
These routines compute the logarithm of n choose m. This is equivalent to the sum \log(n!) - \log(m!) - \log((n-m)!).

Function: double gsl_sf_poch (double a, double x)
Function: int gsl_sf_poch_e (double a, double x, gsl_sf_result * result)
These routines compute the Pochhammer symbol (a)_x := \Gamma(a + x)/\Gamma(x), subject to a and a+x not being negative integers. The Pochhammer symbol is also known as the Apell symbol.

Function: double gsl_sf_lnpoch (double a, double x)
Function: int gsl_sf_lnpoch_e (double a, double x, gsl_sf_result * result)
These routines compute the logarithm of the Pochhammer symbol, \log((a)_x) = \log(\Gamma(a + x)/\Gamma(a)) for a > 0, a+x > 0.

Function: int gsl_sf_lnpoch_sgn_e (double a, double x, gsl_sf_result * result, double * sgn)
These routines compute the sign of the Pochhammer symbol and the logarithm of its magnitude. The computed parameters are result = \log(|(a)_x|) and sgn = sgn((a)_x) where (a)_x := \Gamma(a + x)/\Gamma(a), subject to a, a+x not being negative integers.

Function: double gsl_sf_pochrel (double a, double x)
Function: int gsl_sf_pochrel_e (double a, double x, gsl_sf_result * result)
These routines compute the relative Pochhammer symbol ((a,x) - 1)/x where (a,x) = (a)_x := \Gamma(a + x)/\Gamma(a).

Function: double gsl_sf_gamma_inc_Q (double a, double x)
Function: int gsl_sf_gamma_inc_Q_e (double a, double x, gsl_sf_result * result)
These routines compute the normalized incomplete Gamma Function Q(a,x) = 1/\Gamma(a) \int_x\infty dt t^{a-1} \exp(-t) for a > 0, x >= 0.

Function: double gsl_sf_gamma_inc_P (double a, double x)
Function: int gsl_sf_gamma_inc_P_e (double a, double x, gsl_sf_result * result)
These routines compute the complementary normalized incomplete Gamma Function P(a,x) = 1/\Gamma(a) \int_0^x dt t^{a-1} \exp(-t) for a > 0, x >= 0.

Note that Abramowitz & Stegun call P(a,x) the incomplete gamma function (section 6.5).

Function: double gsl_sf_beta (double a, double b)
Function: int gsl_sf_beta_e (double a, double b, gsl_sf_result * result)
These routines compute the Beta Function, B(a,b) = \Gamma(a)\Gamma(b)/\Gamma(a+b) for a > 0, b > 0.

Function: double gsl_sf_lnbeta (double a, double b)
Function: int gsl_sf_lnbeta_e (double a, double b, gsl_sf_result * result)
These routines compute the logarithm of the Beta Function, \log(B(a,b)) for a > 0, b > 0.

Function: double gsl_sf_beta_inc (double a, double b, double x)
Function: int gsl_sf_beta_inc_e (double a, double b, double x, gsl_sf_result * result)
These routines compute the normalize incomplete Beta function B_x(a,b)/B(a,b) for a > 0, b > 0, and 0 <= x <= 1.