scipy.spatial.cKDTree.query#
- cKDTree.query(self, x, k=1, eps=0, p=2, distance_upper_bound=np.inf, workers=1)#
- Query the kd-tree for nearest neighbors - Parameters:
- xarray_like, last dimension self.m
- An array of points to query. 
- klist of integer or integer
- The list of k-th nearest neighbors to return. If k is an integer it is treated as a list of [1, … k] (range(1, k+1)). Note that the counting starts from 1. 
- epsnon-negative float
- Return approximate nearest neighbors; the k-th returned value is guaranteed to be no further than (1+eps) times the distance to the real k-th nearest neighbor. 
- pfloat, 1<=p<=infinity
- Which Minkowski p-norm to use. 1 is the sum-of-absolute-values “Manhattan” distance 2 is the usual Euclidean distance infinity is the maximum-coordinate-difference distance A finite large p may cause a ValueError if overflow can occur. 
- distance_upper_boundnonnegative float
- 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. - Changed in version 1.9.0: The “n_jobs” argument was renamed “workers”. The old name “n_jobs” was deprecated in SciPy 1.6.0 and was removed in SciPy 1.9.0. 
 
- Returns:
- darray of floats
- The distances to the nearest neighbors. If - xhas shape- tuple+(self.m,), then- dhas shape- tuple+(k,). When k == 1, the last dimension of the output is squeezed. Missing neighbors are indicated with infinite distances.
- indarray of ints
- The index of each neighbor in - self.data. If- xhas shape- tuple+(self.m,), then- ihas shape- tuple+(k,). When k == 1, the last dimension of the output is squeezed. Missing neighbors are indicated with- self.n.
 
 - Notes - If the KD-Tree is periodic, the position - xis wrapped into the box.- When the input k is a list, a query for arange(max(k)) is performed, but only columns that store the requested values of k are preserved. This is implemented in a manner that reduces memory usage. - Examples - >>> import numpy as np >>> from scipy.spatial import cKDTree >>> x, y = np.mgrid[0:5, 2:8] >>> tree = cKDTree(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]]