Zonalized Background States

Zonalization and rolling zonalization of potential vorticity (PV) on isentropes, based on the concept of equivalent latitude.

_images/rzon.png

Example

>>> from rwguide.xarray import zonalization

Start from a dataset of isentropic density and absolute vorticity on isentropic levels:

>>> isen
<xarray.Dataset>
Dimensions:    (isentrope: 3, latitude: 121, longitude: 240, time: 1460)
Coordinates:
  * isentrope  (isentrope) float64 340.0 330.0 320.0
  * latitude   (latitude) float32 90.0 88.5 87.0 85.5 ... -87.0 -88.5 -90.0
  * longitude  (longitude) float32 -180.0 -178.5 -177.0 ... 175.5 177.0 178.5
  * time       (time) datetime64[ns] 2018-01-01 ... 2018-12-31T18:00:00
Data variables:
    sg         (time, isentrope, latitude, longitude) float64 dask.array<...>
    av         (time, isentrope, latitude, longitude) float64 dask.array<...>

Create a weighting window for rolling zonalization:

>>> win = zonalization.fixed_deg_window(isen["longitude"], isen["latitude"], width=60., window="boxcar")
>>> win
<xarray.DataArray 'fixed_deg_window' (latitude: 121, longitude: 41)>
array(...)
Coordinates:
  * latitude   (latitude) float32 90.0 88.5 87.0 85.5 ... -87.0 -88.5 -90.0
  * longitude  (longitude) float32 0.0 1.5 3.0 4.5 6.0 ... 55.5 57.0 58.5 60.0
Attributes:
    width:    60.0

Compute a rolling-mean isentropic density background:

>>> bgsg = zonalization.rolling_mean_background(win, isen["sg"])
>>> bgsg
<xarray.DataArray 'sg_rm' (time: 1460, isentrope: 3, latitude: 121, longitude: 240)>
dask.array<...>
Coordinates:
  * latitude   (latitude) float32 90.0 88.5 87.0 85.5 ... -87.0 -88.5 -90.0
  * isentrope  (isentrope) float64 340.0 330.0 320.0
  * longitude  (longitude) float32 -180.0 -178.5 -177.0 ... 175.5 177.0 178.5
  * time       (time) datetime64[ns] 2018-01-01 ... 2018-12-31T18:00:00

Run the rolling zonalization:

>>> zonalization.zonalize_rolling(win, isen["av"], isen["sg"], bgsg)
<xarray.DataArray 'pv_rz' (time: 1460, isentrope: 3, latitude: 121, longitude: 240)>
dask.array<...>
Coordinates:
  * latitude   (latitude) float32 90.0 88.5 87.0 85.5 ... -87.0 -88.5 -90.0
  * isentrope  (isentrope) float64 340.0 330.0 320.0
  * longitude  (longitude) float32 -180.0 -178.5 -177.0 ... 175.5 177.0 178.5
  * time       (time) datetime64[ns] 2018-01-01 ... 2018-12-31T18:00:00

Xarray Interface

rwguide.xarray.zonalization.area_weights(da_lat, name='area_weights', **kwargs)

Area weights for a lat-lon-gridded sphere.

Added in version 1.2.

Parameters:
Returns:

Meridional profile of weights.

Return type:

xarray.DataArray

rwguide.xarray.zonalization.fixed_deg_window(da_lon, da_lat, name='fixed_deg_window', **kwargs)

Weighting window with a constant-longitude width (in degrees longitude).

Combines a latitude-dependent cosine area weighting and a longitudinally varying weighting window from scipy.signal.

The constant width in terms of degrees longitude should ensure a consistent cut-off wavenumber across the sphere.

Parameters:
  • da_lon (xarray.DataArray) – Longitude coordinate.

  • da_lat (xarray.DataArray) – Latitude coordinate.

  • name (str, optional) – Name assigned to the output.

  • kwargs (any) – Keyword arguments given to rwguide.zonalization.fixed_deg_window().

Returns:

Weighting window.

Return type:

xarray.DataArray

rwguide.xarray.zonalization.fixed_km_window(da_lon, da_lat, name='fixed_km_window', **kwargs)

Weighting window with a constant-distance width (in kilometers).

Combines a latitude-dependent cosine area weighting and a longitudinally varying weighting window from scipy.signal.

The constant width in terms of actual distance should ensure a consistent cut-off wavelength across the sphere. The window is returned folded (see rwguide.zonalization.weighting.fold_periodic()).

Note

We found in testing that a constant-distance width window has rather poor PV conservation properties when used with rolling zonalization. Consider using a constant-longitude window instead (fixed_deg_window()).

