xdem.volume

Volume change calculation tools (aimed for glaciers).

Functions

calculate_hypsometry_area

xdem.volume.calculate_hypsometry_area(ddem_bins, ref_dem, pixel_size, timeframe='reference')[source]

Calculate the associated representative area of the given dDEM bins.

By default, the area bins will be representative of the mean timing between the reference and nonreference DEM: elevations = ref_dem - (h_vs_dh_funcion(ref_dem) / 2) This can be changed to either “reference”: elevations = ref_dem Or “nonreference”: elevations = ref_dem - h_vs_dh_function(ref_dem)

Parameters
  • ddem_bins – A Series or DataFrame of dDEM values. If a DataFrame is given, the column ‘value’ will be used.

  • ref_dem – The reference DEM. This should not have any NaNs.

  • pixel_size – The xy or (x, y) size of the reference DEM pixels in georeferenced coordinates.

  • timeframe – The time at which the area bins are representative. Choices: “reference”, “nonreference”, “mean”

Returns

The representative area within the given dDEM bins.

fit_hypsometric_bins_poly

xdem.volume.fit_hypsometric_bins_poly(hypsometric_bins, value_column='value', degree=3, iterations=1, count_threshold=None)[source]

Fit a polynomial to the hypsometric bins.

Parameters
  • hypsometric_bins (DataFrame) – Bins where nans will be interpolated.

  • value_column (str) – The name of the column in ‘hypsometric_bins’ to use as values.

  • degree (int) – The degree of the polynomial to use.

  • iterations (int) – The number of iterations to run. At each iteration, values with residuals larger than 3 times the residuals’ standard deviation are excluded.

  • count_threshold (Optional[int]) – Optional. A pixel count threshold to exclude during the curve fit (requires a ‘count’ column).

Return type

Series

Returns

Bins replaced by the polynomial fit.

get_regional_hypsometric_signal

xdem.volume.get_regional_hypsometric_signal(ddem, ref_dem, glacier_index_map, n_bins=20, verbose=False, min_coverage=0.05)[source]

Get the normalized regional hypsometric elevation change signal, read “the general shape of it”.

Parameters
  • ddem (Union[ndarray, MaskedArray]) – The dDEM to analyse.

  • ref_dem (Union[ndarray, MaskedArray]) – A void-free reference DEM.

  • glacier_index_map (ndarray) – An array glacier indices of the same shape as the previous inputs.

  • verbose (bool) – Show progress bar.

n_bins = 20 # TODO: This should be an argument. :type n_bins: int :param n_bins: The number of elevation bins to subdivide each glacier in.

Return type

DataFrame

Returns

A DataFrame of bin statistics, scaled by elevation and elevation change.

hypsometric_binning

xdem.volume.hypsometric_binning(ddem, ref_dem, bins=50.0, kind='fixed', aggregation_function=<function median>)[source]

Separate the dDEM in discrete elevation bins. The elevation bins will be calculated based on all ref_dem valid pixels. ddem may contain NaN/masked values over the same area, they will be excluded before the aggregation.

It is assumed that the dDEM is calculated as ‘ref_dem - dem’ (not ‘dem - ref_dem’).

Parameters
  • ddem (ndarray) – The dDEM as a 2D or 1D array.

  • ref_dem (ndarray) – The reference DEM as a 2D or 1D array.

  • bins (Union[float, ndarray]) – The bin size, count, or array, depending on the binning method (‘kind’).

  • kind (str) – The kind of binning to do. Choices: [‘fixed’, ‘count’, ‘quantile’, ‘custom’].

  • aggregation_function (Callable) – The function to aggregate the elevation values within a bin. Defaults to the median.

Return type

DataFrame

Returns

A Pandas DataFrame with elevation bins and dDEM statistics.

hypsometric_interpolation

xdem.volume.hypsometric_interpolation(voided_ddem, ref_dem, mask)[source]

Interpolate a dDEM using hypsometric interpolation within the given mask.

Using ref_dem, elevation bins of constant height (hard-coded to 50 m for now) are created. Gaps in voided-ddem, within the provided mask, are filled with the median dDEM value within that bin.

Parameters
  • voided_ddem (Union[ndarray, MaskedArray]) – A dDEM with voids (either an array with nans or a masked array).

  • ref_dem (Union[ndarray, MaskedArray]) – The reference DEM in the dDEM comparison.

  • mask (ndarray) – A mask to delineate the area that will be interpolated (True means hypsometric will be used).

