Basis Functions
-
basis.
eqvlat_fawa
(ylat, vort, area, n_points, planet_radius=6378000.0, output_fawa=True) → Tuple[numpy.ndarray, Optional[numpy.ndarray]][source] Compute equivalent latitude and finite amplitude wave activity with the full
- Parameters
ylat (sequence or array_like) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.
vort (ndarray) – 2-d numpy array of vorticity values; dimension = (nlat, nlon).
area (ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).
n_points (int) – Analysis resolution to calculate equivalent latitude.
planet_radius (float) – Radius of spherical planet of interest consistent with input area. Default: earth’s radius 6.378e+6
output_fawa (bool) – Whether to output FAWA. If True, return tuple would consist of (qref, fawa). Else, this function will return (qref, None).
- Returns
qref (ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).
fawa (ndarray or None) – 1-d finite-amplitude local wave activity
-
basis.
eqvlat_vgrad
(ylat, vort, area, n_points, planet_radius=6378000.0, vgrad=None)[source] Compute equivalent latitude, and optionally <…>_Q in Nakamura and Zhu (2010).
- Parameters
ylat (sequence or array_like) – 1-d numpy array of latitude (in degree) with equal spacing in ascending order; dimension = nlat.
vort (ndarray) – 2-d numpy array of vorticity values; dimension = (nlat, nlon).
area (ndarray) – 2-d numpy array specifying differential areal element of each grid point; dimension = (nlat, nlon).
n_points (int) – Analysis resolution to calculate equivalent latitude.
planet_radius (float) – Radius of spherical planet of interest consistent with input area. Default: earth’s radius 6.378e+6
vgrad (ndarray, optional) – 2-d numpy array of laplacian (or higher-order laplacian) values; dimension = (nlat, nlon)
- Returns
q_part (ndarray) – 1-d numpy array of value Q(y) where latitude y is given by ylat; dimension = (nlat).
brac (ndarray or None) – 1-d numpy array of averaged vgrad in the square bracket. If vgrad = None, brac = None.
-
basis.
lwa
(nlon, nlat, vort, q_part, dmu, ncforce=None)[source] At each grid point of vorticity q(x,y) and reference state vorticity Q(y), this function calculate the difference between the line integral of [q(x,y+y’)-Q(y)] (and ncforce, if given) over the domain {y+y’>y,q(x,y+y’)<Q(y)} and {y+y’<y,q(x,y+y’)>Q(y)}. See fig. (1) and equation (13) of Huang and Nakamura (2016). dmu is a vector of length nlat: dmu = cos(phi) d(phi) such that phi is the latitude.
- Parameters
nlon (int) – Longitudinal dimension of vort (i.e. vort.shape[1]).
nlat (int) – Latitudinal dimension of vort (i.e. vort.shape[0]).
vort (ndarray) – 2-d numpy array of vorticity values; dimension = (nlat,nlon).
Q_part (sequence or array_like) – 1-d numpy array of Q (vorticity reference state) as a function of latitude. Size = nlat.
dmu (sequence or array_like) – 1-d numpy array of latitudinal differential length element (e.g. dmu = planet_radius * cos(lat) d(lat)). Size = nlat.
ncforce (ndarray or None, optional) – 2-d numpy array of non-conservative force field (i.e. theta in NZ10(a) in equation (23a) and (23b))
- Returns
lwact (ndarray) – 2-d numpy array of local wave activity calculated; dimension = (nlat,nlon).
bigsigma (ndarray or None) – 2-d numpy array of the nonconservative forces acting on local wave activity computed from ncforce. If ncforce = None, bigsigma = None.