Parameters:
  • da_lon (xarray.DataArray) – Longitude coordinate.

  • da_lat (xarray.DataArray) – Latitude coordinate.

  • name (str, optional) – Name assigned to the output.

  • kwargs (any) – Keyword arguments given to rwguide.zonalization.fixed_km_window().

Returns:

Weighting window.

Return type:

xarray.DataArray

rwguide.xarray.zonalization.rolling_mean_background(da_area, da_sg, *, vectorize=True, names=None)

Weighted longitudinal rolling mean.

A simple procedure to approximate the background-state isentropic density field for rolling zonalization.

Parameters:
  • da_area (xarray.DataArray) – Weighting window.

  • da_sg (xarray.DataArray) – Isentropic density field. Core dimensions: latitude, longitude.

  • vectorize (boolean, optional) – Use vectorization of xarray.apply_ufunc().

  • names (dict, optional) – Variable name override.

Returns:

Rolling-mean background isentropic density.

Return type:

xarray.DataArray

rwguide.xarray.zonalization.zonalize(da_area, da_av, da_sg, da_bg=None, *, vectorize=True, names=None)

Zonalization (hemispheric or sectoral/fixed-window).

Integrations are carried out with a conditional boxcounting quadrature scheme. Northern and southern hemispheres are automatically detected and zonalized separately. Regions with 0-valued isentropic density are omitted in the surface integrals.

The PV contours for the zonalizations are automatically determined based on the input data: 10 contours, linearly spaced between 0 and the maximum (NH) or minimum (SH) value of PV found in the domain, are first zonalized to provide an approximation of the zonalized PV profile. The first guess PV profile is then interpolated to the input latitude grid and the interpolated PV values at each grid latitude are then used as contours in a second zonalization pass to obtain a refined zonalized profile. The refined zonalized PV profile is finally interpolated to the input latitude and returned.

E.g., hemispheric zonalization based on a zonal-mean background state isentropic density profile (Ghinassi et al. 2018):

>>> area = area_weights(isen.coords["latitude"])
>>> zonalize(area, isen["av"], isen["sg"], isen["sg"].mean(dim="longitude"))
Parameters:
  • da_area (xarray.DataArray) – Weighting window or meridional profile of weighting coefficients (for “normal” hemispheric zonalization).

  • da_av (xarray.DataArray) – Isentropic absolute vorticity field in 1 / s. Core dimensions: latitude, longitude.

  • da_sg (xarray.DataArray) – Isentropic density field in kg / K / m². Core dimensions: latitude, longitude.

  • da_bg (xarray.DataArray, optional) – Background state isentropic density field in kg / K / m². Core dimensions: latitude, longitude. If none given, da_sg is used for the background too. If no longitude dimension is found, the meridional profile is used at every longitude.

  • vectorize (boolean, optional) – Use vectorization of xarray.apply_ufunc().

  • names (dict, optional) – Variable name override.

Returns:

Zonalized PV in PVU.

Return type:

xarray.DataArray

rwguide.xarray.zonalization.zonalize_rolling(da_area, da_av, da_sg, da_bg=None, *, vectorize=True, names=None)

Rolling zonalization (Polster and Wirth 2023).

The PV contours for the zonalization are automatically determined based on the input data (see zonalize()). Integrations are carried out with a conditional boxcounting quadrature scheme. Northern and southern hemispheres are automatically detected and zonalized separately. Regions with 0-valued isentropic density are omitted in the surface integrals.

Parameters:
  • da_area (xarray.DataArray) – Weighting window.

  • da_av (xarray.DataArray) – Isentropic absolute vorticity field in 1 / s. Core dimensions: latitude, longitude.

  • da_sg (xarray.DataArray) – Isentropic density field in kg / K / m². Core dimensions: latitude, longitude.

  • da_bg (xarray.DataArray, optional) – Background state isentropic density field in kg / K / m². Core dimensions: latitude, longitude. If none given, da_sg is used for the background too.

  • vectorize (boolean, optional) – Use vectorization of xarray.apply_ufunc().

  • names (dict, optional) – Variable name override.

Returns:

Rolling zonalized PV in PVU.

Return type:

xarray.DataArray

Numpy Interface

Computation of zonalized (equivalent-latitude) background states.

rwguide.zonalization.fixed_deg_window(lon, lat, width=60.0, window='boxcar')

Weighting window with a constant-longitude width.

