Package spharm :: Module spharm

Module spharm

source code

Introduction

This module provides a python interface to the NCAR SPHEREPACK library. It is not a one-to-one wrapper for the SPHEREPACK routines, rather it provides a simple interface oriented toward working with atmospheric general circulation model (GCM) data.

Requirements

Installation

Usage

>>> import spharm
>>> x=spharm.Spharmt(144,72,rsphere=8e6,gridtype='gaussian',legfunc='computed')

creates a class instance for spherical harmonic calculations on a 144x72 gaussian grid on a sphere with radius 8000 km. The associated legendre functions are recomputed on the fly (instead of pre-computed and stored). Default values of rsphere, gridtype and legfunc are 6.3712e6, 'regular' and 'stored'. Real-world examples are included in the source distribution.

Class methods

Functions

Conventions

The gridded data is assumed to be oriented such that i=1 is the Greenwich meridian and j=1 is the northernmost point. Grid indices increase eastward and southward. If nlat is odd the equator is included. If nlat is even the equator will lie half way between points nlat/2 and (nlat/2)+1. nlat must be at least 3. For regular grids (gridtype='regular') the poles will be included when nlat is odd. The grid increment in longitude is 2*pi/nlon radians. For example, nlon = 72 for a five degree grid. nlon must be greater than or equal to 4. The efficiency of the computation is improved when nlon is a product of small prime numbers.

The spectral data is assumed to be in a complex array of dimension (ntrunc+1)*(ntrunc+2)/2. ntrunc is the triangular truncation limit (ntrunc = 42 for T42). ntrunc must be <= nlat-1. Coefficients are ordered so that first (nm=0) is m=0,n=0, second is m=0,n=1, nm=ntrunc is m=0,n=ntrunc, nm=ntrunc+1 is m=1,n=1, etc. The values of m (degree) and n (order) as a function of the index nm are given by the arrays indxm, indxn returned by getspecindx.

The associated legendre polynomials are normalized so that the integral (pbar(n,m,theta)**2)*sin(theta) on the interval theta=0 to pi is 1, where pbar(m,n,theta)=sqrt((2*n+1)*factorial(n-m)/(2*factorial(n+m)))* sin(theta)**m/(2**n*factorial(n)) times the (n+m)th derivative of (x**2-1)**n with respect to x=cos(theta). theta = pi/2 - phi, where phi is latitude and theta is colatitude. Therefore, cos(theta) = sin(phi) and sin(theta) = cos(phi). Note that pbar(0,0,theta)=sqrt(2)/2, and pbar(1,0,theta)=.5*sqrt(6)*sin(lat).

The default grid type is regular (equally spaced latitude points). Set gridtype='gaussian' when creating a class instance for gaussian latitude points.

Quantities needed to compute spherical harmonics are precomputed and stored when the class instance is created with legfunc='stored' (the default). If legfunc='computed', they are recomputed on the fly on each method call. The storage requirements for legfunc="stored" increase like nlat**2, while those for legfunc='stored' increase like nlat**3. However, for repeated method invocations on a single class instance, legfunc="stored" will always be faster.


Contact: Jeff Whitaker

Version: 1.0.7

License: Permission to use, copy, modify, and distribute this software and its documentation for any purpose and without fee is hereby granted, provided that the above copyright notice appear in all copies and that both that copyright notice and this permission notice appear in supporting documentation. THE AUTHOR DISCLAIMS ALL WARRANTIES WITH REGARD TO THIS SOFTWARE, INCLUDING ALL IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS. IN NO EVENT SHALL THE AUTHOR BE LIABLE FOR ANY SPECIAL, INDIRECT OR CONSEQUENTIAL DAMAGES OR ANY DAMAGES WHATSOEVER RESULTING FROM LOSS OF USE, DATA OR PROFITS, WHETHER IN AN ACTION OF CONTRACT, NEGLIGENCE OR OTHER TORTIOUS ACTION, ARISING OUT OF OR IN CONNECTION WITH THE USE OR PERFORMANCE OF THIS SOFTWARE.

