xrspatial.proximity.proximity

xrspatial.proximity.proximity(raster: xarray.core.dataarray.DataArray, x: str = 'x', y: str = 'y', target_values: list = [], max_distance: float = inf, distance_metric: str = 'EUCLIDEAN') xarray.core.dataarray.DataArray[source]

Computes the proximity of all pixels in the image to a set of pixels in the source image based on a distance metric.

This function attempts to compute the proximity of all pixels in the image to a set of pixels in the source image. The following options are used to define the behavior of the function. By default all non-zero pixels in raster.values will be considered the “target”, and all proximities will be computed in pixels. Note that target pixels are set to the value corresponding to a distance of zero.

Proximity support NumPy backed, and Dask with NumPy backed xarray DataArray. The return values of proximity are of the same type as the input type. If input raster is a NumPy-backed DataArray, the result is NumPy-backed. If input raster is a Dask-backed DataArray, the result is Dask-backed.

The implementation for NumPy-backed is ported from GDAL, which uses a dynamic programming approach to identify nearest target of a pixel from its surrounding neighborhood in a 3x3 window. The implementation for Dask-backed uses dask.map_overlap to compute proximity chunk by chunk by expanding the chunk’s borders to cover the max_distance.

Parameters
  • raster (xr.DataArray) – 2D array image with raster.shape = (height, width).

  • x (str, default='x') – Name of x-coordinates.

  • y (str, default='y') – Name of y-coordinates.

  • target_values (list) – Target pixel values to measure the distance from. If this option is not provided, proximity will be computed from non-zero pixel values.

  • max_distance (float, default=np.inf) –

    The maximum distance to search. Proximity distances greater than this value will be set to NaN. Should be given in the same distance unit as input. For example, if input raster is in lat-lon and distances between points within the raster is calculated using Euclidean distance metric, max_distance should also be provided in lat-lon unit. If using Great Circle distance metric, and thus all distances is in km, max_distance should also be provided in kilometer unit.

    When scaling with Dask, whether the function scales well depends on the max_distance value. If max_distance is infinite by default, this function only works on a single machine. It should scale well, however, if max_distance is relatively small compared to the maximum possible distance in two arbitrary points in the input raster. Note that if max_distance is equal or larger than the max possible distance between 2 arbitrary points in the input raster, the input data array will be rechunked.

  • distance_metric (str, default='EUCLIDEAN') – The metric for calculating distance between 2 points. Valid distance metrics are: ‘EUCLIDEAN’, ‘GREAT_CIRCLE’, and ‘MANHATTAN’.

Returns

proximity_agg – 2D array of proximity values. All other input attributes are preserved.

Return type

xr.DataArray of same type as raster

References

Examples

import numpy as np
import xarray as xr
from xrspatial import proximity

data = np.array([
    [0., 0., 0., 0., 0.],
    [0., 0., 0., 1., 0.],
    [0., 0., 0., 0., 0.],
    [0., 0., 0., 0., 0.],
    [0., 0., 0., 0., 0.]
])

n, m = data.shape
raster = xr.DataArray(data, dims=['y', 'x'], name='raster')
raster['y'] = np.arange(n)[::-1]
raster['x'] = np.arange(m)

proximity_agg = proximity(raster)
>>> raster
<xarray.DataArray 'raster' (y: 5, x: 5)>
array([[0., 0., 0., 0., 0.],
       [0., 0., 0., 1., 0.],
       [0., 0., 0., 0., 0.],
       [0., 0., 0., 0., 0.],
       [1., 0., 0., 0., 0.]])
Coordinates:
  * y        (y) int64 4 3 2 1 0
  * x        (x) int64 0 1 2 3 4

>>> proximity_agg
<xarray.DataArray (y: 5, x: 5)>
array([[3.1622777, 2.236068 , 1.4142135, 1.       , 1.4142135],
       [3.       , 2.       , 1.       , 0.       , 1.       ],
       [3.1622777, 2.236068 , 1.4142135, 1.       , 1.4142135],
       [3.6055512, 2.828427 , 2.236068 , 2.       , 2.236068 ],
       [4.2426405, 3.6055512, 3.1622777, 3.       , 3.1622777]],
      dtype=float32)
Coordinates:
  * y        (y) int64 4 3 2 1 0
  * x        (x) int64 0 1 2 3 4