blob: 50c93c3ca01e8b09f5f859ed6e425ea22acd42a0 [file] [log] [blame]
@cindex special functions
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 in this chapter are declared in individual header files,
such as @file{gsl_sf_airy.h}, @file{gsl_sf_bessel.h}, etc. The complete
set of header files can be included using the file @file{gsl_sf.h}.
@menu
* Special Function Usage::
* The gsl_sf_result struct::
* Special Function Modes::
* Airy Functions and Derivatives::
* Bessel Functions::
* Clausen Functions::
* Coulomb Functions::
* Coupling Coefficients::
* Dawson Function::
* Debye Functions::
* Dilogarithm::
* Elementary Operations::
* Elliptic Integrals::
* Elliptic Functions (Jacobi)::
* Error Functions::
* Exponential Functions::
* Exponential Integrals::
* Fermi-Dirac Function::
* Gamma and Beta Functions::
* Gegenbauer Functions::
* Hypergeometric Functions::
* Laguerre Functions::
* Lambert W Functions::
* Legendre Functions and Spherical Harmonics::
* Logarithm and Related Functions::
* Mathieu Functions::
* Power Function::
* Psi (Digamma) Function::
* Synchrotron Functions::
* Transport Functions::
* Trigonometric Functions::
* Zeta Functions::
* Special Functions Examples::
* Special Functions References and Further Reading::
@end menu
@node Special Function Usage
@section Usage
The special functions are available in two calling conventions, a
@dfn{natural form} which returns the numerical value of the function and
an @dfn{error-handling form} which returns an error code. The two types
of function provide alternative ways of accessing the same underlying
code.
The @dfn{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
@math{J_0(x)},
@example
double y = gsl_sf_bessel_J0 (x);
@end example
@noindent
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,
@example
gsl_sf_result result;
int status = gsl_sf_bessel_J0_e (x, &result);
@end example
@noindent
The error-handling functions have the suffix @code{_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 @code{GSL_SUCCESS}.
@node The gsl_sf_result struct
@section The gsl_sf_result struct
@cindex gsl_sf_result
@cindex gsl_sf_result_e10
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 @file{gsl_sf_result.h}.
The @code{gsl_sf_result} struct contains value and error fields.
@example
typedef struct
@{
double val;
double err;
@} gsl_sf_result;
@end example
@noindent
The field @var{val} contains the value and the field @var{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
@code{gsl_sf_result_e10} struct contains value and error fields as well
as an exponent field such that the actual result is obtained as
@code{result * 10^(e10)}.
@example
typedef struct
@{
double val;
double err;
int e10;
@} gsl_sf_result_e10;
@end example
@node Special Function Modes
@section Modes
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 @code{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,
@table @code
@item GSL_PREC_DOUBLE
Double-precision, a relative accuracy of approximately @c{$2 \times 10^{-16}$}
@math{2 * 10^-16}.
@item GSL_PREC_SINGLE
Single-precision, a relative accuracy of approximately @c{$1 \times 10^{-7}$}
@math{10^-7}.
@item GSL_PREC_APPROX
Approximate values, a relative accuracy of approximately @c{$5 \times 10^{-4}$}
@math{5 * 10^-4}.
@end table
@noindent
The approximate mode provides the fastest evaluation at the lowest
accuracy.
@node Airy Functions and Derivatives
@section Airy Functions and Derivatives
@include specfunc-airy.texi
@node Bessel Functions
@section Bessel Functions
@include specfunc-bessel.texi
@node Clausen Functions
@section Clausen Functions
@include specfunc-clausen.texi
@node Coulomb Functions
@section Coulomb Functions
@include specfunc-coulomb.texi
@node Coupling Coefficients
@section Coupling Coefficients
@include specfunc-coupling.texi
@node Dawson Function
@section Dawson Function
@include specfunc-dawson.texi
@node Debye Functions
@section Debye Functions
@include specfunc-debye.texi
@node Dilogarithm
@section Dilogarithm
@include specfunc-dilog.texi
@node Elementary Operations
@section Elementary Operations
@include specfunc-elementary.texi
@node Elliptic Integrals
@section Elliptic Integrals
@include specfunc-ellint.texi
@node Elliptic Functions (Jacobi)
@section Elliptic Functions (Jacobi)
@include specfunc-elljac.texi
@node Error Functions
@section Error Functions
@include specfunc-erf.texi
@node Exponential Functions
@section Exponential Functions
@include specfunc-exp.texi
@node Exponential Integrals
@section Exponential Integrals
@include specfunc-expint.texi
@node Fermi-Dirac Function
@section Fermi-Dirac Function
@include specfunc-fermi-dirac.texi
@node Gamma and Beta Functions
@section Gamma and Beta Functions
@include specfunc-gamma.texi
@node Gegenbauer Functions
@section Gegenbauer Functions
@include specfunc-gegenbauer.texi
@node Hypergeometric Functions
@section Hypergeometric Functions
@include specfunc-hyperg.texi
@node Laguerre Functions
@section Laguerre Functions
@include specfunc-laguerre.texi
@node Lambert W Functions
@section Lambert W Functions
@include specfunc-lambert.texi
@node Legendre Functions and Spherical Harmonics
@section Legendre Functions and Spherical Harmonics
@include specfunc-legendre.texi
@node Logarithm and Related Functions
@section Logarithm and Related Functions
@include specfunc-log.texi
@node Mathieu Functions
@section Mathieu Functions
@include specfunc-mathieu.texi
@node Power Function
@section Power Function
@include specfunc-pow-int.texi
@node Psi (Digamma) Function
@section Psi (Digamma) Function
@include specfunc-psi.texi
@node Synchrotron Functions
@section Synchrotron Functions
@include specfunc-synchrotron.texi
@node Transport Functions
@section Transport Functions
@include specfunc-transport.texi
@node Trigonometric Functions
@section Trigonometric Functions
@include specfunc-trig.texi
@node Zeta Functions
@section Zeta Functions
@include specfunc-zeta.texi
@node Special Functions Examples
@section Examples
The following example demonstrates the use of the error handling form of
the special functions, in this case to compute the Bessel function
@math{J_0(5.0)},
@example
@verbatiminclude examples/specfun_e.c
@end example
@noindent
Here are the results of running the program,
@example
$ ./a.out
@verbatiminclude examples/specfun_e.out
@end example
@noindent
The next program computes the same quantity using the natural form of
the function. In this case the error term @var{result.err} and return
status are not accessible.
@example
@verbatiminclude examples/specfun.c
@end example
@noindent
The results of the function are the same,
@example
$ ./a.out
@verbatiminclude examples/specfun.out
@end example
@node Special Functions References and Further Reading
@section References and Further Reading
The library follows the conventions of @cite{Abramowitz & Stegun} where
possible,
@itemize @asis
@item
Abramowitz & Stegun (eds.), @cite{Handbook of Mathematical Functions}
@end itemize
@noindent
The following papers contain information on the algorithms used
to compute the special functions,
@cindex MISCFUN
@itemize @asis
@item
MISCFUN: A software package to compute uncommon special functions.
@cite{ACM Trans.@: Math.@: Soft.}, vol.@: 22, 1996, 288--301
@item
G.N. Watson, A Treatise on the Theory of Bessel Functions,
2nd Edition (Cambridge University Press, 1944).
@item
G. Nemeth, Mathematical Approximations of Special Functions,
Nova Science Publishers, ISBN 1-56072-052-2
@item
B.C. Carlson, Special Functions of Applied Mathematics (1977)
@item
W.J. Thompson, Atlas for Computing Mathematical Functions, John Wiley & Sons,
New York (1997).
@item
Y.Y. Luke, Algorithms for the Computation of Mathematical Functions, Academic
Press, New York (1977).
@comment @item
@comment Fermi-Dirac functions of orders @math{-1/2}, @math{1/2}, @math{3/2}, and
@comment @math{5/2}. @cite{ACM Trans. Math. Soft.}, vol. 24, 1998, 1-12.
@end itemize