HOME | DOWNLOAD | DOCUMENTATION | FAQ |
> Home > Documentation > Python > Spherical Harmonic Transformations
GLQGridCoord
Compute the latitude and longitude coordinates used in Gauss-Legendre quadrature grids.
Usage
latglq
, longlq
= pyshtools.GLQGridCoord (lmax
)
Returns
latglq
: float, dimension (lmax
+1)- The latitude coordinates of a grid, corresponding to the indices [:,0], in DEGREES.
longlq
: float, dimension (2*lmax
+1)- The longitude coordinates of a grid, corresponding to the indices [0,:], in DEGREES. The first node is 0 E.
Parameters
lmax
: integer- The maximum spherical harmonic degree that will be integrated exactly by Gauss-Legendre quadrature.
Description
GLQGridCoord
will compute the latitude and longitude coordinates that are used in Gauss-Legendre quadrature grids for performing spherical harmonic transforms and reconstructions. The latitudinal nodes correspond to the zeros of the Legendre polynomial of degree lmax+1
, and the longitudinal nodes are equally spaced with an interval of 360/(2*lmax+1)
degrees.
See also
shglq
, shexpandglq
, makegridglq
, shexpandglqc
, makegridglqc
, preglq
> Home > Documentation > Python > Spherical Harmonic Transformations
Institut de Physique du Globe de Paris | University of Sorbonne Paris Cité | © 2016 SHTOOLS |