Public Types | Public Member Functions | Private Types | Private Attributes | List of all members
GeographicLib::SphericalHarmonic Class Reference

Spherical harmonic series. More...

#include <SphericalHarmonic.hpp>

Public Types

enum  normalization { FULL, SCHMIDT }
 

Public Member Functions

CircularEngine Circle (real p, real z, bool gradp) const
 
const SphericalEngine::coeffCoefficients () const
 
Math::real operator() (real x, real y, real z) const
 
Math::real operator() (real x, real y, real z, real &gradx, real &grady, real &gradz) const
 
 SphericalHarmonic ()
 
 SphericalHarmonic (const std::vector< real > &C, const std::vector< real > &S, int N, int nmx, int mmx, real a, unsigned norm=FULL)
 
 SphericalHarmonic (const std::vector< real > &C, const std::vector< real > &S, int N, real a, unsigned norm=FULL)
 

Private Types

typedef Math::real real
 

Private Attributes

real _a
 
SphericalEngine::coeff _c [1]
 
unsigned _norm
 

Detailed Description

Spherical harmonic series.

This class evaluates the spherical harmonic sum

V(x, y, z) = sum(n = 0..N)[ q^(n+1) * sum(m = 0..n)[
  (C[n,m] * cos(m*lambda) + S[n,m] * sin(m*lambda)) *
  P[n,m](cos(theta)) ] ]

where

Two normalizations are supported for Pnm

Clenshaw summation is used for the sums over both n and m. This allows the computation to be carried out without the need for any temporary arrays. See SphericalEngine.cpp for more information on the implementation.

References:

Example of use:

// Example of using the GeographicLib::SphericalHarmonic class
#include <iostream>
#include <exception>
#include <vector>
using namespace std;
using namespace GeographicLib;
int main() {
try {
int N = 3; // The maxium degree
double ca[] = {10, 9, 8, 7, 6, 5, 4, 3, 2, 1}; // cosine coefficients
vector<double> C(ca, ca + (N + 1) * (N + 2) / 2);
double sa[] = {6, 5, 4, 3, 2, 1}; // sine coefficients
vector<double> S(sa, sa + N * (N + 1) / 2);
double a = 1;
double x = 2, y = 3, z = 1;
double v, vx, vy, vz;
v = h(x, y, z, vx, vy, vz);
cout << v << " " << vx << " " << vy << " " << vz << "\n";
}
catch (const exception& e) {
cerr << "Caught exception: " << e.what() << "\n";
return 1;
}
}

Definition at line 69 of file SphericalHarmonic.hpp.

Member Typedef Documentation

◆ real

Definition at line 122 of file SphericalHarmonic.hpp.

Member Enumeration Documentation

◆ normalization

Supported normalizations for the associated Legendre polynomials.

Enumerator
FULL 

Fully normalized associated Legendre polynomials.

These are defined by Pnmfull(z) = (−1)m sqrt(k (2n + 1) (nm)! / (n + m)!) Pnm(z), where Pnm(z) is Ferrers function (also known as the Legendre function on the cut or the associated Legendre polynomial) http://dlmf.nist.gov/14.7.E10 and k = 1 for m = 0 and k = 2 otherwise.

The mean squared value of Pnmfull(cosθ) cos(mλ) and Pnmfull(cosθ) sin(mλ) over the sphere is 1.

SCHMIDT 

Schmidt semi-normalized associated Legendre polynomials.

These are defined by Pnmschmidt(z) = (−1)m sqrt(k (nm)! / (n + m)!) Pnm(z), where Pnm(z) is Ferrers function (also known as the Legendre function on the cut or the associated Legendre polynomial) http://dlmf.nist.gov/14.7.E10 and k = 1 for m = 0 and k = 2 otherwise.

The mean squared value of Pnmschmidt(cosθ) cos(mλ) and Pnmschmidt(cosθ) sin(mλ) over the sphere is 1/(2n + 1).

Definition at line 74 of file SphericalHarmonic.hpp.

Constructor & Destructor Documentation

◆ SphericalHarmonic() [1/3]

GeographicLib::SphericalHarmonic::SphericalHarmonic ( const std::vector< real > &  C,
const std::vector< real > &  S,
int  N,
real  a,
unsigned  norm = FULL 
)
inline

Constructor with a full set of coefficients specified.

Parameters
[in]Cthe coefficients Cnm.
[in]Sthe coefficients Snm.
[in]Nthe maximum degree and order of the sum
[in]athe reference radius appearing in the definition of the sum.
[in]normthe normalization for the associated Legendre polynomials, either SphericalHarmonic::FULL (the default) or SphericalHarmonic::SCHMIDT.
Exceptions
GeographicErrif N does not satisfy N ≥ −1.
GeographicErrif C or S is not big enough to hold the coefficients.

The coefficients Cnm and Snm are stored in the one-dimensional vectors C and S which must contain (N + 1)(N + 2)/2 and N (N + 1)/2 elements, respectively, stored in "column-major" order. Thus for N = 3, the order would be: C00, C10, C20, C30, C11, C21, C31, C22, C32, C33. In general the (n,m) element is at index m Nm (m − 1)/2 + n. The layout of S is the same except that the first column is omitted (since the m = 0 terms never contribute to the sum) and the 0th element is S11

The class stores pointers to the first elements of C and S. These arrays should not be altered or destroyed during the lifetime of a SphericalHarmonic object.

Definition at line 167 of file SphericalHarmonic.hpp.

◆ SphericalHarmonic() [2/3]

GeographicLib::SphericalHarmonic::SphericalHarmonic ( const std::vector< real > &  C,
const std::vector< real > &  S,
int  N,
int  nmx,
int  mmx,
real  a,
unsigned  norm = FULL 
)
inline

Constructor with a subset of coefficients specified.

Parameters
[in]Cthe coefficients Cnm.
[in]Sthe coefficients Snm.
[in]Nthe degree used to determine the layout of C and S.
[in]nmxthe maximum degree used in the sum. The sum over n is from 0 thru nmx.
[in]mmxthe maximum order used in the sum. The sum over m is from 0 thru min(n, mmx).
[in]athe reference radius appearing in the definition of the sum.
[in]normthe normalization for the associated Legendre polynomials, either SphericalHarmonic::FULL (the default) or SphericalHarmonic::SCHMIDT.
Exceptions
GeographicErrif N, nmx, and mmx do not satisfy Nnmxmmx ≥ −1.
GeographicErrif C or S is not big enough to hold the coefficients.

The class stores pointers to the first elements of C and S. These arrays should not be altered or destroyed during the lifetime of a SphericalHarmonic object.

Definition at line 198 of file SphericalHarmonic.hpp.

◆ SphericalHarmonic() [3/3]

GeographicLib::SphericalHarmonic::SphericalHarmonic ( )
inline

A default constructor so that the object can be created when the constructor for another object is initialized. This default object can then be reset with the default copy assignment operator.

Definition at line 211 of file SphericalHarmonic.hpp.

Member Function Documentation

◆ Circle()

CircularEngine GeographicLib::SphericalHarmonic::Circle ( real  p,
real  z,
bool  gradp 
) const
inline

Create a CircularEngine to allow the efficient evaluation of several points on a circle of latitude.

Parameters
[in]pthe radius of the circle.
[in]zthe height of the circle above the equatorial plane.
[in]gradpif true the returned object will be able to compute the gradient of the sum.
Exceptions
std::bad_allocif the memory for the CircularEngine can't be allocated.
Returns
the CircularEngine object.

SphericalHarmonic::operator()() exchanges the order of the sums in the definition, i.e., ∑n = 0..Nm = 0..n becomes ∑m = 0..Nn = m..N. SphericalHarmonic::Circle performs the inner sum over degree n (which entails about N2 operations). Calling CircularEngine::operator()() on the returned object performs the outer sum over the order m (about N operations).

Here's an example of computing the spherical sum at a sequence of longitudes without using a CircularEngine object

SphericalHarmonic h(...); // Create the SphericalHarmonic object
double r = 2, lat = 33, lon0 = 44, dlon = 0.01;
double
phi = lat * Math::degree<double>(),
z = r * sin(phi), p = r * cos(phi);
for (int i = 0; i <= 100; ++i) {
lon = lon0 + i * dlon,
lam = lon * Math::degree<double>();
std::cout << lon << " " << h(p * cos(lam), p * sin(lam), z) << "\n";
}

Here is the same calculation done using a CircularEngine object. This will be about N/2 times faster.

SphericalHarmonic h(...); // Create the SphericalHarmonic object
double r = 2, lat = 33, lon0 = 44, dlon = 0.01;
double
phi = lat * Math::degree<double>(),
z = r * sin(phi), p = r * cos(phi);
CircularEngine c(h(p, z, false)); // Create the CircularEngine object
for (int i = 0; i <= 100; ++i) {
lon = lon0 + i * dlon;
std::cout << lon << " " << c(lon) << "\n";
}

