scipy.special.sph_harm_y#
- scipy.special.sph_harm_y(n, m, theta, phi, *, diff_n=0) = <scipy.special._multiufuncs.MultiUFunc object>[source]#
Spherical harmonics. They are defined as
\[Y_n^m(\theta,\phi) = \sqrt{\frac{2 n + 1}{4 \pi} \frac{(n - m)!}{(n + m)!}} P_n^m(\cos(\theta)) e^{i m \phi}\]where \(P_n^m\) are the (unnormalized) associated Legendre polynomials.
- Parameters:
- nArrayLike[int]
Degree of the harmonic. Must have
n >= 0. This is often denoted byl(lower case L) in descriptions of spherical harmonics.- mArrayLike[int]
Order of the harmonic.
- thetaArrayLike[float]
Polar (colatitudinal) coordinate; must be in
[0, pi].- phiArrayLike[float]
Azimuthal (longitudinal) coordinate; must be in
[0, 2*pi].- diff_nOptional[int]
A non-negative integer. Compute and return all derivatives up to order
diff_n. Default is 0.
- Returns:
- yndarray[complex] or tuple[ndarray[complex]]
Spherical harmonics with
diff_nderivatives.
Notes
There are different conventions for the meanings of the input arguments
thetaandphi. In SciPythetais the polar angle andphiis the azimuthal angle. It is common to see the opposite convention, that is,thetaas the azimuthal angle andphias the polar angle.Note that SciPy’s spherical harmonics include the Condon-Shortley phase [2] because it is part of
sph_legendre_p.With SciPy’s conventions, the first several spherical harmonics are
\[\begin{split}Y_0^0(\theta, \phi) &= \frac{1}{2} \sqrt{\frac{1}{\pi}} \\ Y_1^{-1}(\theta, \phi) &= \frac{1}{2} \sqrt{\frac{3}{2\pi}} e^{-i\phi} \sin(\theta) \\ Y_1^0(\theta, \phi) &= \frac{1}{2} \sqrt{\frac{3}{\pi}} \cos(\theta) \\ Y_1^1(\theta, \phi) &= -\frac{1}{2} \sqrt{\frac{3}{2\pi}} e^{i\phi} \sin(\theta).\end{split}\]References
[1]Digital Library of Mathematical Functions, 14.30. https://dlmf.nist.gov/14.30