Blame doc/specfunc-legendre.rst

Packit 67cb25
.. index::
Packit 67cb25
   single: Legendre polynomials
Packit 67cb25
   single: Legendre functions
Packit 67cb25
   single: spherical harmonics
Packit 67cb25
   single: conical functions
Packit 67cb25
   single: hyperbolic space
Packit 67cb25
Packit 67cb25
The Legendre Functions and Legendre Polynomials are described in
Packit 67cb25
Abramowitz & Stegun, Chapter 8.  These functions are declared in 
Packit 67cb25
the header file :file:`gsl_sf_legendre.h`.
Packit 67cb25
Packit 67cb25
Legendre Polynomials
Packit 67cb25
--------------------
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_P1 (double x)
Packit 67cb25
              double gsl_sf_legendre_P2 (double x)
Packit 67cb25
              double gsl_sf_legendre_P3 (double x)
Packit 67cb25
              int gsl_sf_legendre_P1_e (double x, gsl_sf_result * result)
Packit 67cb25
              int gsl_sf_legendre_P2_e (double x, gsl_sf_result * result)
Packit 67cb25
              int gsl_sf_legendre_P3_e (double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These functions evaluate the Legendre polynomials
Packit 67cb25
   :math:`P_l(x)` using explicit representations for :math:`l = 1, 2, 3`.
Packit 67cb25
.. Exceptional Return Values: none
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_Pl (int l, double x)
Packit 67cb25
              int gsl_sf_legendre_Pl_e (int l, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These functions evaluate the Legendre polynomial :math:`P_l(x)`
Packit 67cb25
   for a specific value of :data:`l`, :data:`x` subject to :math:`l \ge 0` and
Packit 67cb25
   :math:`|x| \le 1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_Pl_array (int lmax, double x, double result_array[])
Packit 67cb25
              int gsl_sf_legendre_Pl_deriv_array (int lmax, double x, double result_array[], double result_deriv_array[])
Packit 67cb25
Packit 67cb25
   These functions compute arrays of Legendre polynomials
Packit 67cb25
   :math:`P_l(x)` and derivatives :math:`dP_l(x)/dx`
Packit 67cb25
   for :math:`l = 0, \dots, lmax` and :math:`|x| \le 1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_Q0 (double x)
Packit 67cb25
              int gsl_sf_legendre_Q0_e (double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the Legendre function :math:`Q_0(x)` for
Packit 67cb25
   :math:`x > -1` and :math:`x \ne 1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_Q1 (double x)
Packit 67cb25
              int gsl_sf_legendre_Q1_e (double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the Legendre function :math:`Q_1(x)` for
Packit 67cb25
   :math:`x > -1` and :math:`x \ne 1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_Ql (int l, double x)
Packit 67cb25
              int gsl_sf_legendre_Ql_e (int l, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the Legendre function :math:`Q_l(x)` for
Packit 67cb25
   :math:`x > -1`, :math:`x \ne 1` and :math:`l \ge 0`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
Associated Legendre Polynomials and Spherical Harmonics
Packit 67cb25
-------------------------------------------------------
Packit 67cb25
Packit 67cb25
The following functions compute the associated Legendre polynomials
Packit 67cb25
:math:`P_l^m(x)` which are solutions of the differential equation
Packit 67cb25
Packit 67cb25
.. only:: not texinfo
Packit 67cb25
Packit 67cb25
   .. math:: (1 - x^2) {d^2 \over dx^2} P_l^m(x) - 2x {d \over dx} P_l^m(x) +
Packit 67cb25
             \left( l(l+1) - {m^2 \over 1 - x^2} \right) P_l^m(x) = 0
Packit 67cb25
Packit 67cb25
.. only:: texinfo
Packit 67cb25
Packit 67cb25
   ::
Packit 67cb25
Packit 67cb25
      (1 - x^2) d^2 P_l^m(x) / dx^2 P_l^m(x) - 2x d/dx P_l^m(x) +
Packit 67cb25
      ( l(l+1) - m^2 / (1 - x^2) ) P_l^m(x) = 0
Packit 67cb25
Packit 67cb25
where the degree :math:`l` and order :math:`m` satisfy :math:`0 \le l` and
Packit 67cb25
:math:`0 \le m \le l`.
Packit 67cb25
The functions :math:`P_l^m(x)` grow combinatorially with
Packit 67cb25
:math:`l` and can overflow for :math:`l` larger than about 150.
Packit 67cb25
Alternatively, one may calculate normalized associated Legendre
Packit 67cb25
polynomials. There are a number of different normalization conventions,
Packit 67cb25
and these
Packit 67cb25
functions can be stably computed up to degree and order 2700. The
Packit 67cb25
following normalizations are provided:
Packit 67cb25
Packit 67cb25
* Schmidt semi-normalization
Packit 67cb25
Packit 67cb25
  Schmidt semi-normalized associated Legendre polynomials are often
Packit 67cb25
  used in the magnetics community and are defined as
Packit 67cb25
Packit 67cb25
  .. only:: not texinfo
Packit 67cb25
Packit 67cb25
     .. math::
Packit 67cb25
Packit 67cb25
        S_l^0(x) &= P_l^0(x) \\
Packit 67cb25
        S_l^m(x) &= (-1)^m \sqrt{2 {(l-m)! \over (l+m)!}} P_l^m(x), m > 0 
Packit 67cb25
Packit 67cb25
  .. only:: texinfo
Packit 67cb25
Packit 67cb25
     ::
Packit 67cb25
Packit 67cb25
        S_l^0(x) = P_l^0(x)
Packit 67cb25
        S_l^m(x) = (-1)^m \sqrt((2(l-m)! / (l+m)!)) P_l^m(x), m > 0 
Packit 67cb25
Packit 67cb25
  The factor of :math:`(-1)^m` is called the Condon-Shortley phase
Packit 67cb25
  factor and can be excluded if desired by setting the parameter
Packit 67cb25
  :code:`csphase = 1` in the functions below.
Packit 67cb25
Packit 67cb25
* Spherical Harmonic Normalization
Packit 67cb25
Packit 67cb25
  The associated Legendre polynomials suitable for calculating spherical
Packit 67cb25
  harmonics are defined as
Packit 67cb25
Packit 67cb25
  .. only:: not texinfo
Packit 67cb25
Packit 67cb25
     .. math:: Y_l^m(x) = (-1)^m \sqrt{{2l + 1 \over 4 \pi} {(l-m)! \over (l+m)!}} P_l^m(x)
Packit 67cb25
Packit 67cb25
  .. only:: texinfo
Packit 67cb25
Packit 67cb25
     ::
Packit 67cb25
Packit 67cb25
        Y_l^m(x) = (-1)^m \sqrt((2l + 1) * (l-m)! / (4 \pi) / (l+m)!) P_l^m(x)
Packit 67cb25
Packit 67cb25
  where again the phase factor :math:`(-1)^m` can be included or excluded
Packit 67cb25
  if desired.
Packit 67cb25
Packit 67cb25
* Full Normalization
Packit 67cb25
Packit 67cb25
  The fully normalized associated Legendre polynomials are defined as
Packit 67cb25
Packit 67cb25
  .. only:: not texinfo
Packit 67cb25
Packit 67cb25
     .. math:: N_l^m(x) = (-1)^m \sqrt{(l + {1 \over 2}) {(l-m)! \over (l+m)!}} P_l^m(x)
Packit 67cb25
Packit 67cb25
  .. only:: texinfo
Packit 67cb25
  
Packit 67cb25
     ::
Packit 67cb25
     
Packit 67cb25
        N_l^m(x) = (-1)^m \sqrt((l + 1/2) (l-m)! / (l+m)!) P_l^m(x)
Packit 67cb25
Packit 67cb25
  and have the property
Packit 67cb25
Packit 67cb25
  .. math:: \int_{-1}^1 N_l^m(x)^2 dx = 1
Packit 67cb25
Packit 67cb25
The normalized associated Legendre routines below use a recurrence
Packit 67cb25
relation which is stable up to a degree and order of about 2700.
Packit 67cb25
Beyond this, the computed functions could suffer from underflow
Packit 67cb25
leading to incorrect results. Routines are provided to compute
Packit 67cb25
first and second derivatives
Packit 67cb25
:math:`dP_l^m(x)/dx` and :math:`d^2 P_l^m(x)/dx^2` as well as their alternate
Packit 67cb25
versions :math:`d P_l^m(\cos{\theta})/d\theta` and
Packit 67cb25
:math:`d^2 P_l^m(\cos{\theta})/d\theta^2`. While there is a simple
Packit 67cb25
scaling relationship between the two forms, the derivatives
Packit 67cb25
involving :math:`\theta` are heavily used in spherical harmonic
Packit 67cb25
expansions and so these routines are also provided.
Packit 67cb25
Packit 67cb25
In the functions below, a parameter of type :type:`gsl_sf_legendre_t`
Packit 67cb25
specifies the type of normalization to use. The possible values are
Packit 67cb25
Packit 67cb25
.. type:: gsl_sf_legendre_t
Packit 67cb25
Packit 67cb25
   ================================== ===============================================================================
Packit 67cb25
   Value                              Description
Packit 67cb25
   ================================== ===============================================================================
Packit 67cb25
   :code:`GSL_SF_LEGENDRE_NONE`       The unnormalized associated Legendre polynomials :math:`P_l^m(x)`
Packit 67cb25
   :code:`GSL_SF_LEGENDRE_SCHMIDT`    The Schmidt semi-normalized associated Legendre polynomials :math:`S_l^m(x)`
Packit 67cb25
   :code:`GSL_SF_LEGENDRE_SPHARM`     The spherical harmonic associated Legendre polynomials :math:`Y_l^m(x)`
Packit 67cb25
   :code:`GSL_SF_LEGENDRE_FULL`       The fully normalized associated Legendre polynomials :math:`N_l^m(x)`
Packit 67cb25
   ================================== ===============================================================================
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_array (const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[])
Packit 67cb25
              int gsl_sf_legendre_array_e (const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[])
Packit 67cb25
Packit 67cb25
   These functions calculate all normalized associated Legendre
Packit 67cb25
   polynomials for :math:`0 \le l \le lmax` and
Packit 67cb25
   :math:`0 \le m \le l` for :math:`|x| \le 1`.
Packit 67cb25
   The :data:`norm` parameter specifies which normalization is used.
Packit 67cb25
   The normalized :math:`P_l^m(x)` values are stored in :data:`result_array`, whose
Packit 67cb25
   minimum size can be obtained from calling :func:`gsl_sf_legendre_array_n`.
Packit 67cb25
   The array index of :math:`P_l^m(x)` is obtained from calling
Packit 67cb25
   :code:`gsl_sf_legendre_array_index(l, m)`. To include or exclude
Packit 67cb25
   the Condon-Shortley phase factor of :math:`(-1)^m`, set the parameter
Packit 67cb25
   :data:`csphase` to either :math:`-1` or :math:`1` respectively in the
Packit 67cb25
   :code:`_e` function. This factor is excluded by default.
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_deriv_array (const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[])
Packit 67cb25
              int gsl_sf_legendre_deriv_array_e (const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[])
Packit 67cb25
Packit 67cb25
   These functions calculate all normalized associated Legendre
Packit 67cb25
   functions and their first derivatives up to degree :data:`lmax` for
Packit 67cb25
   :math:`|x| < 1`.
Packit 67cb25
   The parameter :data:`norm` specifies the normalization used. The
Packit 67cb25
   normalized :math:`P_l^m(x)` values and their derivatives
Packit 67cb25
   :math:`dP_l^m(x)/dx` are stored in :data:`result_array` and
Packit 67cb25
   :data:`result_deriv_array` respectively.
Packit 67cb25
   To include or exclude
Packit 67cb25
   the Condon-Shortley phase factor of :math:`(-1)^m`, set the parameter
Packit 67cb25
   :data:`csphase` to either :math:`-1` or :math:`1` respectively in the
Packit 67cb25
   :code:`_e` function. This factor is excluded by default.
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_deriv_alt_array (const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[])
Packit 67cb25
              int gsl_sf_legendre_deriv_alt_array_e (const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[])
Packit 67cb25
Packit 67cb25
   These functions calculate all normalized associated Legendre
Packit 67cb25
   functions and their (alternate) first derivatives up to degree :data:`lmax` for
Packit 67cb25
   :math:`|x| < 1`.
Packit 67cb25
   The normalized :math:`P_l^m(x)` values and their derivatives
Packit 67cb25
   :math:`dP_l^m(\cos{\theta})/d\theta` are stored in :data:`result_array` and
Packit 67cb25
   :data:`result_deriv_array` respectively.
Packit 67cb25
   To include or exclude
Packit 67cb25
   the Condon-Shortley phase factor of :math:`(-1)^m`, set the parameter
Packit 67cb25
   :data:`csphase` to either :math:`-1` or :math:`1` respectively in the
Packit 67cb25
   :code:`_e` function. This factor is excluded by default.
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_deriv2_array (const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[], double result_deriv2_array[])
Packit 67cb25
              int gsl_sf_legendre_deriv2_array_e (const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[], double result_deriv2_array[])
Packit 67cb25
Packit 67cb25
   These functions calculate all normalized associated Legendre
Packit 67cb25
   functions and their first and second derivatives up to degree :data:`lmax` for
Packit 67cb25
   :math:`|x| < 1`.
Packit 67cb25
   The parameter :data:`norm` specifies the normalization used. The
Packit 67cb25
   normalized :math:`P_l^m(x)`, their first derivatives
Packit 67cb25
   :math:`dP_l^m(x)/dx`, and their second derivatives
Packit 67cb25
   :math:`d^2 P_l^m(x)/dx^2` are stored in :data:`result_array`,
Packit 67cb25
   :data:`result_deriv_array`, and :data:`result_deriv2_array` respectively.
Packit 67cb25
   To include or exclude
Packit 67cb25
   the Condon-Shortley phase factor of :math:`(-1)^m`, set the parameter
Packit 67cb25
   :data:`csphase` to either :math:`-1` or :math:`1` respectively in the
Packit 67cb25
   :code:`_e` function. This factor is excluded by default.
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_deriv2_alt_array (const gsl_sf_legendre_t norm, const size_t lmax, const double x, double result_array[], double result_deriv_array[], double result_deriv2_array[])
Packit 67cb25
              int gsl_sf_legendre_deriv2_alt_array_e (const gsl_sf_legendre_t norm, const size_t lmax, const double x, const double csphase, double result_array[], double result_deriv_array[], double result_deriv2_array[])
Packit 67cb25
Packit 67cb25
   These functions calculate all normalized associated Legendre
Packit 67cb25
   functions and their (alternate) first and second derivatives up to degree
Packit 67cb25
   :data:`lmax` for
Packit 67cb25
   :math:`|x| < 1`.
Packit 67cb25
   The parameter :data:`norm` specifies the normalization used. The
Packit 67cb25
   normalized :math:`P_l^m(x)`, their first derivatives
Packit 67cb25
   :math:`dP_l^m(\cos{\theta})/d\theta`, and their second derivatives
Packit 67cb25
   :math:`d^2 P_l^m(\cos{\theta})/d\theta^2` are stored in :data:`result_array`,
Packit 67cb25
   :data:`result_deriv_array`, and :data:`result_deriv2_array` respectively.
Packit 67cb25
   To include or exclude
Packit 67cb25
   the Condon-Shortley phase factor of :math:`(-1)^m`, set the parameter
Packit 67cb25
   :data:`csphase` to either :math:`-1` or :math:`1` respectively in the
Packit 67cb25
   :code:`_e` function. This factor is excluded by default.
Packit 67cb25
Packit 67cb25
.. function:: size_t gsl_sf_legendre_array_n (const size_t lmax)
Packit 67cb25
Packit 67cb25
   This function returns the minimum array size for maximum degree :data:`lmax`
Packit 67cb25
   needed for the array versions of the associated Legendre functions.
Packit 67cb25
   Size is calculated as the total number of :math:`P_l^m(x)` functions,
Packit 67cb25
   plus extra space for precomputing multiplicative factors used in the
Packit 67cb25
   recurrence relations.
Packit 67cb25
Packit 67cb25
.. function:: size_t gsl_sf_legendre_array_index (const size_t l, const size_t m)
Packit 67cb25
Packit 67cb25
   This function returns the index into :data:`result_array`,
Packit 67cb25
   :data:`result_deriv_array`, or :data:`result_deriv2_array` corresponding
Packit 67cb25
   to :math:`P_l^m(x)`, :math:`P_l^{'m}(x)`, or :math:`P_l^{''m}(x)`. The
Packit 67cb25
   index is given by :math:`l(l+1)/2 + m`.
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_Plm (int l, int m, double x)
Packit 67cb25
              int gsl_sf_legendre_Plm_e (int l, int m, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the associated Legendre polynomial
Packit 67cb25
   :math:`P_l^m(x)` for :math:`m \ge 0`,
Packit 67cb25
   :math:`l \ge m`, and :math:`|x| \le 1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM, GSL_EOVRFLW
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_sphPlm (int l, int m, double x)
Packit 67cb25
              int gsl_sf_legendre_sphPlm_e (int l, int m, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the normalized associated Legendre polynomial
Packit 67cb25
   :math:`\sqrt{(2l+1)/(4\pi)} \sqrt{(l-m)!/(l+m)!} P_l^m(x)` suitable
Packit 67cb25
   for use in spherical harmonics.  The parameters must satisfy :math:`m \ge 0`,
Packit 67cb25
   :math:`l \ge m`, and :math:`|x| \le 1`.
Packit 67cb25
   These routines avoid the overflows
Packit 67cb25
   that occur for the standard normalization of :math:`P_l^m(x)`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_Plm_array (int lmax, int m, double x, double result_array[])
Packit 67cb25
              int gsl_sf_legendre_Plm_deriv_array (int lmax, int m, double x, double result_array[], double result_deriv_array[])
Packit 67cb25
Packit 67cb25
   These functions are now deprecated and will be removed in a future
Packit 67cb25
   release; see :func:`gsl_sf_legendre_array` and
Packit 67cb25
   :func:`gsl_sf_legendre_deriv_array`.
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_sphPlm_array (int lmax, int m, double x, double result_array[])
Packit 67cb25
              int gsl_sf_legendre_sphPlm_deriv_array (int lmax, int m, double x, double result_array[], double result_deriv_array[])
Packit 67cb25
Packit 67cb25
   These functions are now deprecated and will be removed in a future
Packit 67cb25
   release; see :func:`gsl_sf_legendre_array` and
Packit 67cb25
   :func:`gsl_sf_legendre_deriv_array`.
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_array_size (const int lmax, const int m)
Packit 67cb25
Packit 67cb25
   This function is now deprecated and will be removed in a future
Packit 67cb25
   release.
Packit 67cb25
Packit 67cb25
Conical Functions
Packit 67cb25
-----------------
Packit 67cb25
Packit 67cb25
The Conical Functions :math:`P^\mu_{-(1/2)+i\lambda}(x)`
Packit 67cb25
and :math:`Q^\mu_{-(1/2)+i\lambda}`
Packit 67cb25
are described in Abramowitz & Stegun, Section 8.12.
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_conicalP_half (double lambda, double x)
Packit 67cb25
              int gsl_sf_conicalP_half_e (double lambda, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the irregular Spherical Conical Function
Packit 67cb25
   :math:`P^{1/2}_{-1/2 + i \lambda}(x)` for :math:`x > -1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_conicalP_mhalf (double lambda, double x)
Packit 67cb25
              int gsl_sf_conicalP_mhalf_e (double lambda, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the regular Spherical Conical Function
Packit 67cb25
   :math:`P^{-1/2}_{-1/2 + i \lambda}(x)` for :math:`x > -1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_conicalP_0 (double lambda, double x)
Packit 67cb25
              int gsl_sf_conicalP_0_e (double lambda, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the conical function
Packit 67cb25
   :math:`P^0_{-1/2 + i \lambda}(x)` for :math:`x > -1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_conicalP_1 (double lambda, double x)
Packit 67cb25
              int gsl_sf_conicalP_1_e (double lambda, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the conical function 
Packit 67cb25
   :math:`P^1_{-1/2 + i \lambda}(x)` for :math:`x > -1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_conicalP_sph_reg (int l, double lambda, double x)
Packit 67cb25
              int gsl_sf_conicalP_sph_reg_e (int l, double lambda, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the Regular Spherical Conical Function
Packit 67cb25
   :math:`P^{-1/2-l}_{-1/2 + i \lambda}(x)`
Packit 67cb25
   for :math:`x > -1` and :math:`l \ge -1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_conicalP_cyl_reg (int m, double lambda, double x)
Packit 67cb25
              int gsl_sf_conicalP_cyl_reg_e (int m, double lambda, double x, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the Regular Cylindrical Conical Function
Packit 67cb25
   :math:`P^{-m}_{-1/2 + i \lambda}(x)`
Packit 67cb25
   for :math:`x > -1` and :math:`m \ge -1`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
Radial Functions for Hyperbolic Space
Packit 67cb25
-------------------------------------
Packit 67cb25
Packit 67cb25
The following spherical functions are specializations of Legendre
Packit 67cb25
functions which give the regular eigenfunctions of the Laplacian on a
Packit 67cb25
3-dimensional hyperbolic space :math:`H^3`.  Of particular interest is
Packit 67cb25
the flat limit, :math:`\lambda \to \infty`, :math:`\eta \to 0`,
Packit 67cb25
:math:`\lambda\eta` fixed.
Packit 67cb25
  
Packit 67cb25
.. function:: double gsl_sf_legendre_H3d_0 (double lambda, double eta)
Packit 67cb25
              int gsl_sf_legendre_H3d_0_e (double lambda, double eta, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the zeroth radial eigenfunction of the Laplacian on the
Packit 67cb25
   3-dimensional hyperbolic space,
Packit 67cb25
Packit 67cb25
   .. math:: L^{H3d}_0(\lambda,\eta) := {\sin(\lambda\eta) \over \lambda\sinh(\eta)}
Packit 67cb25
Packit 67cb25
   for :math:`\eta \ge 0`.
Packit 67cb25
   In the flat limit this takes the form
Packit 67cb25
   :math:`L^{H3d}_0(\lambda,\eta) = j_0(\lambda\eta)`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_H3d_1 (double lambda, double eta)
Packit 67cb25
              int gsl_sf_legendre_H3d_1_e (double lambda, double eta, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the first radial eigenfunction of the Laplacian on
Packit 67cb25
   the 3-dimensional hyperbolic space,
Packit 67cb25
Packit 67cb25
   .. math:: L^{H3d}_1(\lambda,\eta) := {1\over\sqrt{\lambda^2 + 1}} {\left(\sin(\lambda \eta)\over \lambda \sinh(\eta)\right)} \left(\coth(\eta) - \lambda \cot(\lambda\eta)\right)
Packit 67cb25
Packit 67cb25
   for :math:`\eta \ge 0`
Packit 67cb25
   In the flat limit this takes the form 
Packit 67cb25
   :math:`L^{H3d}_1(\lambda,\eta) = j_1(\lambda\eta)`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: double gsl_sf_legendre_H3d (int l, double lambda, double eta)
Packit 67cb25
              int gsl_sf_legendre_H3d_e (int l, double lambda, double eta, gsl_sf_result * result)
Packit 67cb25
Packit 67cb25
   These routines compute the :data:`l`-th radial eigenfunction of the
Packit 67cb25
   Laplacian on the 3-dimensional hyperbolic space :math:`\eta \ge 0` and
Packit 67cb25
   :math:`l \ge 0`.
Packit 67cb25
   In the flat limit this takes the form
Packit 67cb25
   :math:`L^{H3d}_l(\lambda,\eta) = j_l(\lambda\eta)`.
Packit 67cb25
.. Exceptional Return Values: GSL_EDOM
Packit 67cb25
Packit 67cb25
.. function:: int gsl_sf_legendre_H3d_array (int lmax, double lambda, double eta, double result_array[])
Packit 67cb25
Packit 67cb25
   This function computes an array of radial eigenfunctions
Packit 67cb25
   :math:`L^{H3d}_l( \lambda, \eta)`
Packit 67cb25
   for :math:`0 \le l \le lmax`.
Packit 67cb25
.. Exceptional Return Values: