From 7ee6cef5b0a38ae44e982138587faf15ca6adc77 Mon Sep 17 00:00:00 2001 From: Joern Ungermann Date: Wed, 1 Apr 2020 09:45:38 +0200 Subject: [PATCH] Added function to calculate richardson number. Apply suggestions from code review Fixed flake8 errors of pipeline and accomodated further suggestions. Rebased to fix merge conflicts Co-Authored-By: Jon Thielen See #647 and #628 --- docs/_templates/overrides/metpy.calc.rst | 1 + docs/references.rst | 2 ++ src/metpy/calc/thermo.py | 37 ++++++++++++++++++++++++ tests/calc/test_thermo.py | 18 ++++++++++-- 4 files changed, 56 insertions(+), 2 deletions(-) diff --git a/docs/_templates/overrides/metpy.calc.rst b/docs/_templates/overrides/metpy.calc.rst index abd751ead2e..4a05e0d9cf3 100644 --- a/docs/_templates/overrides/metpy.calc.rst +++ b/docs/_templates/overrides/metpy.calc.rst @@ -127,6 +127,7 @@ Boundary Layer/Turbulence brunt_vaisala_frequency_squared brunt_vaisala_period friction_velocity + gradient_richardson_number tke diff --git a/docs/references.rst b/docs/references.rst index 578e473aaa1..ca184da3fa8 100644 --- a/docs/references.rst +++ b/docs/references.rst @@ -92,6 +92,8 @@ References .. [Garratt1994] Garratt, J.R., 1994: *The Atmospheric Boundary Layer*. Cambridge University Press, 316 pp. +.. [Holton2004] Holton, J. R., 2004: *An Introduction to Dynamic Meteorology*. 4th ed. Academic Press, 535 pp. + .. [Koch1983] Koch, S. E., M. DesJardins, and P. J. Kocin, 1983: An interactive Barnes objective map analysis scheme for use with satellite and conventional data. **J. Appl. Meteor. Climatol.**, *22*, 1487-1503, diff --git a/src/metpy/calc/thermo.py b/src/metpy/calc/thermo.py index e7804b1b828..f7007ee1209 100644 --- a/src/metpy/calc/thermo.py +++ b/src/metpy/calc/thermo.py @@ -2617,3 +2617,40 @@ def lifted_index(pressure, temperature, parcel_profile): # calculate the lifted index. lifted_index = T500 - Tp500.to(units.degC) return lifted_index + + +@exporter.export +@preprocess_xarray +@check_units('[length]', '[temperature]', '[speed]', '[speed]') +def gradient_richardson_number(height, potential_temperature, u, v, axis=0): + r"""Calculate the gradient (or flux) Richardson number. + + .. math:: Ri = (g/\theta) * \frac{\left(\partial \theta/\partial z\)} + {[\left(\partial u / \partial z\right)^2 + \left(\partial v / \partial z\right)^2} + + See [Holton2004]_ pg. 121-122. As noted by [Holton2004]_, flux Richardson + number values below 0.25 indicate turbulence. + + Parameters + ---------- + height : `pint.Quantity` + Atmospheric height + potential_temperature : `pint.Quantity` + Atmospheric potential temperature + u : `pint.Quantity` + x component of the wind + v : `pint.Quantity` + y component of the wind + axis : int, optional + The axis corresponding to vertical, defaults to 0. + + Returns + ------- + `pint.Quantity` + Gradient Richardson number + """ + dthetadz = first_derivative(potential_temperature, x=height, axis=axis) + dudz = first_derivative(u, x=height, axis=axis) + dvdz = first_derivative(v, x=height, axis=axis) + + return (mpconsts.g / potential_temperature) * (dthetadz / (dudz ** 2 + dvdz ** 2)) diff --git a/tests/calc/test_thermo.py b/tests/calc/test_thermo.py index a6e07d8c9a8..69ec2511256 100644 --- a/tests/calc/test_thermo.py +++ b/tests/calc/test_thermo.py @@ -11,8 +11,9 @@ brunt_vaisala_period, cape_cin, density, dewpoint, dewpoint_from_relative_humidity, dewpoint_from_specific_humidity, dry_lapse, dry_static_energy, el, equivalent_potential_temperature, - exner_function, isentropic_interpolation, lcl, lfc, lifted_index, - mixed_layer, mixed_layer_cape_cin, mixed_parcel, mixing_ratio, + exner_function, gradient_richardson_number, + isentropic_interpolation, lcl, lfc, lifted_index, mixed_layer, + mixed_layer_cape_cin, mixed_parcel, mixing_ratio, mixing_ratio_from_relative_humidity, mixing_ratio_from_specific_humidity, moist_lapse, moist_static_energy, most_unstable_cape_cin, most_unstable_parcel, @@ -1582,3 +1583,16 @@ def test_lifted_index(): parcel_prof = parcel_profile(pressure, temperature[0], dewpoint[0]) LI = lifted_index(pressure, temperature, parcel_prof) assert_almost_equal(LI, -7.9176350 * units.delta_degree_Celsius, 2) + + +def test_gradient_richardson_number(): + """Test gradient Richardson number calculation.""" + theta = units('K') * np.asarray([254.5, 258.3, 262.2]) + u_wnd = units('m/s') * np.asarray([-2., -1.1, 0.23]) + v_wnd = units('m/s') * np.asarray([3.3, 4.2, 5.2]) + height = units('km') * np.asarray([0.2, 0.4, 0.6]) + + result = gradient_richardson_number(height, theta, u_wnd, v_wnd) + expected = np.asarray([24.2503551, 13.6242603, 8.4673744]) + + assert_array_almost_equal(result, expected, 4)