xdem.spstats

Spatial statistical tools to estimate uncertainties related to DEMs

Functions

cov

xdem.spstats.cov(h, crange, model='Sph', psill=1.0, kappa=0.5, nugget=0)[source]

Covariance function based on variogram function (COV = STD - VGM)

Parameters
  • h (float) – spatial lag

  • crange (float) – correlation range

  • model (str) – model

  • psill (float) – partial sill

  • kappa (float) – smoothing parameter for Exp Class

  • nugget (float) – nugget

Return type

Callable

Returns

covariance function

create_circular_mask

xdem.spstats.create_circular_mask(shape, center=None, radius=None)[source]

Create circular mask on a raster, defaults to the center of the array and it’s half width

Parameters
  • shape – shape of array

  • center – center

  • radius – radius

Returns

create_ring_mask

xdem.spstats.create_ring_mask(shape, center=None, in_radius=0.0, out_radius=None)[source]

Create ring mask on a raster, defaults to the center of the array and a circle mask of half width of the array

Parameters
  • shape – shape of array

  • center – center

  • in_radius – inside radius

  • out_radius – outside radius

Returns

distance_latlon

xdem.spstats.distance_latlon(tup1, tup2, earth_rad=6373000)[source]

Distance between two lat/lon coordinates projected on a spheroid ref: https://stackoverflow.com/questions/19412462/getting-distance-between-two-points-based-on-latitude-longitude :type tup1: tuple :param tup1: lon/lat coordinates of first point :type tup2: tuple :param tup2: lon/lat coordinates of second point :type earth_rad: float :param earth_rad: radius of the earth in meters

Return type

float

Returns

distance

double_sum_covar

xdem.spstats.double_sum_covar(list_tuple_errs, corr_ranges, list_area_tot, list_lat, list_lon, nproc=1)[source]

Double sum of covariances for propagating multi-range correlated errors between disconnected spatial ensembles

Parameters
  • list_tuple_errs – list of tuples of correlated errors by range, by ensemble

  • corr_ranges – list of correlation ranges

  • list_area_tot – list of areas of ensembles

  • list_lat – list of center latitude of ensembles

  • list_lon – list of center longitude of ensembles

  • nproc – number of cores to use for multiprocessing

Returns

sum of covariances

exact_neff_sphsum_circular

xdem.spstats.exact_neff_sphsum_circular(area, crange1, psill1, crange2, psill2)[source]

Number of effective samples derived from exact integration of sum of spherical variogram models over a circular area. The number of effective samples serves to convert between standard deviation/partial sills and standard error over the area. If SE is the standard error, SD the standard deviation and N_eff the number of effective samples, we have: SE = SD / sqrt(N_eff) => N_eff = SD^2 / SE^2 => N_eff = (PS1 + PS2)/SE^2 where PS1 and PS2 are the partial sills estimated from the variogram models, and SE is estimated by integrating the variogram models with parameters PS1/PS2 and R1/R2 where R1/R2 are the correlation ranges. Source: Rolstad et al. (2009), appendix: http://dx.doi.org/10.3189/002214309789470950

Parameters
  • area (float) – circular area

  • crange1 (float) – range of short-range variogram model

  • psill1 (float) – partial sill of short-range variogram model

  • crange2 (float) – range of long-range variogram model

  • psill2 (float) – partial sill of long-range variogram model

Return type

float

Returns

number of effective samples

fit_model_sum_vgm

xdem.spstats.fit_model_sum_vgm(list_model, emp_vgm_df)[source]

Fit a multi-range variogram model to an empirical variogram, weighted based on sampling and elevation errors

Parameters
  • list_model – list of variogram models to sum for the fit: from short-range to long-ranges

  • emp_vgm_df – empirical variogram

Returns

modelled variogram

get_empirical_variogram

xdem.spstats.get_empirical_variogram(dh, coords, **kwargs)[source]

Get empirical variogram from skgstat.Variogram object

Parameters
  • dh (ndarray) – elevation differences

  • coords (ndarray) – coordinates

Return type

DataFrame

Returns

empirical variogram (variance, lags, counts)

integrate_fun

xdem.spstats.integrate_fun(fun, low_b, upp_b)[source]

Numerically integrate function between upper and lower bounds :type fun: Callable :param fun: function :type low_b: float :param low_b: lower bound :type upp_b: float :param upp_b: upper bound

Return type

float

Returns

integral

kernel_sph

xdem.spstats.kernel_sph(xi, x0, a1)[source]

Spherical kernel :type xi: float :param xi: position of first point :type x0: float :param x0: position of second point :type a1: float :param a1: range of kernel

Return type

float

Returns

covariance between the two points

nd_binning

xdem.spstats.nd_binning(values, list_var, list_var_names=typing.Iterable[str], list_var_bins=None, statistics=['count', <function nanmedian>, <function nmad>], list_ranges=None)[source]

N-dimensional binning of values according to one or several explanatory variables. Values input is a (N,) array and variable input is a list of flattened arrays of similar dimensions (N,). For more details on the format of input variables, see documentation of scipy.stats.binned_statistic_dd.

Parameters
  • values (ndarray) – values array (N,)

  • list_var (Iterable[ndarray]) – list (L) of explanatory variables array (N,)

  • list_var_names – list (L) of names of the explanatory variables

  • list_var_bins (Union[int, Iterable[Iterable], None]) – count, or list (L) of counts or custom bin edges for the explanatory variables; defaults to 10 bins

  • statistics (Iterable[Union[str, Callable, None]]) – list (X) of statistics to be computed; defaults to count, median and nmad

  • list_ranges (Optional[Iterable[Sequence]]) – list (L) of minimum and maximum ranges to bin the explanatory variables; defaults to min/max of the data

