Skip to content
Open
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
307 changes: 307 additions & 0 deletions src/metpy/calc/boundarylayer.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,307 @@
# Copyright (c) 2024 MetPy Developers.
# Distributed under the terms of the BSD 3-Clause License.
# SPDX-License-Identifier: BSD-3-Clause
"""
Contains a collection of boundary layer height estimations.

References
----------
[Col14]: Collaud Coen, M., Praz, C., Haefele, A., Ruffieux, D., Kaufmann, P., and Calpini, B. (2014)
Determination and climatology of the planetary boundary layer height above the Swiss plateau by in situ and remote sensing measurements as well as by the COSMO-2 model
Atmos. Chem. Phys., 14, 13205–13221.

[HL06]: Hennemuth, B., & Lammert, A. (2006):
Determination of the atmospheric boundary layer height from radiosonde and lidar backscatter.
Boundary-Layer Meteorology, 120(1), 181-200.

[Guo16]: Guo, J., Miao, Y., Zhang, Y., Liu, H., Li, Z., Zhang, W., ... & Zhai, P. (2016)
The climatology of planetary boundary layer height in China derived from radiosonde and reanalysis data.
Atmos. Chem. Phys, 16(20), 13309-13319.

[Sei00]: Seidel, D. J., Ao, C. O., & Li, K. (2010)
Estimating climatological planetary boundary layer heights from radiosonde observations: Comparison of methods and uncertainty analysis.
Journal of Geophysical Research: Atmospheres, 115(D16).

[VH96]: Vogelezang, D. H. P., & Holtslag, A. A. M. (1996)
Evaluation and model impacts of alternative boundary-layer height formulations.
Boundary-Layer Meteorology, 81(3-4), 245-269.
"""
import numpy as np
from copy import deepcopy

import metpy.calc as mpcalc
import metpy.constants as mpconsts
from metpy.units import units


def smooth(val, span):
Copy link
Contributor

@DWesl DWesl Jul 19, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

XArray calls this a rolling mean. So does pandas.
Bottleneck calls this a moving-window mean.
SciPy appears to call the same thing a uniform filter.

These would likely work better in the See Also section than as a change in the name.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Thanks for the references! I knew some equivalent functions were already existing but they are not quite exactly the same (Xarray works on xarray.Dataset, Scipy has a slightly different strategy at the edges) and, given that the function is simple enough, it was less work to write it than to look for the existing one.

Bottleneck's function seems to do exactly what I want but it is not listed in the Metpy's dependencies. Do you think it's worth adding it so I can use their moving-mean function?

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

You'd have to ask one of the maintainers, but, in the meantime, would this be faster?

cumulative_sums = np.nancumsum(val)
rolling_sums = cumulative_sums[span:] - cumulative_sums[:-span]
valid_index = np.isfinite(val)
cumulative_count = np.cumsum(valid_index)
rolling_count = cumulative_count[span:] - cumulative_count[:-span]
rolling_means = rolling_sums / rolling_count

You'd need to pre-allocate rolling_means and handle the edges still but it should work. (Alternately, use np.lib.stride_tricks.sliding_window_view with np.nanmean, but the note about it being slow is warranted)

Alternately, use SciPy for the bulk of the calculation, then re-do the edges the way you want.

It would probably be a good idea to check whether this takes enough time that it's worth optimizing before going too far, though (as you may have noticed, I am not good at that).

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

With the testing data I have it's almost instantaneous and, as I moved other topics, I don't really have something bigger to quickly try it on. I suggest we leave it that way for now and other users might open another issue if when need to speed it up. Would that be alright?

Copy link
Contributor

@DWesl DWesl Aug 1, 2024

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

You picked the same name as elsewhere in MetPy, though the edge handling is again different from what you do here (they do not smooth close to the edge) and from SciPy, and they do not use nanmean.

