Class providing methods for evaluation of Spherical Harmonics. More...
#include <spherical_harmonics.hpp>
Public Types | |
| using | real_t = double |
| using | uint_t = unsigned int |
Public Member Functions | |
| SphericalHarmonicsTool (uint_t lmax) | |
| ~SphericalHarmonicsTool () | |
| real_t | shconvert_eval (real_t *****sphdata, real_t x, real_t y, real_t z, uint_t ir, uint_t lmax, uint_t kk) |
| real_t | shconvert_eval (uint_t deg, int ord, real_t x, real_t y, real_t z) |
| Evaluate a single spherical harmonics function of given degree and order. | |
| real_t | shconvert_eval (real_t *****sphdata, real_t x, real_t y, real_t z, real_t rmin, real_t rmax, uint_t nr, uint_t lmax, uint_t kk) |
| void | plmbar (real_t *p, uint_t lmax, real_t z) |
| void | dplmbar (real_t *dplm, uint_t lmax, real_t z) |
| Compute derivatives of associated Legendre functions. | |
| real_t | evalVSH (uint_t deg, int ord, real_t x, real_t y, real_t z, uint_t ind, uint_t comp) |
| Method for evaluating vector spherical harmonics. | |
Class providing methods for evaluation of Spherical Harmonics.
This class provides routines for the evaluation and calculation of spherical harmonics. It is a C++ adaptation of old TERRA routines.
| using terra::shell::SphericalHarmonicsTool::real_t = double |
| using terra::shell::SphericalHarmonicsTool::uint_t = unsigned int |
| terra::shell::SphericalHarmonicsTool::SphericalHarmonicsTool | ( | uint_t | lmax | ) |
| terra::shell::SphericalHarmonicsTool::~SphericalHarmonicsTool | ( | ) |
Compute derivatives of associated Legendre functions.
This public method first evaluates the associated Legendre functions up to the specified degree for the given co-latitude and then their derivatives with respect to co-latitude. The actual computations are performed by the private internal methods
| dplm | array to store derivative values |
| lmax | maximal degree up to which to compute |
| z | \(\cos(\theta)\), i.e. cosine of co-latitude |
| SphericalHarmonicsTool::real_t terra::shell::SphericalHarmonicsTool::evalVSH | ( | uint_t | deg, |
| int | ord, | ||
| real_t | x, | ||
| real_t | y, | ||
| real_t | z, | ||
| uint_t | ind, | ||
| uint_t | comp | ||
| ) |
Method for evaluating vector spherical harmonics.
We use vector spherical harmonics given by
\[ \begin{split} \mathcal{Y}_{lm}^0(\theta,\phi) &= \nu_0 \vec{e}_r \\ \mathcal{Y}_{lm}^1(\theta,\phi) &= \nu_1 \vec{e}_\theta + \nu_2 \vec{e}_\phi \\ \mathcal{Y}_{lm}^2(\theta,\phi) &= - \nu_2 \vec{e}_\theta + \nu_1 \vec{e}_\phi \end{split} \]
where \(\{\vec{e}_r,\vec{e}_\theta,\vec{e}_\phi\}\) is the local triad at the point \((\theta,\phi)\) on the unit sphere and the coefficient functions are given by
\[ \begin{split} \nu_0(\theta,\phi) &= Y_{lm}(\theta,\phi) \\ \nu_1(\theta,\phi) &= \frac{1}{\sqrt{l(l+1)}}\cdot \frac{\partial}{\partial \theta} Y_{lm}(\theta,\phi) \\ \nu_2(\theta,\phi) &= \frac{1}{\sqrt{l(l+1)}}\cdot \frac{1}{\sin\theta}\frac{\partial}{\partial \phi} Y_{lm}(\theta,\phi) \end{split} \]
Here \(Y_{lm}\) is the scalar spherical harmonics function.
| deg | degree l |
| ord | order m |
| x | x-coordinate of evaluation node |
| y | y-coordinate of evaluation node |
| z | z-coordinate of evaluation node |
| ind | index for choosing vector from \(\{\mathcal{Y}_{lm}^0, \mathcal{Y}_{lm}^1,\mathcal{Y}_{lm}^2\}\) |
| comp | choosing 0th, 1st or 2nd vector component |
\[ \vec{e}_r = \begin{pmatrix} x \\ y \\ z \end{pmatrix} \enspace,\quad \vec{e}_\theta = \begin{pmatrix} zx/\varrho \\ zy/\varrho \\ -\varrho \end{pmatrix} \enspace,\quad \vec{e}_\phi = \begin{pmatrix} -y/\varrho \\ x/\varrho \\ 0 \end{pmatrix} \enspace,\quad \varrho = \sqrt{x^2 + y^2} \]
Evaluate associated Legendre functions
The subroutine evaluates all associated Legendre functions \( P_l^m \) of the 1st kind for a given argument up to a given degree l.
Normalisation of the associated Legendre functions is done in the following way:
\[ Y_l^m(\theta,\phi) = P_l^m(\cos(\theta)) \,\mbox{trig}\,(m \phi) \]
satisfy\[ \int_S Y_l^m(\theta,\phi)\cdot Y_l^m(\theta,\phi) = 4 \pi \]
Here \((\theta, \phi)\) represent colatitude and longitude. The normalisation is the one preferred in Geodesy. The normalisation factor in detail is given by\[ \sqrt{2-\delta_{0,m}}\sqrt{(2l+1)\frac{(l-m)!}{(l+m)!}} \]
The computation of the associated Legendre functions employs the so called Forward Column Method, see e.g. Holmes & Featherstone, 2002, J. of Geodesy. This is based on the three-term recurrence relation
\[ (l-m)P_l^m(z) = z(2l-1)P_{l-1}^m(z) - (l+m-1)P_{l-2}^m(z) \]
After computing \(P_m^m\) and \(P_{m+1}^m\) to initiate the recurrence we compute \(P_l^m\) by increasing l while keeping m fixed.
Routine is stable in single and real_t precision to l,m = 511 at least; timing is proportional to lmax**2 R.J.O'Connell 7 Sept. 1989; added dp(z) 10 Jan. 1990 Added precalculation and storage of square roots srl(k) 31 Dec 1992
| p | array for storing values \( P_l^m(z) \) |
| lmax | we evaluate all functions up to degree lmax |
| z | argument to \( P_l^m \), i.e. z = cos(colatitude) |
| SphericalHarmonicsTool::real_t terra::shell::SphericalHarmonicsTool::shconvert_eval | ( | real_t ***** | sphdata, |
| real_t | x, | ||
| real_t | y, | ||
| real_t | z, | ||
| real_t | rmin, | ||
| real_t | rmax, | ||
| uint_t | nr, | ||
| uint_t | lmax, | ||
| uint_t | kk | ||
| ) |
Evaluate scalar function described by spherical harmonics coefficients at given cartesian coordinates. In this variant we perform linear interpolation in radial direction between the different radial layers, for which the spherical harmonics expansion was computed.
| sphdata | Storing the data in sph format. Indices as follows: sphdata[kk][ir][ll][mm][cs]
kk: For vector data, kk=0..nk-1
ir: radial layer, counting outwards, 0..nr.
ll: Spherical harmonic coefficients, 0..lmax
mm: Spherical harmonic coefficients, 0..ll
cs: We need two coefficients, for the sine and for
the cosine, 0..1
|
| x | x-coordinate of evaluation node |
| y | y-coordinate of evaluation node |
| z | z-coordinate of evaluation node |
| rmin | radius of innermost shell of expansion |
| rmax | radius of outermost shell of expansion |
| nr | radial index in outermost shell in array sphdata |
| lmax | largest spherical harmonics degree |
| kk | index kk in array sphdata |
| SphericalHarmonicsTool::real_t terra::shell::SphericalHarmonicsTool::shconvert_eval | ( | real_t ***** | sphdata, |
| real_t | x, | ||
| real_t | y, | ||
| real_t | z, | ||
| uint_t | ir, | ||
| uint_t | lmax, | ||
| uint_t | kk | ||
| ) |
Evaluate scalar function described by spherical harmonics coefficients at given cartesian coordinates. The point has to be on the sphere.
| sphdata | Storing the data in sph format. Indices as follows: sphdata[kk][ir][ll][mm][cs]
kk: For vector data, kk=0..nk-1
ir: radial layer, counting outwards, 0..nr-1.
ll: Spherical harmonic coefficients, 0..lmax
mm: Spherical harmonic coefficients, 0..ll
cs: We need two coefficients, for the sine and for
the cosine, 0..1
|
| x | x-coordinate of evaluation node |
| y | y-coordinate of evaluation node |
| z | z-coordinate of evaluation node |
| ir | radial index in array sphdata |
| lmax | largest spherical harmonics degree |
| kk | index kk in array sphdata |
| SphericalHarmonicsTool::real_t terra::shell::SphericalHarmonicsTool::shconvert_eval | ( | uint_t | deg, |
| int | ord, | ||
| real_t | x, | ||
| real_t | y, | ||
| real_t | z | ||
| ) |
Evaluate a single spherical harmonics function of given degree and order.
Evaluate a single spherical harmonics function of specified degree and order (l,m) at given cartesian coordinates (x,y,z). For positive argument ord we use the sine variant, otherwise the cosine variant.
| deg | degree l |
| ord | order m (also decides on sin/cos) |
| x | x-coordinate of evaluation node |
| y | y-coordinate of evaluation node |
| z | z-coordinate of evaluation node |
| real_t* terra::shell::SphericalHarmonicsTool::f1 |
| real_t * terra::shell::SphericalHarmonicsTool::f2 |
| real_t * terra::shell::SphericalHarmonicsTool::fac1 |
| real_t * terra::shell::SphericalHarmonicsTool::fac2 |
| real_t * terra::shell::SphericalHarmonicsTool::srt |