Skip to content

Commit

Permalink
DOC: Improve q_vector documentation (Fixes #3689)
Browse files Browse the repository at this point in the history
Clarify the optional role of static_stability and how that maps to
various GEMPAK calculations. Also ensure that values for
`static_stability` are properly checked for dimensionality.
  • Loading branch information
dopplershift committed Jan 6, 2025
1 parent ea2fcfc commit 407eb0d
Showing 1 changed file with 8 additions and 1 deletion.
9 changes: 8 additions & 1 deletion src/metpy/calc/kinematics.py
Original file line number Diff line number Diff line change
Expand Up @@ -1275,7 +1275,8 @@ def inertial_advective_wind(
broadcast=('u', 'v', 'temperature', 'pressure', 'static_stability', 'parallel_scale',
'meridional_scale')
)
@check_units('[speed]', '[speed]', '[temperature]', '[pressure]', '[length]', '[length]')
@check_units('[speed]', '[speed]', '[temperature]', '[pressure]', '[length]', '[length]',
'[energy] / [mass] / [pressure]**2')
def q_vector(
u,
v,
Expand Down Expand Up @@ -1307,6 +1308,12 @@ def q_vector(
- 2 \nabla_p \cdot \vec{Q} -
\frac{R}{\sigma p} \beta \frac{\partial T}{\partial x}
By default, this function uses a unitless value of 1 for ``static_stability``, which
replicates the functionality of the GEMPAK ``QVEC`` function. If a value is given for
``static_stability``, it should have dimensionality of ``energy / mass / pressure^2``, and
will result in behavior that matches that of GEMPAK's ``QVCL`` function;
`static_stability` can be used to calculate this value if desired.
Parameters
----------
u : (..., M, N) `xarray.DataArray` or `pint.Quantity`
Expand Down

0 comments on commit 407eb0d

Please sign in to comment.