1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159 160 161 162 163 164 165 166 167 168 169 170 171 172 173 174 175 176 177 178 179 180 181 182 183 184 185 186 187 188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204 205 206 207 208 209 210 211 212 213 214 215 216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240 241 242 243 244 245 246 247 248 249 250 251 252 253 254 255 256 257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313 314 315 316 317 318 319 320 321 322 323 324 325 326 327 328 329 330 331 332 333 334 335 336 337 338 339 340 341 342 343 344 345 346 347 348 349 350 351 352 353 354 355 356 357 358 359 360 361 362 363 364 365
|
@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
@tindex gsl_sf_result
@tindex 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.txt
@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.txt
@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 @w{}
@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 @w{}
@item
Allan J. MacLeod, 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
N. M. Temme, Special Functions: An Introduction to the Classical
Functions of Mathematical Physics (1996), ISBN 978-0471113133.
@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).
@item
S. A. Holmes and W. E. Featherstone, A unified approach to the Clenshaw
summation and the recursive computation of very high degree and order
normalised associated Legendre functions, Journal of Geodesy, 76,
pg. 279-299, 2002.
@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
|