scipy.spatial.KDTree.query#
- KDTree.query(x, k=1, eps=0, p=2, distance_upper_bound=inf, workers=1)[source]#
Query the kd-tree for nearest neighbors.
- Parameters:
- xarray_like, last dimension self.m
An array of points to query.
- kint or Sequence[int], optional
Either the number of nearest neighbors to return, or a list of the k-th nearest neighbors to return, starting from 1.
- epsnonnegative float, optional
Return approximate nearest neighbors; the kth returned value is guaranteed to be no further than (1+eps) times the distance to the real kth nearest neighbor.
- pfloat, 1<=p<=infinity, optional
Which Minkowski p-norm to use. 1 is the sum-of-absolute-values distance (“Manhattan” distance). 2 is the usual Euclidean distance. infinity is the maximum-coordinate-difference distance. A large, finite p may cause a ValueError if overflow can occur.
- distance_upper_boundnonnegative float, optional
Return only neighbors within this distance. This is used to prune tree searches, so if you are doing a series of nearest-neighbor queries, it may help to supply the distance to the nearest neighbor of the most recent point.
- workersint, optional
Number of workers to use for parallel processing. If -1 is given all CPU threads are used. Default: 1.
New in version 1.6.0.
- Returns:
- dfloat or array of floats
The distances to the nearest neighbors. If
x
has shapetuple+(self.m,)
, thend
has shapetuple+(k,)
. When k == 1, the last dimension of the output is squeezed. Missing neighbors are indicated with infinite distances. Hits are sorted by distance (nearest first).Changed in version 1.9.0: Previously if
k=None
, then d was an object array of shapetuple
, containing lists of distances. This behavior has been removed, usequery_ball_point
instead.- iinteger or array of integers
The index of each neighbor in
self.data
.i
is the same shape as d. Missing neighbors are indicated withself.n
.
Examples
>>> import numpy as np >>> from scipy.spatial import KDTree >>> x, y = np.mgrid[0:5, 2:8] >>> tree = KDTree(np.c_[x.ravel(), y.ravel()])
To query the nearest neighbours and return squeezed result, use
>>> dd, ii = tree.query([[0, 0], [2.2, 2.9]], k=1) >>> print(dd, ii, sep='\n') [2. 0.2236068] [ 0 13]
To query the nearest neighbours and return unsqueezed result, use
>>> dd, ii = tree.query([[0, 0], [2.2, 2.9]], k=[1]) >>> print(dd, ii, sep='\n') [[2. ] [0.2236068]] [[ 0] [13]]
To query the second nearest neighbours and return unsqueezed result, use
>>> dd, ii = tree.query([[0, 0], [2.2, 2.9]], k=[2]) >>> print(dd, ii, sep='\n') [[2.23606798] [0.80622577]] [[ 6] [19]]
To query the first and second nearest neighbours, use
>>> dd, ii = tree.query([[0, 0], [2.2, 2.9]], k=2) >>> print(dd, ii, sep='\n') [[2. 2.23606798] [0.2236068 0.80622577]] [[ 0 6] [13 19]]
or, be more specific
>>> dd, ii = tree.query([[0, 0], [2.2, 2.9]], k=[1, 2]) >>> print(dd, ii, sep='\n') [[2. 2.23606798] [0.2236068 0.80622577]] [[ 0 6] [13 19]]