The bigger test data would likely be someone trying to find boundary layer height from model data somewhere, and it looks like your code is designed for a profile at a time, not arrays of profiles (either the N x Z that description implies or the Z x Y x X conventional in model output). It might be a simple matter to add an axis keyword argument to most functions and borrow the logic from the derivative functions to create the indexers (and possibly the vertical axis auto-detection for DataArrays), but that should probably be a follow-up PR in case it isn't.

"""Function that calculates the moving average with a given span.
The span is given in number of points on which the average is made.

Parameters
----------
val: array-like
Array of values
span: int
Span of the moving average. The higher the smoother

Returns
-------
smoothed_val: array-like
Array of smoothed values

See also
--------
[`bottleneck.move_mean`](https://bottleneck.readthedocs.io/en/latest/reference.html#bottleneck.move_mean),
[`scipy.ndimage.uniform_filter1d`](https://docs.scipy.org/doc/scipy/reference/generated/scipy.ndimage.uniform_filter1d.html#scipy.ndimage.uniform_filter1d),
[`pandas.DataFrame.rolling`](https://pandas.pydata.org/pandas-docs/stable/reference/api/pandas.DataFrame.rolling.html)
"""
n = len(val)
smoothed_val = deepcopy(val)
for i in range(n):
smoothed_val[i] = np.nanmean(val[i - min(span, i) : i + min(span, n - i)])

return smoothed_val


def bulk_richardson_number(
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Might close #628

height,
potential_temperature,
u,
v,
idxfoot: int = 0,
ustar=0 * units.meter_per_second,
):
r"""Calculate the bulk Richardson number.

See [VH96], eq. (3):

.. math:: Ri = (g/\theta) * \frac{(\Delta z)(\Delta \theta)}
{\left(\Delta u)^2 + (\Delta v)^2 + b(u_*)^2}

Parameters
----------
height : `pint.Quantity`
Altitude (metres above ground) of the points in the profile
potential_temperature : `pint.Quantity`
Potential temperature profile
u : `pint.Quantity`
Zonal wind profile
v : `pint.Quantity`
Meridional wind profile
idxfoot : int, optional
The index of the foot point (first trusted measure), defaults to 0.

Returns
-------
`pint.Quantity`
Bulk Richardson number profile
"""
if idxfoot == 0:
# Force the ground level to have null wind
Du = u
Dv = v
else:
Du = u - u[idxfoot]
Dv = v - v[idxfoot]

Dtheta = potential_temperature - potential_temperature[idxfoot]
Dz = height - height[idxfoot]

idx0 = Du**2 + Dv**2 + ustar**2 == 0
if idx0.sum() > 0:
bRi = np.ones_like(Dtheta) * np.nan * units.dimensionless
bRi[~idx0] = (
(mpconsts.g / potential_temperature[~idx0])
* (Dtheta[~idx0] * Dz[~idx0])
/ (Du[~idx0] ** 2 + Dv[~idx0] ** 2 + ustar**2)
)
else:
bRi = (
(mpconsts.g / potential_temperature)
* (Dtheta * Dz)
/ (Du**2 + Dv**2 + ustar**2)
)

return bRi


def blh_from_richardson_bulk(
height,
potential_temperature,
u,
v,
smoothingspan: int = 10,
idxfoot: int = 0,
bri_threshold=0.25 * units.dimensionless,
ustar=0.1 * units.meter_per_second,
):
"""Calculate atmospheric boundary layer height with the method of
bulk Richardson number.

It is the height where the bulk Richardson number exceeds a given threshold.
Well indicated for unstable boundary layers. See [VH96, Sei00, Col14, Guo16].

Parameters
----------
height : `pint.Quantity`
Altitude (metres above ground) of the points in the profile
potential_temperature : `pint.Quantity`
Potential temperature profile
u : `pint.Quantity`
Zonal wind profile
v : `pint.Quantity`
Meridional wind profile
smoothingspan : int, optional
The amount of smoothing (number of points in moving average)
idxfoot : int, optional
The index of the foot point (first trusted measure), defaults to 0.
bri_threshold : `pint.Quantity`, optional
Threshold to exceed to get boundary layer top. Defaults to 0.25
ustar : `pint.Quantity`, optional
Additional friction term in [VH96]. Defaluts to 0.

Returns
-------
blh : `pint.Quantity`
Boundary layer height estimation
"""
bRi = bulk_richardson_number(
height,
smooth(potential_temperature, smoothingspan),
smooth(u, smoothingspan),
smooth(v, smoothingspan),
idxfoot=idxfoot,
ustar=ustar,
)

