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.
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
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.
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].