spharpy.special#

Module implementing or wrapping special functions required for higher level classes and functions of the spharpy package.

Functions:

chebyshev_coefficients(order)

Calculate the coefficients of a Chebyshev polynomial of the specified order.

legendre_coefficients(order)

Calculate the coefficients of a Legendre polynomial of the specified order.

legendre_function(n, m, z[, cs_phase])

Legendre function of order n and degree m with argument z.

spherical_bessel(n, z[, derivative])

Spherical Bessel function of the first kind and of order n evaluated at z.

spherical_bessel_zeros(n_max, n_zeros)

Compute the zeros of the spherical Bessel function.

spherical_hankel(n, z[, kind, derivative])

Spherical Hankel function of the first kind and order n evaluated at z.

spherical_harmonic(n, m, theta, phi)

The spherical harmonics of order n and degree m.

spherical_harmonic_derivative_phi(n, m, ...)

Calculate the derivative of the spherical harmonics with respect to the azimuth angle phi.

spherical_harmonic_derivative_phi_real(n, m, ...)

The derivative of the real valued spherical harmonics with respect to the azimuth angle \(\phi\).

spherical_harmonic_derivative_theta(n, m, ...)

Calculate the derivative of the spherical harmonics with respect to the elevation angle theta.

spherical_harmonic_derivative_theta_real(n, ...)

The derivative of the real valued spherical harmonics with respect to the elevation angle \(\theta\).

spherical_harmonic_gradient_phi(n, m, theta, phi)

Calculate the derivative of the spherical harmonics with respect to the azimuth angle phi divided by sin(theta).

spherical_harmonic_gradient_phi_real(n, m, ...)

The gradient of the real valued spherical harmonics with respect to the azimuth angle \(\phi\).

spherical_harmonic_normalization(n, m[, norm])

The normalization factor for real valued spherical harmonics.

spherical_harmonic_real(n, m, theta, phi)

Real valued spherical harmonic function of order n and degree m evaluated at the angles theta and phi.

spharpy.special.chebyshev_coefficients(order)[source]#

Calculate the coefficients of a Chebyshev polynomial of the specified order.

Parameters:

order (int) – The order of the polynomial

Returns:

coefficients – The coefficients of the polynomial

Return type:

ndarray, double

Note

This uses numpy.polynomial.chebyshev to compute the coefficients.

spharpy.special.legendre_coefficients(order)[source]#

Calculate the coefficients of a Legendre polynomial of the specified order.

Parameters:

order (int) – The order of the polynomial

Returns:

coefficients – The coefficients of the polynomial

Return type:

ndarray, double

Note

This uses numpy.polynomial.legendre to compute the coefficients.

spharpy.special.legendre_function(n, m, z, cs_phase=True)[source]#

Legendre function of order n and degree m with argument z.

The Legendre function is defined as ([1] Eq. 1.30, 1.33, and 1.34)

\[P_n^m(z) = (-1)^m(1-z^2)^{m/2}\frac{d^m}{dz^m}P_n(z),\]

where the Condon-Shortley phase term \((-1)^m\) is dropped when cs_phase = False and the Legendre Polynomial \(P_n(z)\).

Parameters:
  • n (int) – The order

  • m (int) – The degree

  • z (ndarray, double) – The argument as an array

  • cs_phase (bool, optional) – Whether to include the Condon-Shortley phase term \((-1)^m\) or not. Default is True.

Returns:

legendre – The Legendre function. This will return zeros if \(|m| > n\).

Return type:

ndarray, double

References

Note

This is a wrapper for the Legendre function implementation from scipy. The scipy implementation uses the Condon-Shortley phase. Therefore, the sign needs to be flipped here for uneven degrees when dropping the Condon-Shortley phase.

spharpy.special.spherical_bessel(n, z, derivative=False)[source]#

Spherical Bessel function of the first kind and of order n evaluated at z.

The spherical Bessel function of the first kind is defined as ([2], Eq. 2.29)

\[j_n(z) = \sqrt{\frac{\pi}{2z}} J_{n+\frac{1}{2}} (z)\]

with the Bessel function of the first kind \(J_\alpha\).

Parameters:
  • n (int, ndarray) – Order of the spherical Bessel function

  • z (double, ndarray) – Argument of the spherical Bessel function. Has to be real valued.

  • derivative (bool) – Return the derivative of the spherical Bessel function. Default is False.

Returns:

jn – Spherical Bessel function. Array with dimensions [N x Z], where N is the number of elements in n and Z is the number of elements in z.

Return type:

double, ndarray

Note

This is a wrapper around the Scipy implementation of the spherical Bessel function.

References

spharpy.special.spherical_bessel_zeros(n_max, n_zeros)[source]#

Compute the zeros of the spherical Bessel function.

This function will always start at order zero which is equal to \(\sin(z)/z\) and iteratively compute the roots for higher orders. The roots are computed using Brents algorithm from the scipy package. See spherical_bessel for more information.

Parameters:
  • n_max (int) – The order of the spherical bessel function

  • n_zeros (int) – The number of roots to be computed

Returns:

roots – The roots of the spherical bessel function

Return type:

ndarray, double

spharpy.special.spherical_hankel(n, z, kind=2, derivative=False)[source]#

Spherical Hankel function of the first kind and order n evaluated at z.

The spherical Hankel function of the first kind is defined as ([3], Eq. 2.30)

\[h_n(z) = \sqrt{\frac{\pi}{2z}} H_{n+\frac{1}{2}} (z)\]

with the Hankel function of the first kind \(H_\alpha\).