Return type

MaskedArray

interpolate_hypsometric_bins

xdem.volume.interpolate_hypsometric_bins(hypsometric_bins, value_column='value', method='polynomial', order=3, count_threshold=None)[source]

Interpolate hypsometric bins using any valid Pandas interpolation technique.

NOTE: It will not extrapolate!

Parameters
  • hypsometric_bins (DataFrame) – Bins where nans will be interpolated.

  • value_column – The name of the column in ‘hypsometric_bins’ to use as values.

  • method – Any valid Pandas interpolation technique (e.g. ‘linear’, ‘polynomial’, ‘ffill’, ‘bfill’).

  • order – The polynomial order to use. Only used if method=’polynomial’.

  • count_threshold (Optional[int]) – Optional. A pixel count threshold to exclude during the curve fit (requires a ‘count’ column).

Return type

DataFrame

Returns

Bins interpolated with the chosen interpolation method.

linear_interpolation

xdem.volume.linear_interpolation(array, max_search_distance=10, extrapolate=False, force_fill=False)[source]

Interpolate a 2D array using rasterio’s fillnodata.

Parameters
  • array (Union[ndarray, MaskedArray]) – An array with NaNs or a masked array to interpolate.

  • max_search_distance (int) – The maximum number of pixels to search in all directions to find values to interpolate from. The default is 10.

  • extrapolate (bool) – if False, will remove pixels that have been extrapolated by fillnodata. Default is False.

  • force_fill (bool) – if True, will replace all remaining gaps by the median of all valid values. Default is False.

Return type

ndarray

Returns

A filled array with no NaNs

local_hypsometric_interpolation

xdem.volume.local_hypsometric_interpolation(voided_ddem, ref_dem, mask, min_coverage=0.2, count_threshold=1, nodata=- 9999, plot=False)[source]

Interpolate a dDEM using local hypsometric interpolation. The algorithm loops through each features in the vector file.

The dDEM is assumed to have been created as “voided_ddem = reference_dem - other_dem”.

Parameters
  • voided_ddem (Union[ndarray, MaskedArray]) – A dDEM with voids (either an array with nans or a masked array).

  • ref_dem (Union[ndarray, MaskedArray]) – The reference DEM in the dDEM comparison.

  • mask (ndarray) – A raster of same shape as voided_ddem and ref_dem, containing a diferent non-0 pixel value for each geometry on which to loop.

  • min_coverage (float) – Optional. The minimum coverage fraction to be considered for interpolation.

  • count_threshold (Optional[int]) – Optional. A pixel count threshold to exclude during the hypsometric curve fit.

  • nodata (Union[float, int]) – Optional. No data value to be used for the output masked_array.

  • plot (bool) – Set to True to display intermediate plots.

Return type

MaskedArray

Returns

A dDEM with gaps filled by applying a hypsometric interpolation for each geometry in mask, for areas filling the min_coverage criterion.

norm_regional_hypsometric_interpolation

xdem.volume.norm_regional_hypsometric_interpolation(voided_ddem, ref_dem, glacier_index_map, min_coverage=0.1, regional_signal=None, verbose=False, min_elevation_range=0.33, idealized_ddem=False)[source]

Interpolate missing values by scaling the normalized regional hypsometric signal to each glacier separately.

Only missing values are interpolated. The rest of the glacier’s values are fixed.

Parameters
  • voided_ddem (Union[ndarray, MaskedArray]) – The voided dDEM to fill NaNs in.

  • ref_dem (Union[ndarray, MaskedArray]) – A void-free reference DEM.

  • glacier_index_map (ndarray) – An array glacier indices of the same shape as the previous inputs.

  • min_coverage (float) – The minimum fractional coverage of a glacier to interpolate. Defaults to 10%.

  • regional_signal (Optional[DataFrame]) – A regional signal is already estimate. Otherwise one will be estimated.

  • verbose (bool) – Show progress bars.

  • min_elevation_range (float) – The minimum allowed min/max bin range to scale a signal from. Default: 1/3 of the elevation range needs to be present.

  • idealized_ddem (bool) – Replace observed glacier values with the hypsometric signal. Good for error assessments.

Raises

AssertionError – If ref_dem has voids.

Return type

ndarray

Returns

A dDEM where glacier’s that fit the min_coverage criterion are interpolated.