Parameters:
  • lon (numpy.ndarray) – Longitude coordinates.

  • lat (numpy.ndarray) – Latitude coordinates.

  • width (number) – Window width in °.

  • window (str | tuple) – Window function. Given as the name parameter to weighting.get_window().

Return type:

numpy.ndarray

rwguide.zonalization.fixed_km_window(lon, lat, width=8000.0, window='boxcar')

Weighting window with a constant-distance width.

Parameters:
  • lon (numpy.ndarray) – Longitude coordinates.

  • lat (numpy.ndarray) – Latitude coordinates.

  • width (number) – Window width in km.

  • window (str | tuple) – Window function. Given as the name parameter to weighting.get_window().

Return type:

numpy.ndarray

rwguide.zonalization.rolling_mean_background(window, field)
rwguide.zonalization.zonalize(area, av, sg, bg, nh_last=None, sh_first=None)
rwguide.zonalization.zonalize_rolling(area, av, sg, bg, nh_last=None, sh_first=None)

Weighting Windows

The implementation uses a generalized approach to the area weighting of the surface integrals on the sphere evaluated during the zonalization. This allows for customization of the PV re-arrangement by using arbitrary weighting functions in addition to the normal latitude-dependent cosine-based area weighting.

Regular mass integrals over isentropic density σ on an isentropic surface of the sphere are given by

∫∫ σ(ϕ,λ) a cos(ϕ) ,

where a cos(ϕ) constitutes an area-weighting that accounts for the difference in size of grid boxes of a latitude-longitude grid. We generalize the weighting to

∫∫ σ(ϕ,λ) w(ϕ,λ) ,

where w is now an arbitrary weighting function that can vary both in latitude and longitude. In the limit of a simple boxcar window, we arrive at basic rolling zonalization where PV is zonalized normally in a rolling window. But the weighting window can, e.g., be designed with tapered edges and can widen on the lat-lon grid towards the poles to obtain a fixed window width in terms of actual distance.

We provide a set of functions to build customized weighting windows on the sphere based on cosine area weighting and window functions from scipy.signal:

rwguide.zonalization.weighting.area_weights(lat, radius=1.0, nlon=1)

Area weights for a lat-lon-grid on the sphere.

rwguide.zonalization.weighting.clip_width(width, nmin=3, nmax=3600)

Convert to integer (# of gridpoints) and clip to the given range.

rwguide.zonalization.weighting.earth_circumference(lat)

Circumference of the Earth (in m) at the given latitude(s).

rwguide.zonalization.weighting.fold_periodic(window, n, agg='max')

Fold the window into the periodic domain.

If the window fits into the domain, it is returned unchanged. Otherwise the window is centered in the domain (at n//2). Example:

n = 14, window.size = 18

|<-window.size-->|        |<-window.size-->|
|<----n----->|   |        | |<----n----->| |
|     ..··.. |   |     →  | |   ..··..   | |
...···      ···...        ...···   |  ···...
|------------|----        --|------|-----|--
012345678901234567          01234567890123
                                   ^ n//2
Parameters:
  • window (np.ndarray) – The window to fold.

  • n (number) – The width of the periodic domain (# of gridpoints).

  • agg ("sum" | "clip" | "max", optional) – How to aggregate values in the folded regions. clip restricts values to the maximum of the original window after summation.

Return type:

np.ndarray

rwguide.zonalization.weighting.get_window(name, width, fold=None)

Create a weighting window for the sphere.

Parameters:
  • name (string | tuple) – Which window function to use. Choose a name from scipy.signal.windows or specify a tuple containing the name and additional arguments for the window creation function. Examples: "boxcar", ("tukey", 0.5).

  • width (number | iterable) – Window width in # of gridpoints. Constant or latitude-dependent.

  • fold (None | Tuple | number, optional) – Settings for folding the window into the periodic domain. Provide the n and agg arguments for fold_periodic().

Return type:

numpy.ndarray

rwguide.zonalization.weighting.normalize(window)

Normalize such that window weights add to 1.

rwguide.zonalization.weighting.oddify(width)

Return the next larger or equal odd number.

rwguide.zonalization.weighting.product(window_lon, window_lat)

2D product of a (1D/2D) longitude and a (1D/2D) latitude window.

rwguide.zonalization.weighting.wavelength_to_wavenumber(wavelength, lat)

Convert wavelength to wavenumber (latitude-dependent).

rwguide.zonalization.weighting.wavelength_to_width(wavelength, lat, n)

Convert wavelength to width (latitude-dependent).

rwguide.zonalization.weighting.wavenumber_to_width(wavenumber, n)

Convert wavenumber to width (length or # of gridpoints).