Time stepping

This module contains the function that calls an explicit or an implicit time-stepping scheme, as well as the functions that compute the arrays of coefficients involved in time-stepping.

SWAMPE.time_stepping.RMS_winds(a, I, J, lambdas, mus, U, V)

Computes RMS winds based on the zonal and meridional wind fields.

Parameters:
  • a (float) – planteray radius, m

  • I (int) – number of longitudes

  • J (int) – number of Gaussian latitudes

  • mus (array of float) – Array of Gaussian latitudes of length J

  • lambdas (array of float) – Array of uniformly spaces longitudes of length I.

  • U (array of float) – zonal wind field, JxI

  • V (array of float) – meridional wind field, JxI

Returns:

RMS winds value

Return type:

float

SWAMPE.time_stepping.mJarray(J, M)

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

Parameters:
  • J (int) – number of Gaussian latitudes

  • M (int) – spectral dimension

Return mJarray:

coefficient m in a matrix of size (J, M+1)

Return type:

array of float

SWAMPE.time_stepping.marray(M, N)

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

Parameters:
  • M (int) – highest wavenumber of associated Legendre polynomials

  • N (int) – highest degreee of associated Legendre polynomials

Return marray:

coefficient m in a matrix of size (M+1, N+1)

Return type:

array of float

SWAMPE.time_stepping.narray(M, N)

Computes the array n(n+1).

Parameters:
  • M (int) – spectral dimension

  • N (int) – highest degreee of associated Legendre polynomials

Return narray:

coefficients n(n+1) in a matrix of size (M+1,N+1)

Return type:

array of float

SWAMPE.time_stepping.tstepcoeff(J, M, dt, mus, a)

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

Parameters:
  • J (int) – number of Gaussian latitudes

  • M (int) – spectral dimension

  • mus (array of float) – array of Gaussian latitudes of length J

  • a (float) – planetary radius, m

Return tstepcoeff:

coefficients 2dt/(a(1-mus^2)) in a matrix of size (J,M+1)

Return type:

array of float

SWAMPE.time_stepping.tstepcoeff2(J, M, dt, a)

Computes the time stepping coefficient of the form 2dt/a^2 from Hack and Jakob (1992).

Parameters:
  • J (int) – number of Gaussian latitudes

  • M (int) – spectral dimension

  • dt (float) – time step length, s

  • a (float) – planetary radius, m

Return tstepcoeff2:

time stepping coefficients of size (J,M+1)

Return type:

array of float

SWAMPE.time_stepping.tstepcoeffmn(M, N, a)

Generates the coefficient that multiplies spectral components in invrsUV.

Parameters:
  • M (int) – highest wave number

  • N – highest degree of associated legendre polynomial for m=0

  • a (float) – radius of the planet, in meters

Returns:

an array of coefficients a/(n(n+1))

Return type:

list

SWAMPE.time_stepping.tstepping(etam0, etam1, deltam0, deltam1, Phim0, Phim1, I, J, M, N, Am, Bm, Cm, Dm, Em, Fm, Gm, Um, Vm, fmn, Pmn, Hmn, w, tstepcoeff, tstepcoeff2, tstepcoeffmn, marray, mJarray, narray, PhiFm, dt, a, Phibar, taurad, taudrag, forcflag, diffflag, expflag, sigma, sigmaPhi, test, t)

Calls the timestepping scheme.

Parameters:
  • etam0 (array of float) – Fourier coefficents of absolute vorticity for one time step

  • etam1 (array of float) – Fourier coefficents of absolute vorticity for the following time step

  • deltam0 (array of float) – Fourier coefficents of divergence for one time step

  • deltam1 (array of float) – Fourier coefficents of divergence for the following time step

  • Phim0 (array of float) – Fourier coefficents of geopotential for one time step

  • Phim1 (array of float) – Fourier coefficents of geopotential for the following time step

  • I (int) – number of longitudes

  • J (int) – number of Gaussian latitudes

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

  • Am (array of float) – Fourier coefficients of the nonlinear component A=U*eta

  • Bm (array of float) – Fourier coefficients of the nonlinear component B=V*eta

  • Cm (array of float) – Fourier coefficients of the nonlinear component C=U*Phi

  • Dm (array of float) – Fourier coefficients of the nonlinear component D=V*Phi

  • Em (array of float) – Fourier coefficients of the nonlinear component E=(U^2+V^2)/(2(1-mu^2))

  • Fm (array of float) – Fourier coefficients of the zonal component of wind forcing

  • Gm (array of float) – Fourier coefficients of the meridional component of wind forcing

  • Um (array of float) – Fourier coefficients of the zonal component of wind

  • Vm (array of float) – Fourier coefficients of the meridional component of wind

  • fmn (array of float) – spectral coefficients of the Coriolis force

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

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

  • w (array of float) – Gauss Legendre weights

  • tstepcoeff (array of float) – coefficient for time-stepping of the form 2dt/(a(1-mus^2))

  • tstepcoeff2 (array of float) – time stepping coefficient of the form 2dt/a^2

  • tstepcoeffmn (array of float) – an array of coefficients a/(n(n+1))

  • marray (array of float) – coefficients equal to m=0,1,…,M in a matrix M+1xN+1

  • mJarray (array of float) – coefficients equal to m=0,1,…,M in a matrix M+1xJ

  • narray (array of float) – array n(n+1) in a matrix M+1xN+1

  • PhiFm – Fourier coefficients of the geopotential forcing

  • dt (float) – time step, in seconds

  • a (float) – planetary radius, m

  • Phibar (float) – time-invariant spatial mean geopotential, height of the top layer

  • taurad (float) – radiative timescale

  • taudrag (float) – drag timescale

  • forcflag (float) – forcing flag

  • diffflag (float) – hyperdiffusion filter flag

  • sigma (array of float) – hyperdiffusion filter coefficients for absolute vorticity and divergence

  • sigmaPhi (array of float) – hyperdiffusion filter coefficients for geopotential

  • test (int) – number of test, defaults to None

  • t (int) – number of current time step

Returns:

  • newetamn

    Updated spectral coefficients of absolute vorticity

  • newetatstep

    Updated absolute vorticity

  • newdeltamn

    Updated spectral coefficients of divergence

  • newdeltatstep

    Updated divergence

  • newPhimn

    Updated spectral coefficients of geopotential

  • newPhitstep

    Updated geopotential

  • Unew

    Updated zonal winds

  • Vnew

    Updated meridional winds

Return type:

array of float