Return type

DataFrame

Returns

neff_circ

xdem.spstats.neff_circ(area, list_vgm)[source]

Number of effective samples derived from numerical integration for any sum of variogram models a circular area (generalization of Rolstad et al. (2009): http://dx.doi.org/10.3189/002214309789470950) The number of effective samples N_eff serves to convert between standard deviation/partial sills and standard error over the area: SE = SD / sqrt(N_eff) if SE is the standard error, SD the standard deviation.

Parameters
  • area – area

  • list_vgm – variogram functions to sum

Returns

number of effective samples

neff_rect

xdem.spstats.neff_rect(area, width, crange1, psill1, model1='Sph', crange2=None, psill2=None, model2=None)[source]

Number of effective samples derived from numerical integration for a sum of 2 variogram functions over a rectangular area

Parameters
  • area (float) – area

  • width (float) – width of rectangular area

  • crange1 (float) – correlation range of first variogram

  • psill1 (float) – partial sill of first variogram

  • model1 (str) – model of first variogram

  • crange2 (Optional[float]) – correlation range of second variogram

  • psill2 (Optional[float]) – partial sill of second variogram

  • model2 (Optional[str]) – model of second variogram

Return type

float

Returns

number of effective samples

part_covar_sum

xdem.spstats.part_covar_sum(argsin)[source]

Multiprocessing wrapper for covariance summing :type argsin: tuple :param argsin: Tupled argument for covariance calculation

Return type

float

Returns

covariance sum

patches_method

xdem.spstats.patches_method(dh, mask, gsd, area_size, perc_min_valid=80.0, patch_shape='circular', nmax=1000, verbose=False)[source]

Patches method for empirical estimation of the standard error over an integration area

Parameters
  • dh – elevation differences

  • mask – mask of sampled terrain

  • gsd – ground sampling distance

  • area_size – size of integration area

  • perc_min_valid – minimum valid area in the patch

  • patch_shape – shape of patch [‘circular’ or ‘rectangular’]

  • nmax – maximum number of patch to sample

#TODO: add overlap option?

Returns

tile, mean, median, std and count of each patch

plot_vgm

xdem.spstats.plot_vgm(df, fit_fun=None)[source]

random_subset

xdem.spstats.random_subset(dh, coords, nsamp)[source]

Subsampling of elevation differences with random coordinates

Parameters
  • dh – elevation differences

  • coords – coordinates

  • nsamp – number of sammples for subsampling

Returns

subsets of dh and coords

ring_subset

xdem.spstats.ring_subset(dh, coords, inside_radius=0, outside_radius=0)[source]

Subsampling of elevation differences within a ring/disk (to sample points at similar pairwise distances)

Parameters
  • dh – elevation differences

  • coords – coordinates

  • inside_radius – radius of inside ring disk in pixels

  • outside_radius – radius of outside ring disk in pixels

Returns

subsets of dh and coords

sample_multirange_empirical_variogram

xdem.spstats.sample_multirange_empirical_variogram(dh, gsd=None, coords=None, nsamp=10000, range_list=None, nrun=1, nproc=1, **kwargs)[source]

Wrapper to sample multi-range empirical variograms from the data.

If no option is passed, a varying binning is used with adapted ranges and data subsampling

Parameters
  • dh (ndarray) – elevation differences

  • gsd (Optional[float]) – ground sampling distance (if array is 2D on structured grid)

  • coords (Optional[ndarray]) – coordinates, to be used only with a flattened elevation differences array and passed as an array of he pairs of coordinates: one dimension equal to two and the other to that of the flattened elevation differences

  • range_list (Optional[list]) – successive ranges with even binning

  • nsamp (int) – number of samples to randomly draw from the elevation differences

  • nrun (int) – number of samplings

  • nproc (int) – number of processing cores

Return type

DataFrame

Returns

empirical variogram (variance, lags, counts)

std_err

xdem.spstats.std_err(std, neff)[source]

Standard error

Parameters
  • std (float) – standard deviation

  • neff (float) – number of effective samples

Return type

float

Returns

standard error

std_err_finite

xdem.spstats.std_err_finite(std, neff_tot, neff)[source]

Standard error of subsample of a finite ensemble

Parameters
  • std (float) – standard deviation

  • neff_tot (float) – maximum number of effective samples

  • neff (float) – number of effective samples

Return type

float

Returns

standard error

vgm

xdem.spstats.vgm(h, crange, model='Sph', psill=1.0, kappa=0.5, nugget=0)[source]

Compute variogram model function (Spherical, Exponential, Gaussian or Exponential Class)

Parameters
  • h (float) – spatial lag

  • crange (float) – correlation range

  • model (str) – model

  • psill (float) – partial sill

  • kappa (float) – smoothing parameter for Exp Class

  • nugget (float) – nugget

Returns

variogram function

wrapper_get_empirical_variogram

xdem.spstats.wrapper_get_empirical_variogram(argdict, **kwargs)[source]

Multiprocessing wrapper for get_empirical_variogram

Parameters

argdict (dict) – Keyword argument to pass to get_empirical_variogram()

Return type

DataFrame

Returns

empirical variogram (variance, lags, counts)