inertial_advective_wind#

metpy.calc.inertial_advective_wind(u, v, u_geostrophic, v_geostrophic, dx=None, dy=None, latitude=None, x_dim=-1, y_dim=-2, *, parallel_scale=None, meridional_scale=None, longitude=None, crs=None)[source]#

Calculate the inertial advective wind.

\[\frac{\hat k}{f} \times (\vec V \cdot \nabla)\hat V_g\]
\[\frac{\hat k}{f} \times \left[ \left( u \frac{\partial u_g}{\partial x} + v \frac{\partial u_g}{\partial y} \right) \hat i + \left( u \frac{\partial v_g} {\partial x} + v \frac{\partial v_g}{\partial y} \right) \hat j \right]\]
\[\left[ -\frac{1}{f}\left(u \frac{\partial v_g}{\partial x} + v \frac{\partial v_g}{\partial y} \right) \right] \hat i + \left[ \frac{1}{f} \left( u \frac{\partial u_g}{\partial x} + v \frac{\partial u_g}{\partial y} \right) \right] \hat j\]

This formula is based on equation 27 of [Rochette2006].

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

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

  • u_geostrophic ((…, M, N) xarray.DataArray or pint.Quantity) – x component of the geostrophic (advected) wind

  • v_geostrophic ((…, M, N) xarray.DataArray or pint.Quantity) – y component of the geostrophic (advected) 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 […, Y, X] order). Automatically parsed from input if using xarray.DataArray.

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

  • parallel_scale (pint.Quantity, optional) – Parallel scale of map projection at data coordinate. Optional if xarray.DataArray with latitude/longitude coordinates and MetPy CRS used as input. Also optional if longitude, latitude, and crs are given. If otherwise omitted, calculation will be carried out on a Cartesian, rather than geospatial, grid. Keyword-only argument.

  • meridional_scale (pint.Quantity, optional) – Meridional scale of map projection at data coordinate. Optional if xarray.DataArray with latitude/longitude coordinates and MetPy CRS used as input. Also optional if longitude, latitude, and crs are given. If otherwise omitted, calculation will be carried out on a Cartesian, rather than geospatial, grid. Keyword-only argument.

  • longitude (pint.Quantity, optional) – Longitude of data. Optional if xarray.DataArray with latitude/longitude coordinates used as input. Also optional if parallel_scale and meridional_scale are given. If otherwise omitted, calculation will be carried out on a Cartesian, rather than geospatial, grid. Keyword-only argument.

  • crs (pyproj.crs.CRS, optional) – Coordinate Reference System of data. Optional if xarray.DataArray with MetPy CRS used as input. Also optional if parallel_scale and meridional_scale are given. If otherwise omitted, calculation will be carried out on a Cartesian, rather than geospatial, grid. Keyword-only argument.

Returns:

Notes

Many forms of the inertial advective wind assume the advecting and advected wind to both be the geostrophic wind. To do so, pass the x and y components of the geostrophic wind for u and u_geostrophic/v and v_geostrophic.

Changed in version 1.0: Changed signature from (u, v, u_geostrophic, v_geostrophic, dx, dy, lats)