Isentropic PV Gradient¶
Interpolate to isentropic levels and compute potential vorticity (PV) and its gradient as a waveguide diagnostic.
Example¶
>>> import xarray as xr
>>> from rwguide.xarray import pvgradient
Loading input data (temperature, wind components on pressure levels):
>>> uvt = xr.open_dataset("data/ERA5/ERA5-2018-tuv-1.5.nc", chunks={ "time": 20 })
>>> uvt
<xarray.Dataset>
Dimensions: (longitude: 240, latitude: 121, level: 18, time: 1460)
Coordinates:
* longitude (longitude) float32 -180.0 -178.5 -177.0 ... 175.5 177.0 178.5
* latitude (latitude) float32 90.0 88.5 87.0 85.5 ... -87.0 -88.5 -90.0
* level (level) int32 50 70 100 150 200 250 ... 600 650 700 750 800 850
* time (time) datetime64[ns] 2018-01-01 ... 2018-12-31T18:00:00
Data variables:
u (time, level, latitude, longitude) float32 dask.array<...>
v (time, level, latitude, longitude) float32 dask.array<...>
t (time, level, latitude, longitude) float32 dask.array<...>
Interpolating to isentropes and computing PV and associated variables:
>>> isen = pvgradient.isob_to_isen_all(uvt.assign_coords({ "isentrope": [340., 330., 320.]}))
>>> 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:
pres (time, isentrope, latitude, longitude) float64 dask.array<...>
u (time, isentrope, latitude, longitude) float64 dask.array<...>
v (time, isentrope, latitude, longitude) float64 dask.array<...>
sg (time, isentrope, latitude, longitude) float64 dask.array<...>
av (time, isentrope, latitude, longitude) float64 dask.array<...>
pv (time, isentrope, latitude, longitude) float64 dask.array<...>
Computing the magnitude of the gradient of the logarithm of PV as a waveguidability proxy:
>>> pvgradient.norm_grad_log_abs(isen["pv"])
<xarray.DataArray 'ngl_pv' (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.pvgradient.absolute_vorticity(da_u, da_v, *, names=None)¶
Absolute vorticity of a horizontal wind field.
Coriolis parameter + vertical component of relative vorticity. The planet is Earth with angular frequency 7.292e-5 1/s and radius 6371.2 km.
- Parameters:
da_u (xarray.DataArray) – Zonal wind component field in m / s. Core dimensions: latitude, longitude.
da_v (xarray.DataArray) – Meridional wind component field in m / s. Core dimensions: latitude, longitude.
names (dict, optional) – Variable name override.
- Returns:
Absolute vorticity field in 1 / s.
- Return type:
xarray.DataArray
- rwguide.xarray.pvgradient.curl(da_u, da_v, *, vectorize=True, names=None)¶
Vertical component of rotation ∇ ⨯ · on the sphere.
The sphere has the dimensions of Earth (6371.2 km radius).
- Parameters:
da_u (xarray.DataArray) – Input zonal vector component. Core dimensions: latitude, longitude.
da_v (xarray.DataArray) – Input meridional vector component. Core dimensions: latitude, longitude.
vectorize (boolean, optional) – Use vectorization of
xarray.apply_ufunc()
.names (dict, optional) – Variable name override.
- Returns:
Curl of the vector input field.
- Return type:
xarray.DataArray
- rwguide.xarray.pvgradient.horizontal_gradient(da, *, vectorize=True, names=None)¶
Horizontal gradient ∇· on the sphere.
The sphere has the dimensions of Earth (6371.2 km radius).
- Parameters:
da (xarray.DataArray) – Input field. Core dimensions: latitude, longitude.
vectorize (boolean, optional) – Use vectorization of
xarray.apply_ufunc()
.names (dict, optional) – Variable name override.
- Returns:
Zonal (
gradx_*
) and meridional (grady_*
) gradient of the scalar input field.- Return type:
xarray.Dataset
- rwguide.xarray.pvgradient.interpolate_isob_to_isen(ds, da_isen=None, *, names=None)¶
Interpolation from pressure levels to isentropic levels.
Either assign the target isentropic levels as a coordinate (isentrope) to the input dataset or provide da_isen. The pressure coordinate (level) must be specified in ascending order, the isentropic levels in descending order (both top-down). Out-of-range values are filled with NaN.
- Parameters:
ds (xarray.Dataset) – Input data fields. Core dimensions: level, latitude, longitude.
da_isen (xarray.DataArray, optional) – Override for isentropic levels to interpolate to (in K).
names (dict, optional) – Variable name override.
- Returns:
Dataset with all interpolated fields.
- Return type:
xarray.Dataset
- rwguide.xarray.pvgradient.isentropic_density_isob(da_pt, *, names=None)¶
Isentropic density from potential temperature on pressure levels.
Isentropic density σ = - g⁻¹ ∂p/∂θ computed as -(g ∂θ/∂p)⁻¹ using a 2nd-order finite difference stencil. 0-values are assigned to grid points with unstable stratification.
The pressure coordinate (level) must be specified in hPa.
- Parameters:
da_pt (xarray.DataArray) – Potential temperature field in K. Core dimensions: level, latitude, longitude.
names (dict, optional) – Variable name override.
- Returns:
Isentropic density in kg/K/m².
- Return type:
xarray.DataArray
- rwguide.xarray.pvgradient.isob_to_isen_all(ds, da_isen=None, *, vectorize=True, names=None)¶
All-in-one isentropic PV computation from pressure-level inputs.
Computes pressure (in hPa), isentropic density (in kg / K / m²), absolute vorticity (in 1 / s) and Ertel potential vorticity (in PVU) on isentropes (in K, descending order) from temperature (in K) and zonal and meridional wind components (in m / s) on pressure levels (in hPa, ascending order).
Note
Using this function is generally faster and easier than going through
potential_temperature_isob()
,isentropic_density_isob()
,interpolate_isob_to_isen()
,absolute_vorticity()
andpotential_vorticity_isen()
individually.- Parameters:
ds (xarray.Dataset) – Input dataset with temperature and horizontal wind components as well as isentropic levels to interpolate to. Core dimensions: level, latitude, longitude.
da_isen (xarray.DataArray, optional) – Override for isentropic levels to interpolate to (in K).
vectorize (boolean, optional) – Use vectorization of
xarray.apply_ufunc()
.names (dict, optional) – Variable name override.
- Returns:
Dataset containing all output variables. Core dimensions: isentrope, latitude, longitude.
- Return type:
xarray.Dataset
- rwguide.xarray.pvgradient.norm_grad_log_abs(da, *, threshold=0.0, vectorize=True, names=None)¶
Magnitude of the gradient of the logarithm of the input field.
Used as a waveguide diagnostic of isentropic Ertel PV. To first, order ∇(log(|IPV|)) ≈ 1/f₀ ∇QGPV ~ ∇²u, where IPV is isentropic potential vorticity, QGPV is quasi-geostropic potential vorticity and u is the zonal wind (Martius et al. 2010).
Note
d/dx log(x/a) = d/dx log(x) = 1/x for all a = const, i.e. the result of this function does not depend on the unit of its argument.
- Parameters:
da (xarray.DataArray) – Input field. Core dimensions: latitude, longitude.
threshold (number, optional) – Threshold underneath which values of |·| are ignored, to deal with regions where the input field approaches zero and the logarithm tends toward negative infinity.
vectorize (boolean, optional) – Use vectorization of
xarray.apply_ufunc()
.names (dict, optional) – Variable name override.
- Returns:
‖∇(log(|·|))‖₂ of the input field in 1 / m.
- Return type:
xarray.DataArray
- rwguide.xarray.pvgradient.potential_temperature_isob(da_t, *, names=None)¶
Potential temperature from temperature on pressure levels.
The pressure coordinate (level) must be specified in hPa.
- Parameters:
da_t (xarray.DataArray) – Temperature field in K. Core dimensions: level, latitude, longitude.
names (dict, optional) – Variable name override.
- Returns:
Potential temperature in K.
- Return type:
xarray.DataArray
- rwguide.xarray.pvgradient.potential_vorticity_isen(da_av, da_sg, *, names=None)¶
Isentropic Ertel PV from potential temperature and isentropic density.
Potential vorticity q = ζₐ / σ. NaN-values are assigned where σ = 0.
Note
The output is in potential vorticity units. 1 PVU = K m² / s / kg.
- Parameters:
da_av (xarray.DataArray) – 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.
names (dict, optional) – Variable name override.
- Returns:
Potential vorticity in PVU.
- Return type:
xarray.DataArray
Numpy Interface¶
- rwguide.pvgradient.absolute_vorticity(lat, u, v)¶
- rwguide.pvgradient.curl(lat, u, v)¶
- rwguide.pvgradient.horizontal_gradient(lat, field)¶
- rwguide.pvgradient.interpolate_isob_to_isen(isob, isen, th, *fields)¶
- rwguide.pvgradient.isentropic_density_isob(isob, th)¶
- rwguide.pvgradient.isob_to_isen_all(isob, isen, lat, t, u, v, psfc=None)¶
- rwguide.pvgradient.norm_grad_log_abs(lat, pv, threshold=0.0)¶
- rwguide.pvgradient.potential_temperature_isob(isob, t)¶
- rwguide.pvgradient.potential_vorticity_isen(av, sg)¶