| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
This chapter describes the GSL special function library. The library includes routines for calculating the values of Airy functions, Bessel functions, Clausen functions, Coulomb wave functions, Coupling coefficients, the Dawson function, Debye functions, Dilogarithms, Elliptic integrals, Jacobi elliptic functions, Error functions, Exponential integrals, Fermi-Dirac functions, Gamma functions, Gegenbauer functions, Hypergeometric functions, Laguerre functions, Legendre functions and Spherical Harmonics, the Psi (Digamma) Function, Synchrotron functions, Transport functions, Trigonometric functions and Zeta functions. Each routine also computes an estimate of the numerical error in the calculated value of the function.
The functions are declared in individual header files, such as `gsl_sf_airy.h', `gsl_sf_bessel.h', etc. The complete set of header files can be included using the file `gsl_sf.h'.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The special functions are available in two calling conventions, a natural form which returns the numerical value of the function and an error-handling form which returns an error code. The two types of function provide alternative ways of accessing the same underlying code.
The natural form returns only the value of the function and can be used directly in mathematical expressions. For example, the following function call will compute the value of the Bessel function J_0(x),
double y = gsl_sf_bessel_J0 (x); |
There is no way to access an error code or to estimate the error using this method. To allow access to this information the alternative error-handling form stores the value and error in a modifiable argument,
gsl_sf_result result; int status = gsl_sf_bessel_J0_e (x, &result); |
The error-handling functions have the suffix _e. The returned
status value indicates error conditions such as overflow, underflow or
loss of precision. If there are no errors the error-handling functions
return GSL_SUCCESS.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The error handling form of the special functions always calculate an error estimate along with the value of the result. Therefore, structures are provided for amalgamating a value and error estimate. These structures are declared in the header file `gsl_sf_result.h'.
The gsl_sf_result struct contains value and error fields.
typedef struct
{
double val;
double err;
} gsl_sf_result;
|
The field val contains the value and the field err contains an estimate of the absolute error in the value.
In some cases, an overflow or underflow can be detected and handled by a
function. In this case, it may be possible to return a scaling exponent
as well as an error/value pair in order to save the result from
exceeding the dynamic range of the built-in types. The
gsl_sf_result_e10 struct contains value and error fields as well
as an exponent field such that the actual result is obtained as
result * 10^(e10).
typedef struct
{
double val;
double err;
int e10;
} gsl_sf_result_e10;
|
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The goal of the library is to achieve double precision accuracy wherever
possible. However the cost of evaluating some special functions to
double precision can be significant, particularly where very high order
terms are required. In these cases a mode argument allows the
accuracy of the function to be reduced in order to improve performance.
The following precision levels are available for the mode argument,
GSL_PREC_DOUBLEDouble-precision, a relative accuracy of approximately 2 * 10^-16.
GSL_PREC_SINGLESingle-precision, a relative accuracy of approximately 10^-7.
GSL_PREC_APPROXApproximate values, a relative accuracy of approximately 5 * 10^-4.
The approximate mode provides the fastest evaluation at the lowest accuracy.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Airy functions Ai(x) and Bi(x) are defined by the integral representations,
For further information see Abramowitz & Stegun, Section 10.4. The Airy functions are defined in the header file `gsl_sf_airy.h'.
| 7.4.1 Airy Functions | ||
| 7.4.2 Derivatives of Airy Functions | ||
| 7.4.3 Zeros of Airy Functions | ||
| 7.4.4 Zeros of Derivatives of Airy Functions |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Airy function Ai(x) with an accuracy specified by mode.
These routines compute the Airy function Bi(x) with an accuracy specified by mode.
These routines compute a scaled version of the Airy function S_A(x) Ai(x). For x>0 the scaling factor S_A(x) is \exp(+(2/3) x^(3/2)), and is 1 for x<0.
These routines compute a scaled version of the Airy function S_B(x) Bi(x). For x>0 the scaling factor S_B(x) is exp(-(2/3) x^(3/2)), and is 1 for x<0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Airy function derivative Ai'(x) with an accuracy specified by mode.
These routines compute the Airy function derivative Bi'(x) with an accuracy specified by mode.
These routines compute the derivative of the scaled Airy function S_A(x) Ai(x).
These routines compute the derivative of the scaled Airy function S_B(x) Bi(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the location of the s-th zero of the Airy function Ai(x).
These routines compute the location of the s-th zero of the Airy function Bi(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the location of the s-th zero of the Airy function derivative Ai'(x).
These routines compute the location of the s-th zero of the Airy function derivative Bi'(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The routines described in this section compute the Cylindrical Bessel functions J_n(x), Y_n(x), Modified cylindrical Bessel functions I_n(x), K_n(x), Spherical Bessel functions j_l(x), y_l(x), and Modified Spherical Bessel functions i_l(x), k_l(x). For more information see Abramowitz & Stegun, Chapters 9 and 10. The Bessel functions are defined in the header file `gsl_sf_bessel.h'.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the regular cylindrical Bessel function of zeroth order, J_0(x).
These routines compute the regular cylindrical Bessel function of first order, J_1(x).
These routines compute the regular cylindrical Bessel function of order n, J_n(x).
This routine computes the values of the regular cylindrical Bessel functions J_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the irregular cylindrical Bessel function of zeroth order, Y_0(x), for x>0.
These routines compute the irregular cylindrical Bessel function of first order, Y_1(x), for x>0.
These routines compute the irregular cylindrical Bessel function of order n, Y_n(x), for x>0.
This routine computes the values of the irregular cylindrical Bessel functions Y_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The domain of the function is x>0. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the regular modified cylindrical Bessel function of zeroth order, I_0(x).
These routines compute the regular modified cylindrical Bessel function of first order, I_1(x).
These routines compute the regular modified cylindrical Bessel function of order n, I_n(x).
This routine computes the values of the regular modified cylindrical Bessel functions I_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The start of the range nmin must be positive or zero. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
These routines compute the scaled regular modified cylindrical Bessel function of zeroth order \exp(-|x|) I_0(x).
These routines compute the scaled regular modified cylindrical Bessel function of first order \exp(-|x|) I_1(x).
These routines compute the scaled regular modified cylindrical Bessel function of order n, \exp(-|x|) I_n(x)
This routine computes the values of the scaled regular cylindrical Bessel functions \exp(-|x|) I_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The start of the range nmin must be positive or zero. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the irregular modified cylindrical Bessel function of zeroth order, K_0(x), for x > 0.
These routines compute the irregular modified cylindrical Bessel function of first order, K_1(x), for x > 0.
These routines compute the irregular modified cylindrical Bessel function of order n, K_n(x), for x > 0.
This routine computes the values of the irregular modified cylindrical Bessel functions K_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The start of the range nmin must be positive or zero. The domain of the function is x>0. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
These routines compute the scaled irregular modified cylindrical Bessel function of zeroth order \exp(x) K_0(x) for x>0.
These routines compute the scaled irregular modified cylindrical Bessel function of first order \exp(x) K_1(x) for x>0.
These routines compute the scaled irregular modified cylindrical Bessel function of order n, \exp(x) K_n(x), for x>0.
This routine computes the values of the scaled irregular cylindrical Bessel functions \exp(x) K_n(x) for n from nmin to nmax inclusive, storing the results in the array result_array. The start of the range nmin must be positive or zero. The domain of the function is x>0. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the regular spherical Bessel function of zeroth order, j_0(x) = \sin(x)/x.
These routines compute the regular spherical Bessel function of first order, j_1(x) = (\sin(x)/x - \cos(x))/x.
These routines compute the regular spherical Bessel function of second order, j_2(x) = ((3/x^2 - 1)\sin(x) - 3\cos(x)/x)/x.
These routines compute the regular spherical Bessel function of order l, j_l(x), for l >= 0 and x >= 0.
This routine computes the values of the regular spherical Bessel functions j_l(x) for l from 0 to lmax inclusive for lmax >= 0 and x >= 0, storing the results in the array result_array. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
This routine uses Steed's method to compute the values of the regular spherical Bessel functions j_l(x) for l from 0 to lmax inclusive for lmax >= 0 and x >= 0, storing the results in the array result_array. The Steed/Barnett algorithm is described in Comp. Phys. Comm. 21, 297 (1981). Steed's method is more stable than the recurrence used in the other functions but is also slower.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the irregular spherical Bessel function of zeroth order, y_0(x) = -\cos(x)/x.
These routines compute the irregular spherical Bessel function of first order, y_1(x) = -(\cos(x)/x + \sin(x))/x.
These routines compute the irregular spherical Bessel function of second order, y_2(x) = (-3/x^3 + 1/x)\cos(x) - (3/x^2)\sin(x).
These routines compute the irregular spherical Bessel function of order l, y_l(x), for l >= 0.
This routine computes the values of the irregular spherical Bessel functions y_l(x) for l from 0 to lmax inclusive for lmax >= 0, storing the results in the array result_array. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The regular modified spherical Bessel functions i_l(x) are related to the modified Bessel functions of fractional order, i_l(x) = \sqrt{\pi/(2x)} I_{l+1/2}(x)
These routines compute the scaled regular modified spherical Bessel function of zeroth order, \exp(-|x|) i_0(x).
These routines compute the scaled regular modified spherical Bessel function of first order, \exp(-|x|) i_1(x).
These routines compute the scaled regular modified spherical Bessel function of second order, \exp(-|x|) i_2(x)
These routines compute the scaled regular modified spherical Bessel function of order l, \exp(-|x|) i_l(x)
This routine computes the values of the scaled regular modified cylindrical Bessel functions \exp(-|x|) i_l(x) for l from 0 to lmax inclusive for lmax >= 0, storing the results in the array result_array. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The irregular modified spherical Bessel functions k_l(x) are related to the irregular modified Bessel functions of fractional order, k_l(x) = \sqrt{\pi/(2x)} K_{l+1/2}(x).
These routines compute the scaled irregular modified spherical Bessel function of zeroth order, \exp(x) k_0(x), for x>0.
These routines compute the scaled irregular modified spherical Bessel function of first order, \exp(x) k_1(x), for x>0.
These routines compute the scaled irregular modified spherical Bessel function of second order, \exp(x) k_2(x), for x>0.
These routines compute the scaled irregular modified spherical Bessel function of order l, \exp(x) k_l(x), for x>0.
This routine computes the values of the scaled irregular modified spherical Bessel functions \exp(x) k_l(x) for l from 0 to lmax inclusive for lmax >= 0 and x>0, storing the results in the array result_array. The values are computed using recurrence relations, for efficiency, and therefore may differ slightly from the exact values.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the regular cylindrical Bessel function of fractional order nu, J_\nu(x).
This function computes the regular cylindrical Bessel function of fractional order \nu, J_\nu(x), evaluated at a series of x values. The array v of length size contains the x values. They are assumed to be strictly ordered and positive. The array is over-written with the values of J_\nu(x_i).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the irregular cylindrical Bessel function of fractional order nu, Y_\nu(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the regular modified Bessel function of fractional order nu, I_\nu(x) for x>0, \nu>0.
These routines compute the scaled regular modified Bessel function of fractional order nu, \exp(-|x|)I_\nu(x) for x>0, \nu>0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the irregular modified Bessel function of fractional order nu, K_\nu(x) for x>0, \nu>0.
These routines compute the logarithm of the irregular modified Bessel function of fractional order nu, \ln(K_\nu(x)) for x>0, \nu>0.
These routines compute the scaled irregular modified Bessel function of fractional order nu, \exp(+|x|) K_\nu(x) for x>0, \nu>0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the location of the s-th positive zero of the Bessel function J_0(x).
These routines compute the location of the s-th positive zero of the Bessel function J_1(x).
These routines compute the location of the s-th positive zero of the Bessel function J_\nu(x). The current implementation does not support negative values of nu.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Clausen function is defined by the following integral, It is related to the dilogarithm by Cl_2(\theta) = \Im Li_2(\exp(i \theta)). The Clausen functions are declared in the header file `gsl_sf_clausen.h'.
These routines compute the Clausen integral Cl_2(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Coulomb functions are declared in the header file `gsl_sf_coulomb.h'. Both bound state and scattering solutions are available.
| 7.7.1 Normalized Hydrogenic Bound States | ||
| 7.7.2 Coulomb Wave Functions | ||
| 7.7.3 Coulomb Wave Function Normalization Constant |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the lowest-order normalized hydrogenic bound state radial wavefunction R_1 := 2Z \sqrt{Z} \exp(-Z r).
These routines compute the n-th normalized hydrogenic bound state radial wavefunction, The normalization is chosen such that the wavefunction \psi is given by \psi(n,l,r) = R_n Y_{lm}.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Coulomb wave functions F_L(\eta,x), G_L(\eta,x) are
described in Abramowitz & Stegun, Chapter 14. Because there can be a
large dynamic range of values for these functions, overflows are handled
gracefully. If an overflow occurs, GSL_EOVRFLW is signalled and
exponent(s) are returned through the modifiable parameters exp_F,
exp_G. The full solution can be reconstructed from the following
relations,
This function computes the coulomb wave functions F_L(\eta,x),
G_{L-k}(\eta,x) and their derivatives with respect to x,
F'_L(\eta,x)
G'_{L-k}(\eta,x). The parameters are restricted to L,
L-k > -1/2, x > 0 and integer k. Note that L
itself is not restricted to being an integer. The results are stored in
the parameters F, G for the function values and Fp,
Gp for the derivative values. If an overflow occurs,
GSL_EOVRFLW is returned and scaling exponents are stored in
the modifiable parameters exp_F, exp_G.
This function computes the function F_L(eta,x) for L = Lmin \dots Lmin + kmax storing the results in fc_array. In the case of overflow the exponent is stored in F_exponent.
This function computes the functions F_L(\eta,x), G_L(\eta,x) for L = Lmin \dots Lmin + kmax storing the results in fc_array and gc_array. In the case of overflow the exponents are stored in F_exponent and G_exponent.
This function computes the functions F_L(\eta,x), G_L(\eta,x) and their derivatives F'_L(\eta,x), G'_L(\eta,x) for L = Lmin \dots Lmin + kmax storing the results in fc_array, gc_array, fcp_array and gcp_array. In the case of overflow the exponents are stored in F_exponent and G_exponent.
This function computes the Coulomb wave function divided by the argument F_L(\eta, x)/x for L = Lmin \dots Lmin + kmax, storing the results in fc_array. In the case of overflow the exponent is stored in F_exponent. This function reduces to spherical Bessel functions in the limit \eta \to 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Coulomb wave function normalization constant is defined in Abramowitz 14.1.7.
This function computes the Coulomb wave function normalization constant C_L(\eta) for L > -1.
This function computes the coulomb wave function normalization constant C_L(\eta) for L = Lmin \dots Lmin + kmax, Lmin > -1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Wigner 3-j, 6-j and 9-j symbols give the coupling coefficients for combined angular momentum vectors. Since the arguments of the standard coupling coefficient functions are integer or half-integer, the arguments of the following functions are, by convention, integers equal to twice the actual spin value. For information on the 3-j coefficients see Abramowitz & Stegun, Section 27.9. The functions described in this section are declared in the header file `gsl_sf_coupling.h'.
| 7.8.1 3-j Symbols | ||
| 7.8.2 6-j Symbols | ||
| 7.8.3 9-j Symbols |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Wigner 3-j coefficient,
where the arguments are given in half-integer units, ja = two_ja/2, ma = two_ma/2, etc.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Wigner 6-j coefficient,
where the arguments are given in half-integer units, ja = two_ja/2, ma = two_ma/2, etc.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Wigner 9-j coefficient,
where the arguments are given in half-integer units, ja = two_ja/2, ma = two_ma/2, etc.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Dawson integral is defined by \exp(-x^2) \int_0^x dt \exp(t^2). A table of Dawson's integral can be found in Abramowitz & Stegun, Table 7.5. The Dawson functions are declared in the header file `gsl_sf_dawson.h'.
These routines compute the value of Dawson's integral for x.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Debye functions are defined by the integral D_n(x) = n/x^n \int_0^x dt (t^n/(e^t - 1)). For further information see Abramowitz & Stegun, Section 27.1. The Debye functions are declared in the header file `gsl_sf_debye.h'.
These routines compute the first-order Debye function D_1(x) = (1/x) \int_0^x dt (t/(e^t - 1)).
These routines compute the second-order Debye function D_2(x) = (2/x^2) \int_0^x dt (t^2/(e^t - 1)).
These routines compute the third-order Debye function D_3(x) = (3/x^3) \int_0^x dt (t^3/(e^t - 1)).
These routines compute the fourth-order Debye function D_4(x) = (4/x^4) \int_0^x dt (t^4/(e^t - 1)).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The functions described in this section are declared in the header file `gsl_sf_dilog.h'.
| 7.11.1 Real Argument | ||
| 7.11.2 Complex Argument |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the dilogarithm for a real argument. In Lewin's notation this is Li_2(x), the real part of the dilogarithm of a real x. It is defined by the integral representation Li_2(x) = - \Re \int_0^x ds \log(1-s) / s. Note that \Im(Li_2(x)) = 0 for x <= 1, and -\pi\log(x) for x > 1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
This function computes the full complex-valued dilogarithm for the complex argument z = r \exp(i \theta). The real and imaginary parts of the result are returned in result_re, result_im.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The following functions allow for the propagation of errors when combining quantities by multiplication. The functions are declared in the header file `gsl_sf_elementary.h'.
This function multiplies x and y storing the product and its associated error in result.
This function multiplies x and y with associated absolute errors dx and dy. The product xy +/- xy \sqrt((dx/x)^2 +(dy/y)^2) is stored in result.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The functions described in this section are declared in the header file `gsl_sf_ellint.h'.
| 7.13.1 Definition of Legendre Forms | ||
| 7.13.2 Definition of Carlson Forms | ||
| 7.13.3 Legendre Form of Complete Elliptic Integrals | ||
| 7.13.4 Legendre Form of Incomplete Elliptic Integrals | ||
| 7.13.5 Carlson Forms |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Legendre forms of elliptic integrals F(\phi,k), E(\phi,k) and P(\phi,k,n) are defined by,
The complete Legendre forms are denoted by K(k) = F(\pi/2, k) and E(k) = E(\pi/2, k). Further information on the Legendre forms of elliptic integrals can be found in Abramowitz & Stegun, Chapter 17. The notation used here is based on Carlson, Numerische Mathematik 33 (1979) 1 and differs slightly from that used by Abramowitz & Stegun.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Carlson symmetric forms of elliptical integrals RC(x,y), RD(x,y,z), RF(x,y,z) and RJ(x,y,z,p) are defined by,
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the complete elliptic integral K(k) to the accuracy specified by the mode variable mode.
These routines compute the complete elliptic integral E(k) to the accuracy specified by the mode variable mode.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the incomplete elliptic integral F(\phi,k) to the accuracy specified by the mode variable mode.
These routines compute the incomplete elliptic integral E(\phi,k) to the accuracy specified by the mode variable mode.
These routines compute the incomplete elliptic integral P(\phi,k,n) to the accuracy specified by the mode variable mode.
These functions compute the incomplete elliptic integral D(\phi,k,n) which is defined through the Carlson form RD(x,y,z) by the following relation,
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the incomplete elliptic integral RC(x,y) to the accuracy specified by the mode variable mode.
These routines compute the incomplete elliptic integral RD(x,y,z) to the accuracy specified by the mode variable mode.
These routines compute the incomplete elliptic integral RF(x,y,z) to the accuracy specified by the mode variable mode.
These routines compute the incomplete elliptic integral RJ(x,y,z,p) to the accuracy specified by the mode variable mode.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Jacobian Elliptic functions are defined in Abramowitz & Stegun, Chapter 16. The functions are declared in the header file `gsl_sf_elljac.h'.
This function computes the Jacobian elliptic functions sn(u|m), cn(u|m), dn(u|m) by descending Landen transformations.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The error function is described in Abramowitz & Stegun, Chapter 7. The functions in this section are declared in the header file `gsl_sf_erf.h'.
| 7.15.1 Error Function | ||
| 7.15.2 Complementary Error Function | ||
| 7.15.3 Log Complementary Error Function | ||
| 7.15.4 Probability functions |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the error function erf(x) = (2/\sqrt(\pi)) \int_0^x dt \exp(-t^2).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the complementary error function erfc(x) = 1 - erf(x) = (2/\sqrt(\pi)) \int_x^\infty \exp(-t^2).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the logarithm of the complementary error function \log(\erfc(x)).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The probability functions for the Normal or Gaussian distribution are described in Abramowitz & Stegun, Section 26.2.
These routines compute the Gaussian probability density function Z(x) = (1/\sqrt{2\pi}) \exp(-x^2/2).
These routines compute the upper tail of the Gaussian probability function Q(x) = (1/\sqrt{2\pi}) \int_x^\infty dt \exp(-t^2/2).
The hazard function for the normal distribution, also known as the inverse Mill's ratio, is defined as h(x) = Z(x)/Q(x) = \sqrt{2/\pi \exp(-x^2 / 2) / \erfc(x/\sqrt 2). It decreases rapidly as x approaches -\infty and asymptotes to h(x) \sim x as x approaches +\infty.
These routines compute the hazard function for the normal distribution.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The functions described in this section are declared in the header file `gsl_sf_exp.h'.
| 7.16.1 Exponential Function | ||
| 7.16.2 Relative Exponential Functions | ||
| 7.16.3 Exponentiation With Error Estimate |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines provide an exponential function \exp(x) using GSL semantics and error checking.
This function computes the exponential \exp(x) using the
gsl_sf_result_e10 type to return a result with extended range.
This function may be useful if the value of \exp(x) would
overflow the numeric range of double.
These routines exponentiate x and multiply by the factor y to return the product y \exp(x).
This function computes the product y \exp(x) using the
gsl_sf_result_e10 type to return a result with extended numeric
range.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the quantity \exp(x)-1 using an algorithm that is accurate for small x.
These routines compute the quantity (\exp(x)-1)/x using an algorithm that is accurate for small x. For small x the algorithm is based on the expansion (\exp(x)-1)/x = 1 + x/2 + x^2/(2*3) + x^3/(2*3*4) + \dots.
These routines compute the quantity 2(\exp(x)-1-x)/x^2 using an algorithm that is accurate for small x. For small x the algorithm is based on the expansion 2(\exp(x)-1-x)/x^2 = 1 + x/3 + x^2/(3*4) + x^3/(3*4*5) + \dots.
These routines compute the N-relative exponential, which is the
n-th generalization of the functions gsl_sf_exprel and
gsl_sf_exprel2. The N-relative exponential is given by,
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
This function exponentiates x with an associated absolute error dx.
This function exponentiates a quantity x with an associated absolute
error dx using the gsl_sf_result_e10 type to return a result with
extended range.
This routine computes the product y \exp(x) for the quantities x, y with associated absolute errors dx, dy.
This routine computes the product y \exp(x) for the quantities
x, y with associated absolute errors dx, dy using the
gsl_sf_result_e10 type to return a result with extended range.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
Information on the exponential integrals can be found in Abramowitz & Stegun, Chapter 5. These functions are declared in the header file `gsl_sf_expint.h'.
| 7.17.1 Exponential Integral | ||
| 7.17.2 Ei(x) | ||
| 7.17.3 Hyperbolic Integrals | ||
| 7.17.4 Ei_3(x) | ||
| 7.17.5 Trigonometric Integrals | ||
| 7.17.6 Arctangent Integral |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the exponential integral E_1(x),
These routines compute the second-order exponential integral E_2(x),
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the exponential integral Ei(x), where PV denotes the principal value of the integral.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the integral Shi(x) = \int_0^x dt \sinh(t)/t.
These routines compute the integral Chi(x) := Re[ \gamma_E + \log(x) + \int_0^x dt (\cosh[t]-1)/t] , where \gamma_E is the Euler constant (available as the macro M_EULER).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the exponential integral Ei_3(x) = \int_0^x dt \exp(-t^3) for x >= 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Sine integral Si(x) = \int_0^x dt \sin(t)/t.
These routines compute the Cosine integral Ci(x) = -\int_x^\infty dt \cos(t)/t for x > 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Arctangent integral AtanInt(x) = \int_0^x dt \arctan(t)/t.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The functions described in this section are declared in the header file `gsl_sf_fermi_dirac.h'.
| 7.18.1 Complete Fermi-Dirac Integrals | ||
| 7.18.2 Incomplete Fermi-Dirac Integrals |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The complete Fermi-Dirac integral F_j(x) is given by,
These routines compute the complete Fermi-Dirac integral with an index of -1. This integral is given by F_{-1}(x) = e^x / (1 + e^x).
These routines compute the complete Fermi-Dirac integral with an index of 0. This integral is given by F_0(x) = \ln(1 + e^x).
These routines compute the complete Fermi-Dirac integral with an index of 1, F_1(x) = \int_0^\infty dt (t /(\exp(t-x)+1)).
These routines compute the complete Fermi-Dirac integral with an index of 2, F_2(x) = (1/2) \int_0^\infty dt (t^2 /(\exp(t-x)+1)).
These routines compute the complete Fermi-Dirac integral with an integer index of j, F_j(x) = (1/\Gamma(j+1)) \int_0^\infty dt (t^j /(\exp(t-x)+1)).
These routines compute the complete Fermi-Dirac integral F_{-1/2}(x).
These routines compute the complete Fermi-Dirac integral F_{1/2}(x).
These routines compute the complete Fermi-Dirac integral F_{3/2}(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The incomplete Fermi-Dirac integral F_j(x,b) is given by,
These routines compute the incomplete Fermi-Dirac integral with an index of zero, F_0(x,b) = \ln(1 + e^{b-x}) - (b-x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Gamma function is defined by the following integral,
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'.
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.
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.
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).
These routines compute the regulated Gamma Function \Gamma^*(x) for x > 0. The regulated gamma function is given by,
and is a useful suggestion of Temme.
These routines compute the reciprocal of the gamma function, 1/\Gamma(x) using the real Lanczos method.
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.
These routines compute the Taylor coefficient x^n / n! for x >= 0, n >= 0.
These routines compute the factorial n!. The factorial is related to the Gamma function by n! = \Gamma(n+1).
These routines compute the double factorial n!! = n(n-2)(n-4) \dots.
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.
These routines compute the logarithm of the double factorial of n, \log(n!!).
These routines compute the combinatorial factor n choose m
= n!/(m!(n-m)!)
These routines compute the logarithm of n choose m. This is
equivalent to the sum \log(n!) - \log(m!) - \log((n-m)!).
These routines compute the Pochhammer symbol (a)_x := \Gamma(a + x)/\Gamma(a), subject to a and a+x not being negative integers. The Pochhammer symbol is also known as the Apell symbol.
These routines compute the logarithm of the Pochhammer symbol, \log((a)_x) = \log(\Gamma(a + x)/\Gamma(a)) for a > 0, a+x > 0.
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.
These routines compute the relative Pochhammer symbol ((a,x) - 1)/x where (a,x) = (a)_x := \Gamma(a + x)/\Gamma(a).
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.
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).
These functions compute the incomplete Gamma Function the normalization factor included in the previously defined functions: \Gamma(a,x) = \int_x^\infty dt t^{a-1} \exp(-t) for a real and x >= 0.
These routines compute the Beta Function, B(a,b) = \Gamma(a)\Gamma(b)/\Gamma(a+b) for a > 0, b > 0.
These routines compute the logarithm of the Beta Function, \log(B(a,b)) for a > 0, b > 0.
These routines compute the normalize incomplete Beta function B_x(a,b)/B(a,b) for a > 0, b > 0, and 0 <= x <= 1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Gegenbauer polynomials are defined in Abramowitz & Stegun, Chapter 22, where they are known as Ultraspherical polynomials. The functions described in this section are declared in the header file `gsl_sf_gegenbauer.h'.
These functions evaluate the Gegenbauer polynomials C^{(\lambda)}_n(x) using explicit representations for n =1, 2, 3.
These functions evaluate the Gegenbauer polynomial C^{(\lambda)}_n(x) for a specific value of n, lambda, x subject to \lambda > -1/2, n >= 0.
This function computes an array of Gegenbauer polynomials C^{(\lambda)}_n(x) for n = 0, 1, 2, \dots, nmax, subject to \lambda > -1/2, nmax >= 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
Hypergeometric functions are described in Abramowitz & Stegun, Chapters 13 and 15. These functions are declared in the header file `gsl_sf_hyperg.h'.
These routines compute the hypergeometric function 0F1(c,x).
These routines compute the confluent hypergeometric function 1F1(m,n,x) = M(m,n,x) for integer parameters m, n.
These routines compute the confluent hypergeometric function 1F1(a,b,x) = M(a,b,x) for general parameters a, b.
These routines compute the confluent hypergeometric function U(m,n,x) for integer parameters m, n.
This routine computes the confluent hypergeometric function
U(m,n,x) for integer parameters m, n using the
gsl_sf_result_e10 type to return a result with extended range.
These routines compute the confluent hypergeometric function U(a,b,x).
This routine computes the confluent hypergeometric function
U(a,b,x) using the gsl_sf_result_e10 type to return a
result with extended range.
These routines compute the Gauss hypergeometric function 2F1(a,b,c,x) for |x| < 1.
If the arguments (a,b,c,x) are too close to a singularity then
the function can return the error code GSL_EMAXITER when the
series approximation converges too slowly. This occurs in the region of
x=1, c - a - b = m for integer m.
These routines compute the Gauss hypergeometric function 2F1(a_R + i a_I, a_R - i a_I, c, x) with complex parameters for |x| < 1. exceptions:
These routines compute the renormalized Gauss hypergeometric function 2F1(a,b,c,x) / \Gamma(c) for |x| < 1.
These routines compute the renormalized Gauss hypergeometric function 2F1(a_R + i a_I, a_R - i a_I, c, x) / \Gamma(c) for |x| < 1.
These routines compute the hypergeometric function 2F0(a,b,x). The series representation is a divergent hypergeometric series. However, for x < 0 we have 2F0(a,b,x) = (-1/x)^a U(a,1+a-b,-1/x)
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Laguerre polynomials are defined in terms of confluent hypergeometric functions as L^a_n(x) = ((a+1)_n / n!) 1F1(-n,a+1,x). These functions are declared in the header file `gsl_sf_laguerre.h'.
These routines evaluate the generalized Laguerre polynomials L^a_1(x), L^a_2(x), L^a_3(x) using explicit representations.
These routines evaluate the generalized Laguerre polynomials L^a_n(x) for a > -1, n >= 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
Lambert's W functions, W(x), are defined to be solutions of the equation W(x) \exp(W(x)) = x. This function has multiple branches for x < 0; however, it has only two real-valued branches. We define W_0(x) to be the principal branch, where W > -1 for x < 0, and W_{-1}(x) to be the other real branch, where W < -1 for x < 0. The Lambert functions are declared in the header file `gsl_sf_lambert.h'.
These compute the principal branch of the Lambert W function, W_0(x).
These compute the secondary real-valued branch of the Lambert W function, W_{-1}(x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Legendre Functions and Legendre Polynomials are described in Abramowitz & Stegun, Chapter 8. These functions are declared in the header file `gsl_sf_legendre.h'.
| 7.24.1 Legendre Polynomials | ||
| 7.24.2 Associated Legendre Polynomials and Spherical Harmonics | ||
| 7.24.3 Conical Functions | ||
| 7.24.4 Radial Functions for Hyperbolic Space |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These functions evaluate the Legendre polynomials P_l(x) using explicit representations for l=1, 2, 3.
These functions evaluate the Legendre polynomial P_l(x) for a specific value of l, x subject to l >= 0, |x| <= 1
This function computes an array of Legendre polynomials P_l(x) for l = 0, \dots, lmax, |x| <= 1
These routines compute the Legendre function Q_0(x) for x > -1, x != 1.
These routines compute the Legendre function Q_1(x) for x > -1, x != 1.
These routines compute the Legendre function Q_l(x) for x > -1, x != 1 and l >= 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The following functions compute the associated Legendre Polynomials
P_l^m(x). Note that this function grows combinatorially with
l and can overflow for l larger than about 150. There is
no trouble for small m, but overflow occurs when m and
l are both large. Rather than allow overflows, these functions
refuse to calculate P_l^m(x) and return GSL_EOVRFLW when
they can sense that l and m are too big.
If you want to calculate a spherical harmonic, then do not use
these functions. Instead use gsl_sf_legendre_sphPlm() below,
which uses a similar recursion, but with the normalized functions.
These routines compute the associated Legendre polynomial P_l^m(x) for m >= 0, l >= m, |x| <= 1.
This function computes an array of Legendre polynomials P_l^m(x) for m >= 0, l = |m|, ..., lmax, |x| <= 1.
These routines compute the normalized associated Legendre polynomial $\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x)$ suitable for use in spherical harmonics. The parameters must satisfy m >= 0, l >= m, |x| <= 1. Theses routines avoid the overflows that occur for the standard normalization of P_l^m(x).
This function computes an array of normalized associated Legendre functions $\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x)$ for m >= 0, l = |m|, ..., lmax, |x| <= 1.0
This function returns the size of result_array[] needed for the array versions of P_l^m(x), lmax - m + 1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Conical Functions P^\mu_{-(1/2)+i\lambda}(x), Q^\mu_{-(1/2)+i\lambda} are described in Abramowitz & Stegun, Section 8.12.
These routines compute the irregular Spherical Conical Function P^{1/2}_{-1/2 + i \lambda}(x) for x > -1.
These routines compute the regular Spherical Conical Function P^{-1/2}_{-1/2 + i \lambda}(x) for x > -1.
These routines compute the conical function P^0_{-1/2 + i \lambda}(x) for x > -1.
These routines compute the conical function P^1_{-1/2 + i \lambda}(x) for x > -1.
These routines compute the Regular Spherical Conical Function P^{-1/2-l}_{-1/2 + i \lambda}(x) for x > -1, l >= -1.
These routines compute the Regular Cylindrical Conical Function P^{-m}_{-1/2 + i \lambda}(x) for x > -1, m >= -1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The following spherical functions are specializations of Legendre functions which give the regular eigenfunctions of the Laplacian on a 3-dimensional hyperbolic space H3d. Of particular interest is the flat limit, \lambda \to \infty, \eta \to 0, \lambda\eta fixed.
These routines compute the zeroth radial eigenfunction of the Laplacian on the 3-dimensional hyperbolic space, L^{H3d}_0(\lambda,\eta) := \sin(\lambda\eta)/(\lambda\sinh(\eta)) for \eta >= 0. In the flat limit this takes the form L^{H3d}_0(\lambda,\eta) = j_0(\lambda\eta).
These routines compute the first radial eigenfunction of the Laplacian on the 3-dimensional hyperbolic space, L^{H3d}_1(\lambda,\eta) := 1/\sqrt{\lambda^2 + 1} \sin(\lambda \eta)/(\lambda \sinh(\eta)) (\coth(\eta) - \lambda \cot(\lambda\eta)) for \eta >= 0. In the flat limit this takes the form L^{H3d}_1(\lambda,\eta) = j_1(\lambda\eta).
These routines compute the l'th radial eigenfunction of the Laplacian on the 3-dimensional hyperbolic space \eta >= 0, l >= 0. In the flat limit this takes the form L^{H3d}_l(\lambda,\eta) = j_l(\lambda\eta).
This function computes an array of radial eigenfunctions L^{H3d}_l(\lambda, \eta) for 0 <= l <= lmax.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
Information on the properties of the Logarithm function can be found in Abramowitz & Stegun, Chapter 4. The functions described in this section are declared in the header file `gsl_sf_log.h'.
These routines compute the logarithm of x, \log(x), for x > 0.
These routines compute the logarithm of the magnitude of x, \log(|x|), for x \ne 0.
This routine computes the complex logarithm of z = z_r + i z_i. The results are returned as lnr, theta such that \exp(lnr + i \theta) = z_r + i z_i, where \theta lies in the range [-\pi,\pi].
These routines compute \log(1 + x) for x > -1 using an algorithm that is accurate for small x.
These routines compute \log(1 + x) - x for x > -1 using an algorithm that is accurate for small x.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The following functions are equivalent to the function gsl_pow_int
(see section Small integer powers) with an error estimate. These functions are
declared in the header file `gsl_sf_pow_int.h'.
These routines compute the power x^n for integer n. The power is computed using the minimum number of multiplications. For example, x^8 is computed as ((x^2)^2)^2, requiring only 3 multiplications. For reasons of efficiency, these functions do not check for overflow or underflow conditions.
#include <gsl/gsl_sf_pow_int.h> /* compute 3.0**12 */ double y = gsl_sf_pow_int(3.0, 12); |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The polygamma functions of order m defined by \psi^{(m)}(x) = (d/dx)^m \psi(x) = (d/dx)^{m+1} \log(\Gamma(x)), where \psi(x) = \Gamma'(x)/\Gamma(x) is known as the digamma function. These functions are declared in the header file `gsl_sf_psi.h'.
| 7.27.1 Digamma Function | ||
| 7.27.2 Trigamma Function | ||
| 7.27.3 Polygamma Function |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the digamma function \psi(n) for positive integer n. The digamma function is also called the Psi function.
These routines compute the digamma function \psi(x) for general x, x \ne 0.
These routines compute the real part of the digamma function on the line 1+i y, Re[\psi(1 + i y)].
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the Trigamma function \psi'(n) for positive integer n.
These routines compute the Trigamma function \psi'(x) for general x.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the polygamma function \psi^{(m)}(x) for m >= 0, x > 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The functions described in this section are declared in the header file `gsl_sf_synchrotron.h'.
These routines compute the first synchrotron function x \int_x^\infty dt K_{5/3}(t) for x >= 0.
These routines compute the second synchrotron function x K_{2/3}(x) for x >= 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The transport functions J(n,x) are defined by the integral representations J(n,x) := \int_0^x dt t^n e^t /(e^t - 1)^2. They are declared in the header file `gsl_sf_transport.h'.
These routines compute the transport function J(2,x).
These routines compute the transport function J(3,x).
These routines compute the transport function J(4,x).
These routines compute the transport function J(5,x).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The library includes its own trigonometric functions in order to provide consistency across platforms and reliable error estimates. These functions are declared in the header file `gsl_sf_trig.h'.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the sine function \sin(x).
These routines compute the cosine function \cos(x).
These routines compute the hypotenuse function \sqrt{x^2 + y^2} avoiding overflow and underflow.
These routines compute \sinc(x) = \sin(\pi x) / (\pi x) for any value of x.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
This function computes the complex sine, \sin(z_r + i z_i) storing the real and imaginary parts in szr, szi.
This function computes the complex cosine, \cos(z_r + i z_i) storing the real and imaginary parts in szr, szi.
This function computes the logarithm of the complex sine, \log(\sin(z_r + i z_i)) storing the real and imaginary parts in szr, szi.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute \log(\sinh(x)) for x > 0.
These routines compute \log(\cosh(x)) for any x.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
This function converts the polar coordinates (r,theta) to rectilinear coordinates (x,y), x = r\cos(\theta), y = r\sin(\theta).
This function converts the rectilinear coordinates (x,y) to polar coordinates (r,theta), such that x = r\cos(\theta), y = r\sin(\theta). The argument theta lies in the range [-\pi, \pi].
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines force the angle theta to lie in the range (-\pi,\pi].
These routines force the angle theta to lie in the range [0, 2\pi).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
These routines compute the sine of an angle x with an associated absolute error dx, \sin(x \pm dx).
These routines compute the cosine of an angle x with an associated absolute error dx, \cos(x \pm dx).
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Riemann zeta function is defined in Abramowitz & Stegun, Section 23.2. The functions described in this section are declared in the header file `gsl_sf_zeta.h'.
| 7.31.1 Riemann Zeta Function | ||
| 7.31.2 Riemann Zeta Function Minus One | ||
| 7.31.3 Hurwitz Zeta Function | ||
| 7.31.4 Eta Function |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Riemann zeta function is defined by the infinite sum \zeta(s) = \sum_{k=1}^\infty k^{-s}.
These routines compute the Riemann zeta function \zeta(n) for integer n, n \ne 1.
These routines compute the Riemann zeta function \zeta(s) for arbitrary s, s \ne 1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
For large positive argument, the Riemann zeta function approaches one. In this region the fractional part is interesting, and therefore we need a function to evaluate it explicitly.
These routines compute \zeta(n) - 1 for integer n, n \ne 1.
These routines compute \zeta(s) - 1 for arbitrary s, s \ne 1.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The Hurwitz zeta function is defined by \zeta(s,q) = \sum_0^\infty (k+q)^{-s}.
These routines compute the Hurwitz zeta function \zeta(s,q) for s > 1, q > 0.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The eta function is defined by \eta(s) = (1-2^{1-s}) \zeta(s).
These routines compute the eta function \eta(n) for integer n.
These routines compute the eta function \eta(s) for arbitrary s.
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The following example demonstrates the use of the error handling form of the special functions, in this case to compute the Bessel function J_0(5.0),
|
Here are the results of running the program,
$ ./a.out
|
The next program computes the same quantity using the natural form of the function. In this case the error term result.err and return status are not accessible.
|
The results of the function are the same,
$ ./a.outJ0(5.0) = -0.177596771314338292 exact = -0.177596771314338292 |
| [ < ] | [ > ] | [ << ] | [ Up ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
The library follows the conventions of Abramowitz & Stegun where possible,
The following papers contain information on the algorithms used to compute the special functions,
| [ << ] | [ >> ] | [Top] | [Contents] | [Index] | [ ? ] |
This document was generated by Autobuild on March, 22 2005 using texi2html 1.76.