act.retrievals.generic_sobel_cbh#

act.retrievals.generic_sobel_cbh(ds, variable=None, height_dim=None, var_thresh=None, fill_na=None, return_thresh=False, filter_type='uniform', edge_thresh=5.0)[source]#

Function for calculating cloud base height from lidar/radar data using a basic sobel filter and thresholding. Note, this was not initially based on any published work, but a lit review indicates that there have been similar methods employed to detect boundary layer heights.

NOTE: The returned variable now appends the field name of the data used to generate the CBH as part of the variable name. cbh_sobel_[varname]

Parameters:
  • ds (ACT xarray.Dataset) – ACT xarray dataset where data are stored.

  • variable (string) – Variable on which to process.

  • height_dim (string) – Height variable to use for CBH values.

  • var_thresh (float) – Thresholding for variable if needed.

  • fill_na (float) – Value to fill nans with in DataArray if any.

  • filter_type (string) – Currently the only option is for uniform filtering. uniform: Apply uniform filtering after the sobel filter? Applies a standard area of 3x3 filtering None: Excludes the filtering

  • edge_thresh (float) – Threshold value for finding the edge after the sobel filtering. If the signal is not strong, this may need to be lowered

Returns:

new_ds (ACT xarray.Dataset) – ACT xarray dataset with cbh values included as variable.

Examples

In testing on the ARM KAZR and MPL data, the following methods tended to work best for thresholding/corrections/etc.

kazr = act.retrievals.cbh.generic_sobel_cbh(
    kazr, variable="reflectivity_copol", height_dim="range", var_thresh=-10.0
)

mpl = act.corrections.mpl.correct_mpl(mpl)
mpl.range_bins.values = mpl.height.values[0, :] * 1000.0
mpl.range_bins.attrs["units"] = "m"
mpl["signal_return_co_pol"].values[:, 0:10] = 0.0
mpl = act.retrievals.cbh.generic_sobel_cbh(
    mpl,
    variable="signal_return_co_pol",
    height_dim="range_bins",
    var_thresh=10.0,
    fill_na=0.0,
)

ceil = act.retrievals.cbh.generic_sobel_cbh(
    ceil,
    variable="backscatter",
    height_dim="range",
    var_thresh=1000.0,
    fill_na=0,
)