Package spharm :: Module spharm :: Class Spharmt

Class Spharmt

source code

spherical harmonic transform class.

Instance Methods
 
__delattr__(self, key)
prevent deletion of read-only instance variables.
source code
 
__init__(self, nlon, nlat, rsphere=6371200.0, gridtype='regular', legfunc='stored')
create a Spharmt class instance.
source code
 
__setattr__(self, key, val)
prevent modification of read-only instance variables.
source code
 
getgrad(self, chispec)
compute vector gradient on grid given complex spectral coefficients.
source code
 
getpsichi(self, ugrid, vgrid, ntrunc=None)
compute streamfunction and velocity potential on grid given vector wind.
source code
 
getuv(self, vrtspec, divspec)
compute vector wind on grid given complex spectral coefficients of vorticity and divergence.
source code
 
getvrtdivspec(self, ugrid, vgrid, ntrunc=None)
compute spectral coefficients of vorticity and divergence given vector wind.
source code
 
grdtospec(self, datagrid, ntrunc=None)
grid to spectral transform (spherical harmonic analysis).
source code
 
specsmooth(self, datagrid, smooth)
isotropic spectral smoothing on a sphere.
source code
 
spectogrd(self, dataspec)
spectral to grid transform (spherical harmonic synthesis).
source code
Instance Variables
  gridtype
'regular' (equally spaced in longitude and latitude) or 'gaussian' (equally spaced in longitude, latitudes located at roots of ordinary Legendre polynomial of degree nlat).
  legfunc
'stored' or 'computed'.
  nlat
number of latitudes (set when class instance is created, cannot be changed).
  nlon
number of longitudes (set when class instance is created, cannot be changed).
  rsphere
The radius of the sphere in meters (set when class instance is created, cannot be changed).
Method Details

__init__(self, nlon, nlat, rsphere=6371200.0, gridtype='regular', legfunc='stored')
(Constructor)

source code 

create a Spharmt class instance.

Parameters:
  • nlon - Number of longitudes. The grid must be oriented from east to west, with the first point at the Greenwich meridian and the last point at 360-delta degrees east (where delta = 360/nlon degrees). Must be >= 4. Transforms will be faster when nlon is the product of small primes.
  • nlat - Number of latitudes. The grid must be oriented from north to south. If nlat is odd the equator is included. If nlat is even the equator will lie half way between points points nlat/2 and (nlat/2)+1. Must be >=3.
  • rsphere - The radius of the sphere in meters. Default 6371200 (the value for Earth).
  • legfunc - 'stored' (default) or 'computed'. If 'stored', associated legendre functions are precomputed and stored when the class instance is created. This uses O(nlat**3) memory, but speeds up the spectral transforms. If 'computed', associated legendre functions are computed on the fly when transforms are requested. This uses O(nlat**2) memory, but slows down the spectral transforms a bit.
  • gridtype - 'regular' (default) or 'gaussian'. Regular grids will include the poles and equator if nlat is odd. Gaussian grids never include the poles, but will include the equator if nlat is odd.

getgrad(self, chispec)

source code 

compute vector gradient on grid given complex spectral coefficients.

Parameters:
  • chispec - rank 1 or 2 numpy complex array with shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) containing complex spherical harmonic coefficients (where ntrunc is the triangular truncation limit and nt is the number of spectral arrays to be transformed). If chispec is rank 1, nt is assumed to be 1.
Returns:
uchi, vchi - rank 2 or 3 numpy float32 arrays containing gridded zonal and meridional components of the vector gradient. Shapes are either (nlat,nlon) or (nlat,nlon,nt).

getpsichi(self, ugrid, vgrid, ntrunc=None)

source code 

compute streamfunction and velocity potential on grid given vector wind.

Parameters:
  • ugrid - rank 2 or 3 numpy float32 array containing grid of zonal winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number of grids to be transformed. If ugrid is rank 2, nt is assumed to be 1.
  • vgrid - rank 2 or 3 numpy float32 array containing grid of meridional winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number of grids to be transformed. Both ugrid and vgrid must have the same shape.
  • ntrunc - optional spectral truncation limit. (default self.nlat-1)