height = height[~np.isnan(bRi)]
bRi = bRi[~np.isnan(bRi)]

if any(bRi > bri_threshold):
iblh = np.where(bRi > bri_threshold)[0][0]
blh = height[iblh]
else:
blh = np.nan * units.meter

return blh


def blh_from_parcel(
height,
potential_temperature,
smoothingspan: int = 5,
theta0=None,
):
"""Calculate atmospheric boundary layer height with the "parcel method"
(or "potential temperature threshold method").

It is the height where the potential temperature profile exceeds its
foot value. Well indicated for unstable boundary layers. See [Sei00, HL06, Col14].

Parameters
----------
height : `pint.Quantity`
Altitude (metres above ground) of the points in the profile
potential_temperature : `pint.Quantity`
Potential temperature profile
smoothingspan : int, optional
The amount of smoothing (number of points in moving average)
theta0 : `pint.Quantity`, optional
Value of theta at the foot point (skip unstruted points or add extra term). If not provided, theta[0] is taken.

Returns
-------
blh : `pint.Quantity`
Boundary layer height estimation
"""
potential_temperature = smooth(potential_temperature, smoothingspan)

if theta0 is None:
theta0 = potential_temperature[0]

if any(potential_temperature > theta0):
iblh = np.where(potential_temperature > theta0)[0][0]
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

So this looks like a potential temperature threshold method. I would prefer "exceeds" over "reaches" in the documentation, given the usual description of the convective boundary layer.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Fine by me. I will also add that it's only suited for unstable boundary layer, as suggested in the main comment. The name of the method might vary with the authors, "parcel method" is the one I have seen the most, but I can include other names (e.g. "potential temperature threshold method") in the doc.

Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

There might be a difference in expectations if the boundary layer is saturated (i.e. fog or fair-weather cumulus), but describing alternate names should avoid that.

blh = height[iblh]
else:
blh = np.nan * units.meter

return blh


def blh_from_concentration_gradient(
height,
concentration_profile,
smoothingspan: int = 5,
idxfoot: int = 0,
):
"""Calculate atmospheric boundary layer height from a concentration
profile (specific/relative humidity, aerosol backscatter, TKE..)

It is the height where the gradient of the concentration profile reaches a minimum.
Well indicated for stable boundary layers. See [Sei00, HL06, Col14].

Parameters
----------
height : `pint.Quantity`
Altitude (metres above ground) of the points in the profile
concentration_profile : `pint.Quantity`
Concentration profile (specific/relative humidity, aerosol backscatter, TKE..)
smoothingspan : int, optional
The amount of smoothing (number of points in moving average)
idxfoot : int, optional
The index of the foot point (first trusted measure), defaults to 0.

Returns
-------
blh : `pint.Quantity`
Boundary layer height estimation
"""
dcdz = mpcalc.first_derivative(smooth(concentration_profile, smoothingspan), x=height)
dcdz = dcdz[idxfoot:]
height = height[idxfoot:]
iblh = np.argmin(dcdz)

return height[iblh]


