|
Packit |
67cb25 |
.. index:: elliptic integrals
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
The functions described in this section are declared in the header
|
|
Packit |
67cb25 |
file :file:`gsl_sf_ellint.h`. Further information about the elliptic
|
|
Packit |
67cb25 |
integrals can be found in Abramowitz & Stegun, Chapter 17.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Definition of Legendre Forms
|
|
Packit |
67cb25 |
----------------------------
|
|
Packit |
67cb25 |
.. index:: Legendre forms of elliptic integrals
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
The Legendre forms of elliptic integrals :math:`F(\phi,k)`,
|
|
Packit |
67cb25 |
:math:`E(\phi,k)` and :math:`\Pi(\phi,k,n)` are defined by,
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. only:: not texinfo
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. math::
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
F(\phi,k) &= \int_0^\phi dt {1 \over \sqrt{(1 - k^2 \sin^2(t))}} \\
|
|
Packit |
67cb25 |
E(\phi,k) &= \int_0^\phi dt \sqrt{(1 - k^2 \sin^2(t))} \\
|
|
Packit |
67cb25 |
\Pi(\phi,k,n) &= \int_0^\phi dt {1 \over (1 + n \sin^2(t)) \sqrt{1 - k^2 \sin^2(t)}}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. only:: texinfo
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
::
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
F(\phi,k) = \int_0^\phi dt 1/\sqrt((1 - k^2 \sin^2(t)))
|
|
Packit |
67cb25 |
E(\phi,k) = \int_0^\phi dt \sqrt((1 - k^2 \sin^2(t)))
|
|
Packit |
67cb25 |
Pi(\phi,k,n) = \int_0^\phi dt 1/((1 + n \sin^2(t))\sqrt(1 - k^2 \sin^2(t)))
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
The complete Legendre forms are denoted by :math:`K(k) = F(\pi/2, k)` and
|
|
Packit |
67cb25 |
:math:`E(k) = E(\pi/2, k)`.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
The notation used here is based on Carlson, "Numerische
|
|
Packit |
67cb25 |
Mathematik" 33 (1979) 1 and differs slightly from that used by
|
|
Packit |
67cb25 |
Abramowitz & Stegun, where the functions are given in terms of the
|
|
Packit |
67cb25 |
parameter :math:`m = k^2` and :math:`n` is replaced by :math:`-n`.
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Definition of Carlson Forms
|
|
Packit |
67cb25 |
---------------------------
|
|
Packit |
67cb25 |
.. index:: Carlson forms of Elliptic integrals
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
The Carlson symmetric forms of elliptical integrals :math:`RC(x,y)`,
|
|
Packit |
67cb25 |
:math:`RD(x,y,z)`, :math:`RF(x,y,z)` and :math:`RJ(x,y,z,p)` are defined
|
|
Packit |
67cb25 |
by,
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. only:: not texinfo
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. math::
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
RC(x,y) &= 1/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1} \\
|
|
Packit |
67cb25 |
RD(x,y,z) &= 3/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1/2} (t+z)^{-3/2} \\
|
|
Packit |
67cb25 |
RF(x,y,z) &= 1/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1/2} (t+z)^{-1/2} \\
|
|
Packit |
67cb25 |
RJ(x,y,z,p) &= 3/2 \int_0^\infty dt (t+x)^{-1/2} (t+y)^{-1/2} (t+z)^{-1/2} (t+p)^{-1}
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. only:: texinfo
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
::
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
RC(x,y) = 1/2 \int_0^\infty dt (t+x)^(-1/2) (t+y)^(-1)
|
|
Packit |
67cb25 |
RD(x,y,z) = 3/2 \int_0^\infty dt (t+x)^(-1/2) (t+y)^(-1/2) (t+z)^(-3/2)
|
|
Packit |
67cb25 |
RF(x,y,z) = 1/2 \int_0^\infty dt (t+x)^(-1/2) (t+y)^(-1/2) (t+z)^(-1/2)
|
|
Packit |
67cb25 |
RJ(x,y,z,p) = 3/2 \int_0^\infty dt
|
|
Packit |
67cb25 |
(t+x)^(-1/2) (t+y)^(-1/2) (t+z)^(-1/2) (t+p)^(-1)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Legendre Form of Complete Elliptic Integrals
|
|
Packit |
67cb25 |
--------------------------------------------
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_Kcomp (double k, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_Kcomp_e (double k, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the complete elliptic integral :math:`K(k)` to
|
|
Packit |
67cb25 |
the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
Note that Abramowitz & Stegun define this function in terms of the
|
|
Packit |
67cb25 |
parameter :math:`m = k^2`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_Ecomp (double k, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_Ecomp_e (double k, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the complete elliptic integral :math:`E(k)` to the
|
|
Packit |
67cb25 |
accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
Note that Abramowitz & Stegun define this function in terms of the
|
|
Packit |
67cb25 |
parameter :math:`m = k^2`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_Pcomp (double k, double n, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_Pcomp_e (double k, double n, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the complete elliptic integral :math:`\Pi(k,n)` to the
|
|
Packit |
67cb25 |
accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
Note that Abramowitz & Stegun define this function in terms of the
|
|
Packit |
67cb25 |
parameters :math:`m = k^2` and :math:`\sin^2(\alpha) = k^2`, with the
|
|
Packit |
67cb25 |
change of sign :math:`n \to -n`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Legendre Form of Incomplete Elliptic Integrals
|
|
Packit |
67cb25 |
----------------------------------------------
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_F (double phi, double k, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_F_e (double phi, double k, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`F(\phi,k)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
Note that Abramowitz & Stegun define this function in terms of the
|
|
Packit |
67cb25 |
parameter :math:`m = k^2`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_E (double phi, double k, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_E_e (double phi, double k, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`E(\phi,k)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
Note that Abramowitz & Stegun define this function in terms of the
|
|
Packit |
67cb25 |
parameter :math:`m = k^2`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_P (double phi, double k, double n, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_P_e (double phi, double k, double n, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`\Pi(\phi,k,n)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
Note that Abramowitz & Stegun define this function in terms of the
|
|
Packit |
67cb25 |
parameters :math:`m = k^2` and :math:`\sin^2(\alpha) = k^2`, with the
|
|
Packit |
67cb25 |
change of sign :math:`n \to -n`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_D (double phi, double k, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_D_e (double phi, double k, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These functions compute the incomplete elliptic integral
|
|
Packit |
67cb25 |
:math:`D(\phi,k)` which is defined through the Carlson form :math:`RD(x,y,z)`
|
|
Packit |
67cb25 |
by the following relation,
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. only:: not texinfo
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. math:: D(\phi,k) = {1 \over 3} (\sin \phi)^3 RD (1-\sin^2(\phi), 1-k^2 \sin^2(\phi), 1)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. only:: texinfo
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
::
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
D(\phi,k) = (1/3)(\sin(\phi))^3 RD (1-\sin^2(\phi), 1-k^2 \sin^2(\phi), 1).
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
Carlson Forms
|
|
Packit |
67cb25 |
-------------
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_RC (double x, double y, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_RC_e (double x, double y, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`RC(x,y)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_RD (double x, double y, double z, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_RD_e (double x, double y, double z, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`RD(x,y,z)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_RF (double x, double y, double z, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_RF_e (double x, double y, double z, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`RF(x,y,z)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
.. function:: double gsl_sf_ellint_RJ (double x, double y, double z, double p, gsl_mode_t mode)
|
|
Packit |
67cb25 |
int gsl_sf_ellint_RJ_e (double x, double y, double z, double p, gsl_mode_t mode, gsl_sf_result * result)
|
|
Packit |
67cb25 |
|
|
Packit |
67cb25 |
These routines compute the incomplete elliptic integral :math:`RJ(x,y,z,p)`
|
|
Packit |
67cb25 |
to the accuracy specified by the mode variable :data:`mode`.
|
|
Packit |
67cb25 |
.. Exceptional Return Values: GSL_EDOM
|