Distance#

Point#

pyinterpolate.distance.point.point_distance(points: ArrayLike, other: ArrayLike, metrics: str = 'euclidean') ndarray[source]

Calculates the Euclidean distance from set of points to other set of points.

Parameters:
pointsarray

Spatial coordinates.

otherarray

Other array with spatial coordinates.

metricsstr, default = ‘euclidean’

Metrics used to calculate distance. See scipy.spatial.distance.cdist for more details.

Returns:
distancesarray

Distances matrix. Row index = points point index, and column index = other point index.

Notes

The function creates array of size MxN, where M = number of points and N = number of other. Large arrays may cause memory errors.

Examples

>>> points = [(0, 0), (0, 1), (0, 2)]
>>> other = [(2, 2), (3, 3)]
>>> distances = point_distance(points=points, other=other)
>>> print(distances)
[[2.82842712 4.24264069]
 [2.23606798 3.60555128]
 [2.         3.16227766]]

Block#

pyinterpolate.distance.block.calc_block_to_block_distance(ps_blocks: DataFrame | PointSupport, lon_col_name=None, lat_col_name=None, val_col_name=None, block_id_col_name=None) DataFrame[source]

Function calculates distances between the blocks’ point supports.

Parameters:
ps_blocksUnion[pd.DataFrame, PointSupport]
The point support of polygons.
  • DataFrame, then lon_col_name, lat_col_name, val_col_name, block_id_col_name parameters must be provided,

  • PointSupport.

lon_col_nameoptional

Longitude or x coordinate.

lat_col_nameoptional

Latitude or y coordinate.

val_col_nameoptional

The point support values column.

block_id_col_nameoptional

Column with block names / indexes.

Returns:
block_distancesDataFrame

Indexes and columns representing unique blocks, values in cells are weighted distances between blocks.

Raises:
AttributeError

Blocks are provided as DataFrame but column names were not given.

Notes

The weighted distance between blocks is derived from the equation 3 given in publication [1] from References.

\[d(v_{a}, v_{b})=\frac{1}{\sum_{s=1}^{P_{a}} \sum_{s'=1}^{P_{b}} n(u_{s}) n(u_{s'})} * \sum_{s=1}^{P_{a}} \sum_{s'=1}^{P_{b}} n(u_{s})n(u_{s'})||u_{s}-u_{s'}||\]
where:
  • \(P_{a}\) and \(P_{b}\): number of points \(u_{s}\) and \(u_{s'}\) used to discretize the two units \(v_{a}\) and \(v_{b}\)

  • \(n(u_{s})\) and \(n(u_{s'})\) - population size in the cells \(u_{s}\) and \(u_{s'}\)

References

[1]

Goovaerts, P. Kriging and Semivariogram Deconvolution in the Presence of Irregular Geographical Units. Math Geosci 40, 101–128 (2008). https://doi.org/10.1007/s11004-007-9129-1

Examples

>>> import os
>>> import geopandas as gpd
>>> from pyinterpolate import (
...     calc_block_to_block_distance,
...     Blocks,
...     ExperimentalVariogram,
...     PointSupport,
...     TheoreticalVariogram
... )
>>>
>>>
>>> FILENAME = 'cancer_data.gpkg'
>>> LAYER_NAME = 'areas'
>>> DS = gpd.read_file(FILENAME, layer=LAYER_NAME)
>>> AREA_VALUES = 'rate'
>>> AREA_INDEX = 'FIPS'
>>> AREA_GEOMETRY = 'geometry'
>>> PS_LAYER_NAME = 'points'
>>> PS_VALUES = 'POP10'
>>> PS_GEOMETRY = 'geometry'
>>> PS = gpd.read_file(FILENAME, layer=PS_LAYER_NAME)
>>>
>>> CANCER_DATA = {
...    'ds': DS,
...    'index_column_name': AREA_INDEX,
...    'value_column_name': AREA_VALUES,
...    'geometry_column_name': AREA_GEOMETRY
... }
>>> POINT_SUPPORT_DATA = {
...     'ps': PS,
...     'value_column_name': PS_VALUES,
...     'geometry_column_name': PS_GEOMETRY
... }
>>> BLOCKS = Blocks(**CANCER_DATA)
>>> indexes = BLOCKS.block_indexes
>>>
>>> PS = PointSupport(
...     points=POINT_SUPPORT_DATA['ps'],
...     ps_blocks=BLOCKS,
...     points_value_column=POINT_SUPPORT_DATA['value_column_name'],
...     points_geometry_column=POINT_SUPPORT_DATA['geometry_column_name']
... )
>>> distances = calc_block_to_block_distance(PS)
>>> print(distances.head())
               42049.0        42039.0  ...       23009.0       23029.0
42049.0       0.000000   48885.422652  ...  9.780781e+05  1.061496e+06
42039.0   48885.422652       0.000000  ...  9.952927e+05  1.079528e+06
42085.0   91350.566468   52240.566481  ...  1.030860e+06  1.115642e+06
42073.0  120259.614440   78853.981609  ...  1.042438e+06  1.127643e+06
42007.0  152110.979345  109955.408914  ...  1.055898e+06  1.141472e+06