potential_vorticity_baroclinic

metpy.calc.potential_vorticity_baroclinic(potential_temperature, pressure, u, v, dx=None, dy=None, latitude=None, x_dim=- 1, y_dim=- 2, vertical_dim=- 3)[source]

Calculate the baroclinic potential vorticity.

\[PV = -g \left(\frac{\partial u}{\partial p}\frac{\partial \theta}{\partial y} - \frac{\partial v}{\partial p}\frac{\partial \theta}{\partial x} + \frac{\partial \theta}{\partial p}(\zeta + f) \right)\]

This formula is based on equation 4.5.93 [Bluestein1993]

Parameters
  • potential_temperature ((…, P, M, N) xarray.DataArray or pint.Quantity) – potential temperature

  • pressure ((…, P, M, N) xarray.DataArray or pint.Quantity) – vertical pressures

  • u ((…, P, M, N) xarray.DataArray or pint.Quantity) – x component of the wind

  • v ((…, P, M, N) xarray.DataArray or pint.Quantity) – y component of the wind

  • dx (pint.Quantity, optional) – The grid spacing(s) in the x-direction. If an array, there should be one item less than the size of u along the applicable axis. Optional if xarray.DataArray with latitude/longitude coordinates used as input.

  • dy (pint.Quantity, optional) – The grid spacing(s) in the y-direction. If an array, there should be one item less than the size of u along the applicable axis. Optional if xarray.DataArray with latitude/longitude coordinates used as input.

  • latitude (pint.Quantity, optional) – Latitude of the wind data. Optional if xarray.DataArray with latitude/longitude coordinates used as input. Note that an argument without units is treated as dimensionless, which translates to radians.

  • x_dim (int, optional) – Axis number of x dimension. Defaults to -1 (implying […, Z, Y, X] order). Automatically parsed from input if using xarray.DataArray.

  • y_dim (int, optional) – Axis number of y dimension. Defaults to -2 (implying […, Z, Y, X] order). Automatically parsed from input if using xarray.DataArray.

  • vertical_dim (int, optional) – Axis number of vertical dimension. Defaults to -3 (implying […, Z, Y, X] order). Automatically parsed from input if using xarray.DataArray.

Returns

(…, P, M, N) xarray.DataArray or pint.Quantity – baroclinic potential vorticity

Notes

The same function can be used for isobaric and isentropic PV analysis. Provide winds for vorticity calculations on the desired isobaric or isentropic surface. At least three layers of pressure/potential temperature are required in order to calculate the vertical derivative (one above and below the desired surface). The first two terms will be zero if isentropic level data is used due to the gradient of theta in both the x and y-directions will be zero since you are on an isentropic surface.

This function expects pressure/isentropic level to increase with increasing array element (e.g., from higher in the atmosphere to closer to the surface. If the pressure array is one-dimensional, and not given as xarray.DataArray, p[:, None, None] can be used to make it appear multi-dimensional.)

Changed in version 1.0: Changed signature from (potential_temperature, pressure, u, v, dx, dy, lats)