Diagnostics

Waveguide Diagnostics

barotropic.diagnostics.extract_ks_waveguides(ks_or_state, k, grid=None)

Extract waveguide boundaries based on the stationary wavenumber.

Parameters:
  • ks_or_state (State | array) – Input state or stationary wavenumber profile (complex variant).

  • k (int) – Wavenumber for which waveguides are extracted.

  • grid (None | Grid) – Grid information only required if ks_or_state is not a State.

Returns:

List of tuples that contain the boundary latitudes of the detected waveguides (one tuple per waveguide) in degrees.

WKB and ray-tracing theory say that Rossby waves are refracted towards latitudes of higher stationary wavenumber. A local maximum in the meridional profile of stationary zonal wavenumber Ks can therefore trap waves and constitutes a waveguide. See e.g. Petoukhov et al. (2013), Wirth (2020).

barotropic.diagnostics.stationary_wavenumber(u_or_state, grid=None, order=4, min_u=0.001, kind='complex')

Non-dimensionalised stationary (zonal) wavenumber.

Parameters:
  • u_or_state (State | array) – Input state or 2D field of zonal wind component.

  • grid (None | Grid) – Grid information only required if u_or_state is not a State.

  • order (2 | 4) – Order of the derivative used in the calculation. Parameter of Grid.derivative_meridional().

  • min_u (number) – Threshold for the absolute value of the zonal wind below which it is considered to be zero.

  • kind ("complex" | "real" | "squared") –

    Ks² can be negative. Determine how Ks is returned:

    • ”complex”: return complex number √(Ks²).

    • ”real”: return real number Re(Ks) - Im(Ks) where Ks is the complex number √(Ks²). Since Ks² is real, Im(Ks) is always zero when Re(Ks) is non-zero and vice versa, rendering this expression of Ks unique.

    • ”squared”: return Ks².

Returns:

Meridional profile of stationary wavenumber.

Stationary zonal wavenumber:

Ks² = r²·βM/uM,

where r is the radius of the sphere, βM is the Mercator projection gradient of zonal-mean PV and uM is the Mercator projection zonal-mean wind. See e.g. Hoskins and Karoly (1981), Wirth (2020).

Wave Activity Diagnostics

barotropic.diagnostics.falwa(pv_or_state, grid=None, levels=None, interpolate=True, quad='sptrapz')

Finite-Amplitude Local Wave Activity according to Huang and Nakamura (2016).

Parameters:
  • pv_or_state (State | array) – Input state or 2D PV field.

  • grid (None | Grid) – Grid information only required if pv_or_state is not a State.

  • levels (None | array | int) – Parameter of Grid.zonalize().

  • interpolate (bool) – Interpolate output onto the regular latitudes of the grid. If False, FALWA is returned on the equivalent latitudes obtained from the zonalization procedure. These may be irregular and unordered.

  • quad (str) – Quadrature rule used in the surface integrals of the computation. See Grid.quad().

Returns:

If interpolate is True, return FALWA interpolated onto the regular latitudes as an array. Otherwise, return a tuple of FALWA and the associated equivalent latitudes (both arrays).

To improve numerical stability if a boxcounting quadrature is selected, the zonalized PV profile is always interpolated to the dual latitude grid. Therefore, levels does not affect the output latitudes, even if interpolate is set to False.

barotropic.diagnostics.fawa(pv_or_state, grid=None, levels=None, interpolate=True, quad='sptrapz')

Finite-Amplitude Wave Activity according to Nakamura and Zhu (2010).

Computes FAWA as the zonal mean of falwa().

Wave Packet Diagnostics

barotropic.diagnostics.envelope_hilbert(field, wavenumber_range=(2, 10))

Compute the envelope of wave packets with the Hilbert transform.

Parameters:
  • field (array) – 2D input field.

  • wavenumber_range ((number, number)) – Restrict the interval in which the dominant wavenumber is determined.

Returns:

Filtered output field.

Applied to the merdional wind for RWP detection by Zimin et al. (2003).

Other

barotropic.diagnostics.dominant_wavenumber_fourier(field, grid, smooth=('boxcar', 10), wavenumber_range=(1, 20))

Dominant zonal wavenumber based on the Fourier power spectrum.

Parameters:
  • field (array) – 2D input field.

  • grid (Grid) – Grid associated with input field.

  • smooth (None | tuple) – If not None, the output wavenumber profile is filtered with Grid.filter_meridional(), whose window and width arguments are set by smooth.

  • wavenumber_range ((number, number)) – Restrict the interval in which the dominant wavenumber is determined.

Returns:

Returns the dominant zonal wavenumber as a function of latitude.

barotropic.diagnostics.dominant_wavenumber_wavelet(field, grid, smooth=('hann', 40, 10))

Dominant zonal wavenumber based on Wavelet Analysis.

Parameters:
  • field (array) – 2D input field.

  • grid (Grid) – Grid associated with input field.

  • smooth (None | tuple) – The smoothing applied to the dominant wavenumber field. Set to None if no smoothing is desired. Otherwise provide a 3-tuple (window, width_lon, width_lat) used as input to Grid.filter_meridional() and Grid.filter_zonal().

Returns:

Gridded dominant zonal wavenumber.

Implements the procedure of Ghinassi et al. (2018) that determines a local dominant wavenumber at every gridpoint of the input field.

Requires pywt (version >= 1.1.0).

barotropic.diagnostics.filter_by_wavenumber(field, wavenumber)

Zonal Hann smoothing based on a space-dependent wavenumber field.

Parameters:
  • field (array) – 2D input field.

  • wavenumber (array) – nlat-sized 1D array or field-shaped 2D array of wavenumbers.

Returns:

Filtered field.

Apply Hann smoothing in zonal direction with the Hann window width governed at every gridpoint by the wavenumber given for the same location. The wavelength corresponding to each wavenumber determines the full width at half maximum of the Hann window.

Used by Ghinassi et al. (2018) to filter the FALWA field, discounting phase information in order to diagnose wave packets as a whole.