Parameters:
  • n (int, ndarray) – Order of the spherical bessel function

  • z (double, ndarray) – Argument of the spherical bessel function. Has to be real valued.

  • kind (int, optional) – The kind of the spherical Hankel function. Must be 1 or 2. The default is 2.

  • derivative (bool, optional) – Flag to indicate if the derivative os the spherical Hankel function should be returned. The default is False.

Returns:

hn – Spherical bessel function. Array with dimensions [N x Z], where N is the number of elements in n and Z is the number of elements in z.

Return type:

double, ndarray

References

Note

This is based on the Hankel functions implemented in the scipy package.

spharpy.special.spherical_harmonic(n, m, theta, phi)[source]#

The spherical harmonics of order n and degree m. The spherical harmonic functions are fully normalized (N3D) and include the Condon-Shortley phase according to [4].

nunsigned int

The spherical harmonic order

mint

The spherical harmonic degree

thetandarray, double

The colatitude angle

phindarray, double

The azimuth angle

Returns:

spherical_harmonic – The complex valued spherical harmonic of order n and degree m

Return type:

ndarray, double

Note

This function wraps the spherical harmonic implementation from scipy.

References

spharpy.special.spherical_harmonic_derivative_phi(n, m, theta, phi)[source]#

Calculate the derivative of the spherical harmonics with respect to the azimuth angle phi.

Parameters:
  • n (int) – Spherical harmonic order

  • m (int) – Spherical harmonic degree

  • theta (double) – Elevation angle 0 < theta < pi

  • phi (double) – Azimuth angle 0 < phi < 2*pi

Returns:

sh_diff – Spherical harmonic derivative

Return type:

complex double

spharpy.special.spherical_harmonic_derivative_phi_real(n, m, theta, phi)[source]#

The derivative of the real valued spherical harmonics with respect to the azimuth angle \(\phi\).

Parameters:
  • n (int) – The spherical harmonic order.

  • m (int) – The spherical harmonic degree.

  • theta (ndarray, double) – The elevation angle

  • phi (ndarray, double) – The azimuth angle

Returns:

derivative – The derivative

Return type:

ndarray, double

Note

This implementation does not include the Condon-Shortley phase term.

spharpy.special.spherical_harmonic_derivative_theta(n, m, theta, phi)[source]#

Calculate the derivative of the spherical harmonics with respect to the elevation angle theta.

Parameters:
  • n (int) – Spherical harmonic order

  • m (int) – Spherical harmonic degree

  • theta (double) – Elevation angle 0 < theta < pi

  • phi (double) – Azimuth angle 0 < phi < 2*pi

Returns:

sh_diff – Spherical harmonic derivative

Return type:

complex double

spharpy.special.spherical_harmonic_derivative_theta_real(n, m, theta, phi)[source]#

The derivative of the real valued spherical harmonics with respect to the elevation angle \(\theta\).

Parameters:
  • n (int) – The spherical harmonic order.

  • m (int) – The spherical harmonic degree.

  • theta (ndarray, double) – The elevation angle

  • phi (ndarray, double) – The azimuth angle

Returns:

derivative – The derivative

Return type:

ndarray, double

Note

This implementation does not include the Condon-Shortley phase term.

spharpy.special.spherical_harmonic_gradient_phi(n, m, theta, phi)[source]#

Calculate the derivative of the spherical harmonics with respect to the azimuth angle phi divided by sin(theta).

Parameters:
  • n (int) – Spherical harmonic order

  • m (int) – Spherical harmonic degree

  • theta (double) – Elevation angle 0 < theta < pi

  • phi (double) – Azimuth angle 0 < phi < 2*pi

Returns:

sh_diff – Spherical harmonic derivative

Return type:

complex double

spharpy.special.spherical_harmonic_gradient_phi_real(n, m, theta, phi)[source]#

The gradient of the real valued spherical harmonics with respect to the azimuth angle \(\phi\). The singularities at the poles are avoided using the formulation proposed in [5].

Parameters:
  • n (int) – The spherical harmonic order.

  • m (int) – The spherical harmonic degree.

  • theta (ndarray, double) – The elevation angle

  • phi (ndarray, double) – The azimuth angle

Returns:

derivative – The derivative

Return type:

ndarray, double

Note

This implementation does not include the Condon-Shortley phase term.

References

spharpy.special.spherical_harmonic_normalization(n, m, norm='full')[source]#

The normalization factor for real valued spherical harmonics.

\[N_n^m = \sqrt{\frac{2n+1}{4\pi}\frac{(n-m)!}{(n+m)!}}\]
Parameters:
  • n (int) – The spherical harmonic order.

  • m (int) – The spherical harmonic degree.

  • norm ('full', 'semi', optional) – Normalization to use. Can be either fully normalized on the sphere or semi-normalized. Default is 'full'.

Returns:

norm – The normalization factor.

Return type:

double

spharpy.special.spherical_harmonic_real(n, m, theta, phi)[source]#

Real valued spherical harmonic function of order n and degree m evaluated at the angles theta and phi. The spherical harmonic functions are fully normalized (N3D) and follow the AmbiX phase convention, which does not include the Condon-Shortley phase [6].

\[Y_n^m(\theta, \phi) = \sqrt{\frac{2n+1}{4\pi} \frac{(n-|m|)!}{(n+|m|)!}} P_n^{|m|}(\cos \theta) \begin{cases} \displaystyle \cos(|m|\phi), & \text{if $m \ge 0$} \newline \displaystyle \sin(|m|\phi) , & \text{if $m < 0$} \end{cases}\]

References

Parameters:
  • n (unsigned int) – The spherical harmonic order

  • m (int) – The spherical harmonic degree

  • theta (ndarray, double) – The colatitude angle

  • phi (ndarray, double) – The azimuth angle

Returns:

spherical_harmonic – The real valued spherial harmonic of order n and degree m

Return type:

ndarray, double