xrspatial.slope.slope#

xrspatial.slope.slope(agg: DataArray, name: str | None = 'slope', method: str = 'planar', z_unit: str = 'meter', boundary: str = 'nan') DataArray[source]#

Returns slope of input aggregate in degrees.

Parameters:
  • agg (xr.DataArray or xr.Dataset) – 2D array of elevation data. If a Dataset is passed, the operation is applied to each data variable independently.

  • name (str, default='slope') – Name of output DataArray.

  • method (str, default='planar') – 'planar' uses the classic Horn algorithm with uniform cell size. 'geodesic' converts cells to Earth-Centered Earth-Fixed (ECEF) coordinates and fits a 3D plane, yielding accurate results for geographic (lat/lon) coordinate systems.

  • z_unit (str, default='meter') – Unit of the elevation values. Only used when method='geodesic'. Accepted values: 'meter', 'foot', 'kilometer', 'mile' (and common aliases).

  • boundary (str, default='nan') – How to handle edges where the kernel extends beyond the raster. 'nan' — fill missing neighbours with NaN (default). 'nearest' — repeat edge values. 'reflect' — mirror at boundary. 'wrap' — periodic / toroidal.

Returns:

slope_agg – If agg is a DataArray, returns a DataArray of the same type. If agg is a Dataset, returns a Dataset with slope computed for each data variable. 2D array of slope values in degrees. The output attrs['units'] is set to 'degrees'; all other input attributes are preserved.

Return type:

xr.DataArray or xr.Dataset

Notes

The 'planar' method uses the coordinate spacing directly as the cell size. If the coordinates are in degrees (lat/lon) but the elevation values are in meters, the result is wrong by orders of magnitude. When this mismatch is detected, a UserWarning is emitted suggesting you reproject to a projected CRS or use method='geodesic'. The 'planar' method also raises a ValueError if the cell size on either axis is zero, negative, or non-finite, since those values cannot produce a valid slope.

References

Examples

>>> import numpy as np
>>> import xarray as xr
>>> from xrspatial import slope
>>> data = np.array([
...     [0, 0, 0, 0, 0],
...     [0, 0, 0, -1, 2],
...     [0, 0, 0, 0, 1],
...     [0, 0, 0, 5, 0]])
>>> agg = xr.DataArray(data)
>>> slope_agg = slope(agg)
>>> slope_agg
<xarray.DataArray 'slope' (dim_0: 4, dim_1: 5)>
array([[      nan,       nan,       nan,       nan,       nan],
       [      nan,  0.      , 14.036243, 32.512516,       nan],
       [      nan,  0.      , 42.031113, 53.395725,       nan],
       [      nan,       nan,       nan,       nan,       nan]],
      dtype=float32)
Dimensions without coordinates: dim_0, dim_1