xrspatial.aspect.aspect#

xrspatial.aspect.aspect(agg: xarray.core.dataarray.DataArray, name: Optional[str] = 'aspect') xarray.core.dataarray.DataArray[source]#

Calculates the aspect value of an elevation aggregate.

Calculates, for all cells in the array, the downward slope direction of each cell based on the elevation of its neighbors in a 3x3 grid. The value is measured clockwise in degrees with 0 (due north), and 360 (again due north). Values along the edges are not calculated.

Direction of the aspect can be determined by its value: From 0 to 22.5: North From 22.5 to 67.5: Northeast From 67.5 to 112.5: East From 112.5 to 157.5: Southeast From 157.5 to 202.5: South From 202.5 to 247.5: West From 247.5 to 292.5: Northwest From 337.5 to 360: North

Note that values of -1 denote flat areas.

Parameters
  • agg (xarray.DataArray) – 2D NumPy, CuPy, or Dask with NumPy-backed xarray DataArray of elevation values.

  • name (str, default='aspect') – Name of ouput DataArray.

Returns

aspect_agg – 2D aggregate array of calculated aspect values. All other input attributes are preserved.

Return type

xarray.DataArray of the same type as agg

References

Examples

Aspect works with NumPy backed xarray DataArray .. sourcecode:: python

>>> import numpy as np
>>> import xarray as xr
>>> from xrspatial import aspect
>>> data = np.array([
    [1, 1, 1, 1, 1],
    [1, 1, 1, 2, 0],
    [1, 1, 1, 0, 0],
    [4, 4, 9, 2, 4],
    [1, 5, 0, 1, 4],
    [1, 5, 0, 5, 5]
], dtype=np.float32)
>>> raster = xr.DataArray(data, dims=['y', 'x'], name='raster')
>>> print(raster)
<xarray.DataArray 'raster' (y: 6, x: 5)>
array([[1., 1., 1., 1., 1.],
       [1., 1., 1., 2., 0.],
       [1., 1., 1., 0., 0.],
       [4., 4., 9., 2., 4.],
       [1., 5., 0., 1., 4.],
       [1., 5., 0., 5., 5.]])
Dimensions without coordinates: y, x
>>> aspect_agg = aspect(raster)
>>> print(aspect_agg)
<xarray.DataArray 'aspect' (y: 6, x: 5)>
array([[ nan,  nan        ,   nan       ,   nan       , nan],
       [ nan,  -1.        ,   225.      ,   135.      , nan],
       [ nan, 343.61045967,   8.97262661,  33.69006753, nan],
       [ nan, 307.87498365,  71.56505118,  54.46232221, nan],
       [ nan, 191.30993247, 144.46232221, 255.96375653, nan],
       [ nan,  nan        ,   nan       ,   nan       , nan]])
Dimensions without coordinates: y, x

Aspect works with Dask with NumPy backed xarray DataArray .. sourcecode:: python

>>> import dask.array as da
>>> data_da = da.from_array(data, chunks=(3, 3))
>>> raster_da = xr.DataArray(data_da, dims=['y', 'x'], name='raster_da')
>>> print(raster_da)
<xarray.DataArray 'raster' (y: 6, x: 5)>
dask.array<array, shape=(6, 5), dtype=int64, chunksize=(3, 3), chunktype=numpy.ndarray>
Dimensions without coordinates: y, x
>>> aspect_da = aspect(raster_da)
>>> print(aspect_da)
<xarray.DataArray 'aspect' (y: 6, x: 5)>
dask.array<_trim, shape=(6, 5), dtype=float32, chunksize=(3, 3), chunktype=numpy.ndarray>
Dimensions without coordinates: y, x
>>> print(aspect_da.compute())  # compute the results
<xarray.DataArray 'aspect' (y: 6, x: 5)>
array([[ nan,  nan        ,   nan       ,   nan       , nan],
       [ nan,  -1.        ,   225.      ,   135.      , nan],
       [ nan, 343.61045967,   8.97262661,  33.69006753, nan],
       [ nan, 307.87498365,  71.56505118,  54.46232221, nan],
       [ nan, 191.30993247, 144.46232221, 255.96375653, nan],
       [ nan,  nan        ,   nan       ,   nan       , nan]])
Dimensions without coordinates: y, x

Aspect works with CuPy backed xarray DataArray. Make sure you have a GPU and CuPy installed to run this example. .. sourcecode:: python

>>> import cupy
>>> data_cupy = cupy.asarray(data)
>>> raster_cupy = xr.DataArray(data_cupy, dims=['y', 'x'])
>>> aspect_cupy = aspect(raster_cupy)
>>> print(type(aspect_cupy.data))
<class 'cupy.core.core.ndarray'>
>>> print(aspect_cupy)
<xarray.DataArray 'aspect' (y: 6, x: 5)>
array([[       nan,       nan,        nan,        nan,        nan],
       [       nan,       -1.,       225.,       135.,        nan],
       [       nan, 343.61047,   8.972626,  33.690067,        nan],
       [       nan, 307.87497,  71.56505 ,  54.462322,        nan],
       [       nan, 191.30994, 144.46233 ,  255.96376,        nan],
       [       nan,       nan,        nan,        nan,        nan]],
      dtype=float32)
Dimensions without coordinates: y, x