def blh_from_temperature_inversion(
height,
temperature,
smoothingspan: int = 5,
idxfoot: int = 0,
):
"""Calculate atmospheric boundary layer height from the inversion of
absolute temperature gradient

It is the height where the temperature gradient (absolute or potential) changes of sign.
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Out of curiosity, how well does this work for the convective boundary layer with potential temperature? Or, for that matter, with the nocturnal stable boundary layer with either?

I was expecting to see a threshold method on $\frac{d\theta}{dz}$.

Copy link
Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

From what I remember (my experience with this is now a bit old), this method is more suited for nocturnal stable boundary layers as it will track the end of the stable layer at the surface.

For convective boundary layer, the parcel method should be preferred to this method, as this method would gives underestimated height with a large variability.

The threshold of $\frac{d\theta}{dz}$ would be interesting to try too. I think I have seen it used in some works, just not in the ones I mention in this code.

Well indicated for stable boundary layers. See [Col14].

Parameters
----------
height : `pint.Quantity`
Altitude (metres above ground) of the points in the profile
humidity : `pint.Quantity`
Temperature (absolute or potential) profile
smoothingspan : int, optional
The amount of smoothing (number of points in moving average)
idxfoot : int, optional
The index of the foot point (first trusted measure), defaults to 0.

Returns
-------
blh : `pint.Quantity`
Boundary layer height estimation
"""
dTdz = mpcalc.first_derivative(smooth(temperature, smoothingspan), x=height)

dTdz = dTdz[idxfoot:]
height = height[idxfoot:]

if any(dTdz * dTdz[0] < 0):
iblh = np.where(dTdz * dTdz[0] < 0)[0][0]
blh = height[iblh]
else:
blh = np.nan * units.meter

return blh
66 changes: 66 additions & 0 deletions tests/test_boundarylayer.py
Original file line number Diff line number Diff line change
@@ -0,0 +1,66 @@
#!/usr/bin/python
# -*-coding:utf-8 -*-
"""Testing program for the MetPy boundary layer module"""

import numpy as np

Check notice

Code scanning / CodeQL

Unused import

Import of 'np' is not used.
import pandas as pd

import metpy.calc as mpcalc

Check notice

Code scanning / CodeQL

Module is imported with 'import' and 'import from'

Module 'metpy.calc' is imported with both 'import' and 'import from'.
from metpy.calc import boundarylayer
from metpy.cbook import get_test_data
from metpy.units import units

# SAMPLE DATA
# ===========
col_names = ["pressure", "height", "temperature", "dewpoint", "direction", "speed"]

df = pd.read_fwf(
get_test_data("may4_sounding.txt", as_file_obj=False),

Check warning

Code scanning / CodeQL

File is not always closed

File is opened but is not closed.
skiprows=5,
usecols=[0, 1, 2, 3, 6, 7],
names=col_names,
)

# Drop any rows with all NaN values for T, Td, winds
df = df.dropna(
subset=("temperature", "dewpoint", "direction", "speed"), how="all"
).reset_index(drop=True)

height = df["height"].values * units.metres
pressure = df["pressure"].values * units.hPa
temperature = df["temperature"].values * units.degC
dewpoint = df["dewpoint"].values * units.degC
wind_speed = df["speed"].values * units.knots
wind_dir = df["direction"].values * units.degrees

u, v = mpcalc.wind_components(wind_speed, wind_dir)
relative_humidity = mpcalc.relative_humidity_from_dewpoint(temperature, dewpoint)
potential_temperature = mpcalc.potential_temperature(pressure, temperature)
specific_humidity = mpcalc.specific_humidity_from_dewpoint(pressure, dewpoint)


# BOUNDARY LAYER HEIGHT ESTIMATIONS
# =================================

def test_blh_from_richardson_bulk():
blh = boundarylayer.blh_from_richardson_bulk(height, potential_temperature, u, v)
blh_true = 1397 * units.meter
assert blh == blh_true


def test_blh_from_parcel():
blh = boundarylayer.blh_from_parcel(height, potential_temperature)
blh_true = 610 * units.meter
assert blh == blh_true


def test_blh_from_concentration_gradient():
blh = boundarylayer.blh_from_concentration_gradient(height, specific_humidity)
blh_true = 1766 * units.meter
assert blh == blh_true


def test_blh_from_temperature_inversion():
blh = boundarylayer.blh_from_temperature_inversion(height, potential_temperature)
blh_true = 610 * units.meter
assert blh == blh_true