Spherical Harmonic Transforms

SWAMPE.spectral_transform.PmnHmn(J, M, N, mus)

Calculates the values of associated Legendre polynomials and their derivatives evaluated at Gaussian latitudes (mus) up to wavenumber M.

Parameters
  • J (int) – number of latitudes

  • M (int) – highest wavenumber for associated Legendre polynomials

  • N (int) – highest degree of the Legendre functions for m=0

  • mus (array of floats) – Gaussian latitudes

Returns

Pmn array (associated legendre polynomials), Hmn array (derivatives of Pmn*(1-x^2)), both evaluated at the Gaussian latitudes mus

Return type

array of float

SWAMPE.spectral_transform.diagnostic_eta_delta(Um, Vm, fmn, I, J, M, N, Pmn, Hmn, w, tstepcoeff, mJarray, dt)

Computes vorticity and divergence from zonal and meridional wind fields. This is a diagnostic relationship. For details, see Hack and Jakob (1992) equations (5.26)-(5.27).

Parameters

param Um

Fourier coefficient of zonal winds

type Um

array of float

param Vm

Fourier coefficient of meridional winds

type Vm

array of float

param fmn

spectal coefficients of the Coriolic force

type fmn

array of float

param I

number of longitudes

type I

int

param J

number of Gaussian latitudes

type J

int

param M

highest wavenumber for associated Legendre polynomials

type M

int

param N

highest degree of associated Legendre polynomials

type N

int

param Pmn

values of the associated Legendre polynomials at Gaussian latitudes mus up to wavenumber M

type Pmn

array of float64

param Hmn

values of the associated Legendre polynomial derivatives at Gaussian latitudes up to wavenumber M

type Hmn

array of float

param w

Gauss Legendre weights

type w

array of float

param tstepcoeff

a coefficient for time-stepping of the form 2dt/(a(1-mus^2) from Hack and Jakob (1992)

type tstepcoeff

array of float

param mJarray

coefficients equal to m=0,1,…,M

type mJarray

array of float

param dt

time step, in seconds

type dt

float Returns

Absolute vorticity

  • newdelta

    Divergence

  • etamn

    Spectral coefficients of absolute vorticity

  • deltamn

    Spectral coefficients of divergence

rtype

array of float

SWAMPE.spectral_transform.fwd_fft_trunc(data, I, M)

Calculates and truncates the fast forward Fourier transform of the input.

Parameters
  • data (array of float) – array of dimension IxJ (usually the values of state variables at lat-long coordinates)

  • I (int) – number of longitudes

  • M (int) – highest wavenumber for associated Legendre polynomials

Return datam

Fourier coefficients 0 through M

Rtype datam

array of complex

SWAMPE.spectral_transform.fwd_leg(data, J, M, N, Pmn, w)

Calculates the forward legendre transform.

Parameters
  • data (array of float or array of complex) – input to be transformed (usually output of fft)

  • J (int) – number of latitudes

  • M (int) – highest wavenumber for associated Legendre polynomials

  • N (int) – highest degree of the Legendre functions for m=0

  • Pmn (array of float) – associated legendre functions evaluated at the Gaussian latitudes mus up to wavenumber M

  • w (array of float) – Gauss Legendre weights

Return legcoeff

Legendre coefficients (if data was output from FFT, then legcoeff are the spectral coefficients)

Rtype legcoeff

array of complex

SWAMPE.spectral_transform.invrsUV(deltamn, etamn, fmn, I, J, M, N, Pmn, Hmn, tstepcoeffmn, marray)

Computes the wind velocity from the values of vorticity and divergence. This is a diagnostic relationship. For details, see Hack and Jakob (1992) equations (5.24)-(5.25).

Parameters

param deltamn

Fourier coefficients of divergence

type deltamn

array of complex

param etamn

Fourier coefficients of vorticity

type etamn

array of complex

param fmn

spectral coefficients of the Coriolis force

type etamn

array of float

param I

number of longitudes

type I

int

param J

number of latitudes

type J

int

param M

highest wavenumber for associated Legendre polynomials

type M

int

param N

highest degree of associated Legendre polynomials

type N

int

param Pmn

values of the associated Legendre polynomials at Gaussian latitudes mus up to wavenumber M

type Pmn

array of float

param Hmn

values of the associated Legendre polynomial derivatives at Gaussian latitudes up to wavenumber M

type Hmn

array of float

param tstepcoeffmn

coefficient to scale spectral components

type tstepcoeffmn

array of float

param marray

array to multiply a quantity by a factor of m ranging from 0 through M.

type marray

array of float

Returns

return

  • Unew

    Zonal velocity component

  • Vnew

    Meridional velocity component

rtype

array of float

SWAMPE.spectral_transform.invrs_fft(approxXim, I)

Calculates the inverse Fourier transform.

Parameters
  • approxXim (array of complex) – Fourier coefficients

  • I (int) – number of longitudes

Returns

long-lat coefficients

Return type

array of complex

SWAMPE.spectral_transform.invrs_leg(legcoeff, I, J, M, N, Pmn)

Calculates the inverse Legendre transform.

Parameters
  • legcoeff (array of complex) – Legendre coefficients

  • J (int) – number of latitudes

  • M (int) – highest wavenumber for associated Legendre polynomials

  • Pmn (array of float) – associated legendre functions evaluated at the Gaussian latitudes mus up to wavenumber M

Returns

transformed spectral coefficients

Return type

array of complex