Definition at line 324 of file SphericalHarmonic.hpp.

◆ Coefficients()

const SphericalEngine::coeff& GeographicLib::SphericalHarmonic::Coefficients ( ) const
inline
Returns
the zeroth SphericalEngine::coeff object.

Definition at line 348 of file SphericalHarmonic.hpp.

◆ operator()() [1/2]

Math::real GeographicLib::SphericalHarmonic::operator() ( real  x,
real  y,
real  z 
) const
inline

Compute the spherical harmonic sum.

Parameters
[in]xcartesian coordinate.
[in]ycartesian coordinate.
[in]zcartesian coordinate.
Returns
V the spherical harmonic sum.

This routine requires constant memory and thus never throws an exception.

Definition at line 224 of file SphericalHarmonic.hpp.

◆ operator()() [2/2]

Math::real GeographicLib::SphericalHarmonic::operator() ( real  x,
real  y,
real  z,
real gradx,
real grady,
real gradz 
) const
inline

Compute a spherical harmonic sum and its gradient.

Parameters
[in]xcartesian coordinate.
[in]ycartesian coordinate.
[in]zcartesian coordinate.
[out]gradxx component of the gradient
[out]gradyy component of the gradient
[out]gradzz component of the gradient
Returns
V the spherical harmonic sum.

This is the same as the previous function, except that the components of the gradients of the sum in the x, y, and z directions are computed. This routine requires constant memory and thus never throws an exception.

Definition at line 257 of file SphericalHarmonic.hpp.

Member Data Documentation

◆ _a

real GeographicLib::SphericalHarmonic::_a
private

Definition at line 124 of file SphericalHarmonic.hpp.

◆ _c

SphericalEngine::coeff GeographicLib::SphericalHarmonic::_c[1]
private

Definition at line 123 of file SphericalHarmonic.hpp.

◆ _norm

unsigned GeographicLib::SphericalHarmonic::_norm
private

Definition at line 125 of file SphericalHarmonic.hpp.


The documentation for this class was generated from the following file:
e
Array< double, 1, 3 > e(1./3., 0.5, 2.)
ceres::sin
Jet< T, N > sin(const Jet< T, N > &f)
Definition: jet.h:439
GeographicLib
Namespace for GeographicLib.
Definition: JacobiConformal.hpp:15
c
Scalar Scalar * c
Definition: benchVecAdd.cpp:17
x
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition: gnuplot_common_settings.hh:12
GeographicLib::SphericalHarmonic::SphericalHarmonic
SphericalHarmonic()
Definition: SphericalHarmonic.hpp:211
h
const double h
Definition: testSimpleHelicopter.cpp:19
main
int main(int argc, char **argv)
Definition: cmake/example_cmake_find_gtsam/main.cpp:63
ceres::cos
Jet< T, N > cos(const Jet< T, N > &f)
Definition: jet.h:426
vy
StridedVectorType vy(make_vector(y, *n, std::abs(*incy)))
example::lon0
const double lon0
Definition: testGPSFactor.cpp:41
GeographicLib::SphericalHarmonic
Spherical harmonic series.
Definition: SphericalHarmonic.hpp:69
pybind_wrapper_test_script.z
z
Definition: pybind_wrapper_test_script.py:61
y
Scalar * y
Definition: level1_cplx_impl.h:124
a
ArrayXXi a
Definition: Array_initializer_list_23_cxx11.cpp:1
C
Matrix< Scalar, Dynamic, Dynamic > C
Definition: bench_gemm.cpp:50
std
Definition: BFloat16.h:88
p
float * p
Definition: Tutorial_Map_using.cpp:9
v
Array< int, Dynamic, 1 > v
Definition: Array_initializer_list_vector_cxx11.cpp:1
N
#define N
Definition: igam.h:9
lon
static const double lon
Definition: testGeographicLib.cpp:34
real
Definition: main.h:100
i
int i
Definition: BiCGSTAB_step_by_step.cpp:9
S
DiscreteKey S(1, 2)
SphericalHarmonic.hpp
Header for GeographicLib::SphericalHarmonic class.
vx
StridedVectorType vx(make_vector(x, *n, std::abs(*incx)))
lat
static const double lat
Definition: testGeographicLib.cpp:34


gtsam
Author(s):
autogenerated on Sat Nov 16 2024 04:15:07