Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Check that pressure decreases monotonically in parcel calculations #1967

Merged
merged 4 commits into from
Aug 2, 2021
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions src/metpy/calc/__init__.py
Original file line number Diff line number Diff line change
Expand Up @@ -5,6 +5,7 @@

from .basic import * # noqa: F403
from .cross_sections import * # noqa: F403
from .exceptions import * # noqa: F403
from .indices import * # noqa: F403
from .kinematics import * # noqa: F403
from .thermo import * # noqa: F403
Expand Down
8 changes: 8 additions & 0 deletions src/metpy/calc/exceptions.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,8 @@
# Copyright (c) 2021 MetPy Developers.
# Distributed under the terms of the BSD 3-Clause License.
# SPDX-License-Identifier: BSD-3-Clause
"""Contains specific exceptions raised by calculations."""


class InvalidSoundingError(ValueError):
"""Raise when a sounding does not meet thermodynamic calculation expectations."""
16 changes: 16 additions & 0 deletions src/metpy/calc/thermo.py
Original file line number Diff line number Diff line change
Expand Up @@ -10,6 +10,7 @@
import scipy.optimize as so
import xarray as xr

from .exceptions import InvalidSoundingError
from .tools import (_greater_or_close, _less_or_close, _remove_nans, find_bounding_indices,
find_intersections, first_derivative, get_layer)
from .. import constants as mpconsts
Expand Down Expand Up @@ -873,13 +874,28 @@ def parcel_profile_with_lcl_as_dataset(pressure, temperature, dewpoint):
)


def _check_pressure(pressure):
"""Check that pressure decreases monotonically.

Returns True if the pressure decreases monotonically; otherwise, returns False.
"""
return np.all(pressure[:-1] > pressure[1:])


def _parcel_profile_helper(pressure, temperature, dewpoint):
"""Help calculate parcel profiles.

Returns the temperature and pressure, above, below, and including the LCL. The
other calculation functions decide what to do with the pieces.

"""
# Check that pressure is monotonic decreasing.
if not _check_pressure(pressure):
msg = """
Pressure does not decrease monotonically in your sounding.
Using scipy.signal.medfilt may fix this."""
raise InvalidSoundingError(msg)

# Find the LCL
press_lcl, temp_lcl = lcl(pressure[0], temperature, dewpoint)
press_lcl = press_lcl.to(pressure.units)
Expand Down
18 changes: 14 additions & 4 deletions tests/calc/test_thermo.py
Original file line number Diff line number Diff line change
Expand Up @@ -11,10 +11,10 @@
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, gradient_richardson_number, isentropic_interpolation,
isentropic_interpolation_as_dataset, lcl, lfc, lifted_index,
mixed_layer, mixed_layer_cape_cin, mixed_parcel, mixing_ratio,
mixing_ratio_from_relative_humidity,
exner_function, gradient_richardson_number, InvalidSoundingError,
isentropic_interpolation, isentropic_interpolation_as_dataset, 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, parcel_profile,
parcel_profile_with_lcl, parcel_profile_with_lcl_as_dataset,
Expand Down Expand Up @@ -190,6 +190,16 @@ def test_parcel_profile_lcl():
assert_array_almost_equal(prof, true_prof, 2)


def test_parcel_profile_lcl_not_monotonic():
"""Test parcel profile with lcl calculation."""
with pytest.raises(InvalidSoundingError):
p = np.array([1004., 1000., 943., 925., 928., 850., 839., 749., 700.]) * units.hPa
t = np.array([24.2, 24., 20.2, 21.6, 21.4, 20.4, 20.2, 14.4, 13.2]) * units.degC
td = np.array([21.9, 22.1, 19.2, 20.5, 20.4, 18.4, 17.4, 8.4, -2.8]) * units.degC

_ = parcel_profile_with_lcl(p, t, td)


def test_parcel_profile_with_lcl_as_dataset():
"""Test parcel profile with lcl calculation with xarray."""
p = np.array([1004., 1000., 943., 928., 925., 850., 839., 749., 700., 699.]) * units.hPa
Expand Down
13 changes: 8 additions & 5 deletions tutorials/xarray_tutorial.py
Original file line number Diff line number Diff line change
Expand Up @@ -260,21 +260,24 @@
#########################################################################
# For profile-based calculations (and most remaining calculations in the ``metpy.calc``
# module), xarray ``DataArray``\s are accepted as inputs, but the outputs remain Pint
# Quantities (typically scalars)
# Quantities (typically scalars). Note that MetPy's profile calculations (such as CAPE and
# CIN) require the sounding to be ordered from highest to lowest pressure. As seen earlier
# in this tutorial, this data is ordered the other way, so we need to reverse the inputs
# to ``mpcalc.surface_based_cape_cin``.

data_at_point = data.metpy.sel(
time1='2017-09-05 12:00',
latitude=40 * units.degrees_north,
latitude=30 * units.degrees_north,
dopplershift marked this conversation as resolved.
Show resolved Hide resolved
longitude=260 * units.degrees_east
)
dewpoint = mpcalc.dewpoint_from_relative_humidity(
data_at_point['Temperature_isobaric'],
data_at_point['Relative_humidity_isobaric']
)
cape, cin = mpcalc.surface_based_cape_cin(
data_at_point['isobaric3'],
data_at_point['Temperature_isobaric'],
dewpoint
data_at_point['isobaric3'][::-1],
data_at_point['Temperature_isobaric'][::-1],
dewpoint[::-1]
)
cape

Expand Down