Returns:
psigrid, chigrid - rank 2 or 3 numpy float32 arrays of gridded streamfunction and velocity potential. Shapes are either (nlat,nlon) or (nlat,nlon,nt).

getuv(self, vrtspec, divspec)

source code 

compute vector wind on grid given complex spectral coefficients of vorticity and divergence.

Parameters:
  • vrtspec - rank 1 or 2 numpy complex array of vorticity spectral coefficients, with shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) (where ntrunc is the triangular truncation and nt is the number of spectral arrays to be transformed). If vrtspec is rank 1, nt is assumed to be 1.
  • divspec - rank 1 or 2 numpy complex array of divergence spectral coefficients, with shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) (where ntrunc is the triangular truncation and nt is the number of spectral arrays to be transformed). Both vrtspec and divspec must have the same shape.
Returns:
ugrid, vgrid - rank 2 or 3 numpy float32 arrays containing gridded zonal and meridional winds. Shapes are either (nlat,nlon) or (nlat,nlon,nt).

getvrtdivspec(self, ugrid, vgrid, ntrunc=None)

source code 

compute spectral coefficients of vorticity and divergence given vector wind.

Parameters:
  • ugrid - rank 2 or 3 numpy float32 array containing grid of zonal winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number of grids to be transformed. If ugrid is rank 2, nt is assumed to be 1.
  • vgrid - rank 2 or 3 numpy float32 array containing grid of meridional winds. Must have shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number of grids to be transformed. Both ugrid and vgrid must have the same shape.
  • ntrunc - optional spectral truncation limit. (default self.nlat-1)
Returns:
vrtspec, divspec - rank 1 or 2 numpy complex arrays of vorticity and divergence spherical harmonic coefficients with shape shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt).

grdtospec(self, datagrid, ntrunc=None)

source code 

grid to spectral transform (spherical harmonic analysis).

Parameters:
  • datagrid - rank 2 or 3 numpy float32 array with shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number of grids to be transformed. If datagrid is rank 2, nt is assumed to be 1.
  • ntrunc - optional spectral truncation limit. (default self.nlat-1)
Returns:
dataspec - rank 1 or 2 numpy complex array with shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) containing complex spherical harmonic coefficients resulting from the spherical harmonic analysis of datagrid.

specsmooth(self, datagrid, smooth)

source code 

isotropic spectral smoothing on a sphere.

Parameters:
  • datagrid - rank 2 or 3 numpy float32 array with shape (nlat,nlon) or (nlat,nlon,nt), where nt is the number of grids to be smoothed. If datagrid is rank 2, nt is assumed to be 1.
  • smooth - rank 1 array of length nlat containing smoothing factors as a function of total wavenumber.
Returns:
datagrid - rank 2 or 3 numpy float32 array with shape (nlat,nlon) or (nlat,nlon,nt) containing the smoothed grids.

spectogrd(self, dataspec)

source code 

spectral to grid transform (spherical harmonic synthesis).

Parameters:
  • dataspec - rank 1 or 2 numpy complex array with shape (ntrunc+1)*(ntrunc+2)/2 or ((ntrunc+1)*(ntrunc+2)/2,nt) containing complex spherical harmonic coefficients (where ntrunc is the triangular truncation limit and nt is the number of spectral arrays to be transformed). If dataspec is rank 1, nt is assumed to be 1.
Returns:
datagrid - rank 2 or 3 numpy float32 array with shape (nlat,nlon) or (nlat,nlon,nt) containing the gridded data resulting from the spherical harmonic synthesis of dataspec.

Instance Variable Details

gridtype

'regular' (equally spaced in longitude and latitude) or 'gaussian' (equally spaced in longitude, latitudes located at roots of ordinary Legendre polynomial of degree nlat). Set when class instance is created, cannot be changed.

legfunc

'stored' or 'computed'. If 'stored', associated legendre functions are precomputed and stored when the class instance is created. If 'computed', associated legendre functions are computed on the fly when transforms are requested. Set when class instance is created, cannot be changed.