Classes
  Spharmt
spherical harmonic transform class.
Functions
 
gaussian_lats_wts(nlat)
compute the gaussian latitudes (in degrees) and quadrature weights.
source code
 
getgeodesicpts(m)
computes the lat/lon values of the points on the surface of the sphere corresponding to a twenty-sided (icosahedral) geodesic.
source code
 
getspecindx(ntrunc)
compute indices of zonal wavenumber (indxm) and degree (indxn) for complex spherical harmonic coefficients.
source code
 
legendre(lat, ntrunc)
calculate associated legendre functions for triangular truncation T(ntrunc), at a given latitude.
source code
 
regrid(grdin, grdout, datagrid, ntrunc=None, smooth=None)
regrid data using spectral interpolation, while performing optional spectral smoothing and/or truncation.
source code
 
specintrp(lon, dataspec, legfuncs)
spectral interpolation given spherical harmonic coefficients.
source code
Variables
  __package__ = 'spharm'
  __version__ = '1.0.9'
Function Details

gaussian_lats_wts(nlat)

source code 

compute the gaussian latitudes (in degrees) and quadrature weights.

Parameters:
  • nlat - number of gaussian latitudes desired.
Returns:
lats, wts - rank 1 numpy float64 arrays containing gaussian latitudes (in degrees north) and gaussian quadrature weights.

getgeodesicpts(m)

source code 

computes the lat/lon values of the points on the surface of the sphere corresponding to a twenty-sided (icosahedral) geodesic.

Parameters:
  • m - the number of points on the edge of a single geodesic triangle. There are 10*(m-1)**2+2 total geodesic points, including the poles.
Returns:
lats, lons - rank 1 numpy float32 arrays containing the latitudes and longitudes of the geodesic points (in degrees). These points are nearly evenly distributed on the surface of the sphere.

getspecindx(ntrunc)

source code 

compute indices of zonal wavenumber (indxm) and degree (indxn) for complex spherical harmonic coefficients.

Parameters:
  • ntrunc - spherical harmonic triangular truncation limit.
Returns:
indxm, indxn - rank 1 numpy Int32 arrays containing zonal wavenumber (indxm) and degree (indxn) of spherical harmonic coefficients.

legendre(lat, ntrunc)

source code 

calculate associated legendre functions for triangular truncation T(ntrunc), at a given latitude.

Parameters:
  • lat - the latitude (in degrees) to compute the associate legendre functions.
  • ntrunc - the triangular truncation limit.
Returns:
pnm - rank 1 numpy float32 array containing the (ntrunc+1)*(ntrunc+2)/2 associated legendre functions at latitude lat.

regrid(grdin, grdout, datagrid, ntrunc=None, smooth=None)

source code 

regrid data using spectral interpolation, while performing optional spectral smoothing and/or truncation.

Parameters:
  • grdin - Spharmt class instance describing input grid.
  • grdout - Spharmt class instance describing output grid.
  • datagrid - data on input grid (grdin.nlat x grdin.nlon). If datagrid is rank 3, last dimension is the number of grids to interpolate.
  • ntrunc - optional spectral truncation limit for datagrid (default min(grdin.nlat-1,grdout.nlat-1)).
  • smooth - rank 1 array of length grdout.nlat containing smoothing factors as a function of total wavenumber (default is no smoothing).
Returns:
datagrid - interpolated (and optionally smoothed) array(s) on grdout.nlon x grdout.nlat grid.

specintrp(lon, dataspec, legfuncs)

source code 

spectral interpolation given spherical harmonic coefficients.

Parameters:
  • lon - longitude (in degrees) of point on a sphere to interpolate to.
  • dataspec - spectral coefficients of function to interpolate.
  • legfuncs - associated legendre functions with same triangular truncation as dataspec (computed using legendre), computed at latitude of interpolation point.
Returns:
ob - interpolated value.