1
|
|
|
"""Implements core function nearest_neighbors used for AMD and PDD |
2
|
|
|
calculations. |
3
|
|
|
""" |
4
|
|
|
|
5
|
|
|
from typing import Tuple, Iterable, Iterator, Callable, List |
6
|
|
|
from itertools import product, tee, accumulate |
7
|
|
|
import functools |
8
|
|
|
|
9
|
|
|
import numba |
10
|
|
|
import numpy as np |
11
|
|
|
import numpy.typing as npt |
12
|
|
|
from scipy.spatial import KDTree |
13
|
|
|
from scipy.spatial.distance import cdist |
14
|
|
|
|
15
|
|
|
FloatArray = npt.NDArray[np.floating] |
16
|
|
|
IntArray = npt.NDArray[np.integer] |
17
|
|
|
|
18
|
|
|
__all__ = [ |
19
|
|
|
'nearest_neighbors', |
20
|
|
|
'nearest_neighbors_data', |
21
|
|
|
'nearest_neighbors_minval', |
22
|
|
|
'generate_concentric_cloud' |
23
|
|
|
] |
24
|
|
|
|
25
|
|
|
|
26
|
|
|
def nearest_neighbors( |
27
|
|
|
motif: FloatArray, cell: FloatArray, x: FloatArray, k: int |
28
|
|
|
) -> FloatArray: |
29
|
|
|
"""Find distances to ``k`` nearest neighbors in a periodic set for |
30
|
|
|
each point in ``x``. |
31
|
|
|
|
32
|
|
|
Given a periodic set described by ``motif`` and ``cell``, a query |
33
|
|
|
set of points ``x`` and an integer ``k``, find distances to the |
34
|
|
|
``k`` nearest neighbors in the periodic set for all points in |
35
|
|
|
``x``. Returns an array with shape (x.shape[0], k) of distances to |
36
|
|
|
the neighbors. This function only returns distances, see the |
37
|
|
|
function nearest_neighbors_data() to also get the point cloud and |
38
|
|
|
indices of the points which are neighbors. |
39
|
|
|
|
40
|
|
|
Parameters |
41
|
|
|
---------- |
42
|
|
|
motif : :class:`numpy.ndarray` |
43
|
|
|
Cartesian coordinates of the motif, shape (no points, dims). |
44
|
|
|
cell : :class:`numpy.ndarray` |
45
|
|
|
The unit cell as a square array, shape (dims, dims). |
46
|
|
|
x : :class:`numpy.ndarray` |
47
|
|
|
Array of points to query for neighbors. For AMD/PDD invariants |
48
|
|
|
this is the motif, or more commonly an asymmetric unit of it. |
49
|
|
|
k : int |
50
|
|
|
Number of nearest neighbors to find for each point in ``x``. |
51
|
|
|
|
52
|
|
|
Returns |
53
|
|
|
------- |
54
|
|
|
dists : numpy.ndarray |
55
|
|
|
Array shape ``(x.shape[0], k)`` of distances from points in |
56
|
|
|
``x`` to their ``k`` nearest neighbors in the periodic set in |
57
|
|
|
order, e.g. ``dists[m][n]`` is the distance from ``x[m]`` to its |
58
|
|
|
n-th nearest neighbor in the periodic set. |
59
|
|
|
""" |
60
|
|
|
|
61
|
|
|
m, dims = motif.shape |
62
|
|
|
|
63
|
|
|
# Get an initial collection of points and a generator for more |
64
|
|
|
int_lattice, int_lat_generator = _integer_lattice_batches(dims, m, k) |
65
|
|
|
cloud = _lattice_to_cloud(motif, int_lattice @ cell) |
66
|
|
|
|
67
|
|
|
# Squared distances to k nearest neighbors |
68
|
|
|
sqdists = _cdist_sqeulcidean(x, cloud) |
69
|
|
|
motif_diam = np.sqrt(sqdists[:, :m].max()) |
70
|
|
|
sqdists.partition(k - 1) |
71
|
|
|
sqdists = sqdists[:, :k] |
72
|
|
|
sqdists.sort() |
73
|
|
|
|
74
|
|
|
# Generate layers of the lattice until they are too far away to give nearer |
75
|
|
|
# neighbors. If a lattice point l has |l| >= max|p-p'| + max_d, (p in x, |
76
|
|
|
# p' in motif) then |p-(l+p')| >= max_d for all p in x, p' in motif, so all |
77
|
|
|
# points l+p' are too far away to be nearer neighbors and l is discarded. |
78
|
|
|
max_sqd = sqdists[:, -1].max() |
79
|
|
|
bound = (np.sqrt(max_sqd) + motif_diam) ** 2 |
80
|
|
|
|
81
|
|
|
while True: |
82
|
|
|
|
83
|
|
|
# Get next layer of lattice and prune distant lattice points |
84
|
|
|
lattice = next(int_lat_generator) @ cell |
85
|
|
|
# The einsum call is equivalent to (lattice ** 2).sum(axis=-1) |
86
|
|
|
lattice = lattice[np.einsum('ij,ij->i', lattice, lattice) < bound] |
87
|
|
|
if lattice.size == 0: # None are close enough |
88
|
|
|
break |
89
|
|
|
|
90
|
|
|
# Squared distances to new points |
91
|
|
|
cloud = _lattice_to_cloud(motif, lattice) |
92
|
|
|
sqdists_ = _cdist_sqeulcidean(x, cloud) |
93
|
|
|
|
94
|
|
|
where_close = (sqdists_ < max_sqd).any(axis=0) |
95
|
|
|
if not np.any(where_close): # None are close enough |
96
|
|
|
break |
97
|
|
|
|
98
|
|
|
# Squared distances to up to k nearest new points |
99
|
|
|
sqdists_ = sqdists_[:, where_close] |
100
|
|
|
if sqdists_.shape[-1] > k: |
101
|
|
|
sqdists_.partition(k - 1) |
102
|
|
|
sqdists_ = sqdists_[:, :k] |
103
|
|
|
sqdists_.sort() |
104
|
|
|
|
105
|
|
|
# Merge existing and new distances |
106
|
|
|
sqdists = _merge_sorted_arrays(sqdists, sqdists_) |
107
|
|
|
max_sqd = sqdists[:, -1].max() |
108
|
|
|
bound = (np.sqrt(max_sqd) + motif_diam) ** 2 |
109
|
|
|
|
110
|
|
|
return np.sqrt(sqdists) |
111
|
|
|
|
112
|
|
|
|
113
|
|
|
def nearest_neighbors_data( |
114
|
|
|
motif: FloatArray, cell: FloatArray, x: FloatArray, k: int |
115
|
|
|
) -> Tuple[FloatArray, FloatArray, IntArray]: |
116
|
|
|
"""Find the ``k`` nearest neighbors in a periodic set for each |
117
|
|
|
point in ``x``, along with the point cloud generated during the |
118
|
|
|
search and indices of the nearest neighbors in the cloud. |
119
|
|
|
|
120
|
|
|
Note: the points in ``x`` are |
121
|
|
|
|
122
|
|
|
Given a periodic set described by ``motif`` and ``cell``, a query |
123
|
|
|
set of points ``x`` and an integer ``k``, find the ``k`` nearest |
124
|
|
|
neighbors in the periodic set for all points in ``x``. Return |
125
|
|
|
an array of distances to neighbors, the point cloud generated |
126
|
|
|
during the search and the indices of which points in the cloud are |
127
|
|
|
the neighbors of points in ``x``. |
128
|
|
|
|
129
|
|
|
Note: the point ``cloud[i]`` in the periodic set comes from (is |
130
|
|
|
translationally equivalent to) the motif point |
131
|
|
|
``motif[i % len(motif)]``, because points are added to ``cloud`` in |
132
|
|
|
batches of whole unit cells and not rearranged. |
133
|
|
|
|
134
|
|
|
Parameters |
135
|
|
|
---------- |
136
|
|
|
motif : :class:`numpy.ndarray` |
137
|
|
|
Cartesian coordinates of the motif, shape (no points, dims). |
138
|
|
|
cell : :class:`numpy.ndarray` |
139
|
|
|
The unit cell as a square array, shape (dims, dims). |
140
|
|
|
x : :class:`numpy.ndarray` |
141
|
|
|
Array of points to query for neighbors. For AMD/PDD invariants |
142
|
|
|
this is the motif, or more commonly an asymmetric unit of it. |
143
|
|
|
k : int |
144
|
|
|
Number of nearest neighbors to find for each point in ``x``. |
145
|
|
|
|
146
|
|
|
Returns |
147
|
|
|
------- |
148
|
|
|
dists : numpy.ndarray |
149
|
|
|
Array shape ``(x.shape[0], k)`` of distances from points in |
150
|
|
|
``x`` to their ``k`` nearest neighbors in the periodic set in |
151
|
|
|
order, e.g. ``dists[m][n]`` is the distance from ``x[m]`` to its |
152
|
|
|
n-th nearest neighbor in the periodic set. |
153
|
|
|
cloud : numpy.ndarray |
154
|
|
|
Collection of points in the periodic set that were generated |
155
|
|
|
during the search. Arranged such that cloud[i] comes from the |
156
|
|
|
motif point motif[i % len(motif)] by translation. |
157
|
|
|
inds : numpy.ndarray |
158
|
|
|
Array shape ``(x.shape[0], k)`` containing the indices of |
159
|
|
|
nearest neighbors in ``cloud``, e.g. the n-th nearest neighbor |
160
|
|
|
to ``x[m]`` is ``cloud[inds[m][n]]``. |
161
|
|
|
""" |
162
|
|
|
|
163
|
|
|
full_cloud = [] |
164
|
|
|
m, dims = motif.shape |
165
|
|
|
|
166
|
|
|
# Get an initial collection of lattice points + a generator for more |
167
|
|
|
int_lattice, int_lat_generator = _integer_lattice_batches(dims, m, k) |
168
|
|
|
cloud = _lattice_to_cloud(motif, int_lattice @ cell) |
169
|
|
|
full_cloud.append(cloud) |
170
|
|
|
cloud_ind_offset = len(cloud) |
171
|
|
|
|
172
|
|
|
# Squared distances to k nearest neighbors + inds of neighbors in cloud |
173
|
|
|
sqdists = _cdist_sqeulcidean(x, cloud) |
174
|
|
|
motif_diam = np.sqrt(sqdists[:, :m].max()) |
175
|
|
|
part_inds = np.argpartition(sqdists, k - 1)[:, :k] |
176
|
|
|
part_sqdists = np.take_along_axis(sqdists, part_inds, axis=-1)[:, :k] |
177
|
|
|
part_sort_inds = np.argsort(part_sqdists) |
178
|
|
|
inds = np.take_along_axis(part_inds, part_sort_inds, axis=-1) |
179
|
|
|
sqdists = np.take_along_axis(part_sqdists, part_sort_inds, axis=-1) |
180
|
|
|
|
181
|
|
|
# Generate layers of lattice until they are too far away to give nearer |
182
|
|
|
# neighbors. For a lattice point l, points in l + motif are further away |
183
|
|
|
# from x than |l| - max|p-p'| (p in x, p' in motif), giving a bound we can |
184
|
|
|
# use to rule out distant lattice points. |
185
|
|
|
max_sqd = sqdists[:, -1].max() |
186
|
|
|
bound = (np.sqrt(max_sqd) + motif_diam) ** 2 |
187
|
|
|
|
188
|
|
|
while True: |
189
|
|
|
|
190
|
|
|
# Get next layer of lattice and prune distant lattice points |
191
|
|
|
lattice = next(int_lat_generator) @ cell |
192
|
|
|
lattice = lattice[np.einsum('ij,ij->i', lattice, lattice) < bound] |
193
|
|
|
if lattice.size == 0: # None are close enough |
194
|
|
|
break |
195
|
|
|
|
196
|
|
|
cloud = _lattice_to_cloud(motif, lattice) |
197
|
|
|
full_cloud.append(cloud) |
198
|
|
|
# Squared distances to new points |
199
|
|
|
sqdists_ = _cdist_sqeulcidean(x, cloud) |
200
|
|
|
close = sqdists_ < max_sqd |
201
|
|
|
if not np.any(close): # None are close enough |
202
|
|
|
break |
203
|
|
|
|
204
|
|
|
# Squared distances to up to k nearest new points + inds |
205
|
|
|
argpart_upto = min(k - 1, sqdists_.shape[-1] - 1) |
206
|
|
|
part_inds = np.argpartition(sqdists_, argpart_upto)[:, :k] |
207
|
|
|
part_sqdists = np.take_along_axis(sqdists_, part_inds, axis=-1)[:, :k] |
208
|
|
|
part_sort_inds = np.argsort(part_sqdists) |
209
|
|
|
inds_ = np.take_along_axis(part_inds, part_sort_inds, axis=-1) |
210
|
|
|
sqdists_ = np.take_along_axis(part_sqdists, part_sort_inds, axis=-1) |
211
|
|
|
|
212
|
|
|
# Offset inds_ so they point to full_cloud instead of cloud |
213
|
|
|
inds_ += cloud_ind_offset |
214
|
|
|
cloud_ind_offset += len(cloud) |
215
|
|
|
|
216
|
|
|
# Merge sqdists & sqdists_, inds & inds_ |
217
|
|
|
sqdists, inds = _merge_sorted_arrays_inds(sqdists, sqdists_, inds, inds_) |
218
|
|
|
max_sqd = sqdists[:, -1].max() |
219
|
|
|
bound = (np.sqrt(max_sqd) + motif_diam) ** 2 |
220
|
|
|
|
221
|
|
|
return np.sqrt(sqdists), np.concatenate(full_cloud), inds |
222
|
|
|
|
223
|
|
|
|
224
|
|
|
def _integer_lattice_batches_cache(func: Callable) -> Callable: |
225
|
|
|
"""Specialised cache for ``_integer_lattice_batches``. Stores layers |
226
|
|
|
of integer lattice points from ``_integer_lattice_generator``. |
227
|
|
|
How many layers are needed depends on the ratio of k to m, so this |
228
|
|
|
is passed to the cache. One cache is kept for each dimension. |
229
|
|
|
""" |
230
|
|
|
|
231
|
|
|
cache = {} # (dims, n_layers): (layers, generator) |
232
|
|
|
npoints_cache = {} # dims: [num points accumulated in layers 1,2,...] |
233
|
|
|
|
234
|
|
|
@functools.wraps(func) |
235
|
|
|
def wrapper( |
236
|
|
|
dims: int, m: int, k: int |
237
|
|
|
) -> Tuple[List[FloatArray], Iterator[FloatArray]]: |
238
|
|
|
|
239
|
|
|
if dims not in npoints_cache: |
|
|
|
|
240
|
|
|
npoints_cache[dims] = [] |
241
|
|
|
|
242
|
|
|
# Use npoints_cache to get how many layers of the lattice are needed |
243
|
|
|
# Delicate code; change in accordance with _integer_lattice_batches |
244
|
|
|
n_layers = 1 |
245
|
|
|
for n_points in npoints_cache[dims]: |
246
|
|
|
if n_points * m > k: |
247
|
|
|
break |
248
|
|
|
n_layers += 1 |
249
|
|
|
n_layers += 1 |
250
|
|
|
|
251
|
|
|
# Update cache and possibly npoints_cache |
252
|
|
|
if (dims, n_layers) not in cache: |
|
|
|
|
253
|
|
|
layers, generator = func(dims, m, k) |
254
|
|
|
n_layers = len(layers) |
255
|
|
|
# If more layers were made than seen so far, update npoints_cache |
256
|
|
|
if len(npoints_cache[dims]) < n_layers: |
257
|
|
|
npoints_cache[dims] = list(accumulate(len(i) for i in layers)) |
258
|
|
|
cache[(dims, n_layers)] = [np.concatenate(layers), generator] |
259
|
|
|
|
260
|
|
|
layers, generator = cache[(dims, n_layers)] |
261
|
|
|
cache[(dims, n_layers)][1], ret_generator = tee(generator) |
262
|
|
|
return layers, ret_generator |
263
|
|
|
|
264
|
|
|
return wrapper |
265
|
|
|
|
266
|
|
|
|
267
|
|
|
@_integer_lattice_batches_cache |
268
|
|
|
def _integer_lattice_batches( |
269
|
|
|
dims: int, m: int, k: int |
270
|
|
|
) -> Tuple[List[FloatArray], Iterator[FloatArray]]: |
271
|
|
|
"""Return an initial batch of integer lattice points (number |
272
|
|
|
according to k & m) and a generator for more distant points. |
273
|
|
|
|
274
|
|
|
Parameters |
275
|
|
|
---------- |
276
|
|
|
dims : int |
277
|
|
|
The dimension of Euclidean space the lattice is in. |
278
|
|
|
m : int |
279
|
|
|
Number of motif points. Used to determine how many layers of |
280
|
|
|
lattice points to put into the initial batch. |
281
|
|
|
k : int |
282
|
|
|
Parameter of ``nearest_neighbors``, the number of neighbors to |
283
|
|
|
find for each point. Used to determine how many layers of |
284
|
|
|
lattice points to put into the initial batch. |
285
|
|
|
|
286
|
|
|
Returns |
287
|
|
|
------- |
288
|
|
|
initial_integer_lattice : :class:`numpy.ndarray` |
289
|
|
|
A collection of integer lattice points. Consists of the first |
290
|
|
|
few layers generated by ``integer_lattice_generator`` (number of |
291
|
|
|
layers depends on m, k). |
292
|
|
|
integer_lattice_generator |
293
|
|
|
A generator for integer lattice points more distant than those |
294
|
|
|
in ``initial_integer_lattice``. |
295
|
|
|
""" |
296
|
|
|
|
297
|
|
|
int_lattice_generator = iter(_integer_lattice_generator(dims)) |
298
|
|
|
layers = [next(int_lattice_generator)] |
299
|
|
|
n_points = 1 |
300
|
|
|
while n_points * m <= k: |
301
|
|
|
layer = next(int_lattice_generator) |
302
|
|
|
n_points += layer.shape[0] |
303
|
|
|
layers.append(layer) |
304
|
|
|
layers.append(next(int_lattice_generator)) |
305
|
|
|
return layers, int_lattice_generator |
306
|
|
|
|
307
|
|
|
|
308
|
|
|
def memoized_generator(generator_function: Callable) -> Callable: |
309
|
|
|
"""Caches results of a generator.""" |
310
|
|
|
cache = {} |
311
|
|
|
@functools.wraps(generator_function) |
312
|
|
|
def wrapper(*args) -> Iterable: |
313
|
|
|
if args not in cache: |
|
|
|
|
314
|
|
|
cache[args] = generator_function(*args) |
315
|
|
|
cache[args], r = tee(cache[args]) |
316
|
|
|
return r |
317
|
|
|
return wrapper |
318
|
|
|
|
319
|
|
|
|
320
|
|
|
@memoized_generator |
321
|
|
|
def _integer_lattice_generator(dims: int) -> Iterable[FloatArray]: |
322
|
|
|
"""Generate batches of integer lattice points. Each yield gives all |
323
|
|
|
points (that have not already been yielded) inside a sphere centered |
324
|
|
|
at the origin with radius d; d starts at 0 and increments by 1 on |
325
|
|
|
each loop. |
326
|
|
|
|
327
|
|
|
Parameters |
328
|
|
|
---------- |
329
|
|
|
dims : int |
330
|
|
|
The dimension of Euclidean space the lattice is in. |
331
|
|
|
|
332
|
|
|
Yields |
333
|
|
|
------- |
334
|
|
|
:class:`numpy.ndarray` |
335
|
|
|
Yields arrays of integer lattice points in `dims`-dimensional |
336
|
|
|
Euclidean space. |
337
|
|
|
""" |
338
|
|
|
|
339
|
|
|
d = 0 |
340
|
|
|
if dims == 1: |
341
|
|
|
yield np.zeros((1, 1), dtype=np.float64) |
342
|
|
|
while True: |
343
|
|
|
d += 1 |
344
|
|
|
yield np.array([[-d], [d]], dtype=np.float64) |
345
|
|
|
|
346
|
|
|
ymax = {} |
347
|
|
|
while True: |
348
|
|
|
positive_int_lattice = [] |
349
|
|
|
while True: |
350
|
|
|
batch = False |
351
|
|
|
for xy in product(range(d + 1), repeat=dims-1): |
352
|
|
|
if xy not in ymax: |
353
|
|
|
ymax[xy] = 0 |
354
|
|
|
if sum(i**2 for i in xy) + ymax[xy]**2 <= d**2: |
355
|
|
|
positive_int_lattice.append((*xy, ymax[xy])) |
356
|
|
|
batch = True |
357
|
|
|
ymax[xy] += 1 |
358
|
|
|
if not batch: |
359
|
|
|
break |
360
|
|
|
pos_int_lat = np.array(positive_int_lattice, dtype=np.float64) |
361
|
|
|
yield _reflect_positive_integer_lattice(pos_int_lat) |
362
|
|
|
d += 1 |
363
|
|
|
|
364
|
|
|
|
365
|
|
|
@numba.njit(cache=True, fastmath=True) |
366
|
|
|
def _reflect_positive_integer_lattice( |
367
|
|
|
positive_int_lattice: FloatArray |
368
|
|
|
) -> FloatArray: |
369
|
|
|
"""Reflect points in the positive quadrant across all combinations |
370
|
|
|
of axes without duplicating points that are invariant under |
371
|
|
|
reflections. |
372
|
|
|
""" |
373
|
|
|
|
374
|
|
|
dims = positive_int_lattice.shape[-1] |
375
|
|
|
batches = [] |
376
|
|
|
batches.extend(positive_int_lattice) |
377
|
|
|
|
378
|
|
|
for n_reflections in range(1, dims + 1): |
379
|
|
|
|
380
|
|
|
axes = np.arange(n_reflections, dtype=np.int64) |
381
|
|
|
off_axes = (positive_int_lattice[:, axes] == 0).sum(axis=-1) == 0 |
382
|
|
|
int_lattice = positive_int_lattice[off_axes] |
383
|
|
|
int_lattice[:, axes] *= -1 |
384
|
|
|
batches.extend(int_lattice) |
385
|
|
|
|
386
|
|
|
while True: |
387
|
|
|
i = n_reflections - 1 |
388
|
|
|
for _ in range(n_reflections): |
389
|
|
|
if axes[i] != i + dims - n_reflections: |
390
|
|
|
break |
391
|
|
|
i -= 1 |
392
|
|
|
else: |
393
|
|
|
break |
394
|
|
|
axes[i] += 1 |
395
|
|
|
for j in range(i + 1, n_reflections): |
396
|
|
|
axes[j] = axes[j-1] + 1 |
397
|
|
|
|
398
|
|
|
off_axes = (positive_int_lattice[:, axes] == 0).sum(axis=-1) == 0 |
399
|
|
|
int_lattice = positive_int_lattice[off_axes] |
400
|
|
|
int_lattice[:, axes] *= -1 |
401
|
|
|
batches.extend(int_lattice) |
402
|
|
|
|
403
|
|
|
int_lattice = np.empty(shape=(len(batches), dims), dtype=np.float64) |
404
|
|
|
for i in range(len(batches)): |
405
|
|
|
int_lattice[i] = batches[i] |
406
|
|
|
return int_lattice |
407
|
|
|
|
408
|
|
|
|
409
|
|
|
@numba.njit(cache=True, fastmath=True) |
410
|
|
|
def _lattice_to_cloud(motif: FloatArray, lattice: FloatArray) -> FloatArray: |
411
|
|
|
"""Create a cloud of points from a periodic set with ``motif``, |
412
|
|
|
and a collection of lattice points ``lattice``. |
413
|
|
|
""" |
414
|
|
|
m = len(motif) |
415
|
|
|
layer = np.empty((m * len(lattice), motif.shape[-1]), dtype=np.float64) |
416
|
|
|
i1 = 0 |
417
|
|
|
for translation in lattice: |
418
|
|
|
i2 = i1 + m |
419
|
|
|
layer[i1:i2] = motif + translation |
420
|
|
|
i1 = i2 |
421
|
|
|
return layer |
422
|
|
|
|
423
|
|
|
|
424
|
|
|
@numba.njit(fastmath=True, cache=True) |
425
|
|
|
def _cdist_sqeulcidean(XA, XB): |
426
|
|
|
mA, mB = XA.shape[0], XB.shape[0] |
427
|
|
|
dims = XA.shape[-1] |
428
|
|
|
dm = np.empty((mA, mB), np.float64) |
429
|
|
|
for i in range(mA): |
430
|
|
|
for j in range(mB): |
431
|
|
|
v = 0 |
432
|
|
|
for d in range(dims): |
433
|
|
|
v += (XA[i, d] - XB[j, d]) ** 2 |
434
|
|
|
dm[i, j] = v |
435
|
|
|
return dm |
436
|
|
|
|
437
|
|
|
|
438
|
|
|
@numba.njit(cache=True, fastmath=True) |
439
|
|
|
def _merge_sorted_arrays(dists: FloatArray, dists_: FloatArray) -> FloatArray: |
440
|
|
|
"""Merge two 2D arrays with sorted rows into one sorted array with |
441
|
|
|
same number of columns as ``dists``. |
442
|
|
|
""" |
443
|
|
|
|
444
|
|
|
m, n_new_points = dists_.shape |
445
|
|
|
ret = np.copy(dists) |
446
|
|
|
|
447
|
|
|
for i in range(m): |
448
|
|
|
# Traverse dists[i] backwards until value smaller than dists_[i, 0] |
449
|
|
|
j = -1 |
450
|
|
|
dp_ = 0 |
451
|
|
|
dist_ = dists_[i, dp_] |
452
|
|
|
while True: |
453
|
|
|
if dists[i, j] <= dist_: |
454
|
|
|
j += 1 |
455
|
|
|
break |
456
|
|
|
j -= 1 |
457
|
|
|
|
458
|
|
|
if j == 0: # If dists_[i, 0] >= dists[i, -1], no need to insert |
459
|
|
|
continue |
460
|
|
|
|
461
|
|
|
# dp points to dists[i], dp_ points to dists_[i], j points to ret[i] |
462
|
|
|
# Fill ret with the larger value, increment pointers and repeat |
463
|
|
|
dp = j |
464
|
|
|
dist = dists[i, dp] |
465
|
|
|
|
466
|
|
|
while j < 0: |
467
|
|
|
if dist <= dist_: |
468
|
|
|
ret[i, j] = dist |
469
|
|
|
dp += 1 |
470
|
|
|
dist = dists[i, dp] |
471
|
|
|
else: |
472
|
|
|
ret[i, j] = dist_ |
473
|
|
|
dp_ += 1 |
474
|
|
|
if dp_ < n_new_points: |
475
|
|
|
dist_ = dists_[i, dp_] |
476
|
|
|
else: # ran out of points in dists_ |
477
|
|
|
dist_ = np.inf |
478
|
|
|
j += 1 |
479
|
|
|
|
480
|
|
|
return ret |
481
|
|
|
|
482
|
|
|
|
483
|
|
|
@numba.njit(cache=True, fastmath=True) |
484
|
|
|
def _merge_sorted_arrays_inds( |
485
|
|
|
dists: FloatArray, |
486
|
|
|
dists_: FloatArray, |
487
|
|
|
inds: IntArray, |
488
|
|
|
inds_: IntArray |
489
|
|
|
) -> Tuple[FloatArray, IntArray]: |
490
|
|
|
"""Similar to _merge_sorted_arrays, but also merges two arrays |
491
|
|
|
``inds`` and ``inds_`` with the same pattern for |
492
|
|
|
``nearest_neighbors_data``. |
493
|
|
|
""" |
494
|
|
|
|
495
|
|
|
m, n_new_points = dists_.shape |
496
|
|
|
ret_dists = np.copy(dists) |
497
|
|
|
ret_inds = np.copy(inds) |
498
|
|
|
|
499
|
|
|
for i in range(m): |
500
|
|
|
|
501
|
|
|
j = -1 |
502
|
|
|
dp_ = 0 |
503
|
|
|
dist_ = dists_[i, dp_] |
504
|
|
|
p_ = inds_[i, dp_] |
505
|
|
|
|
506
|
|
|
while True: |
507
|
|
|
if dists[i, j] <= dist_: |
508
|
|
|
j += 1 |
509
|
|
|
break |
510
|
|
|
j -= 1 |
511
|
|
|
|
512
|
|
|
if j == 0: |
513
|
|
|
continue |
514
|
|
|
|
515
|
|
|
dp = j |
516
|
|
|
dist = dists[i, dp] |
517
|
|
|
p = inds[i, dp] |
518
|
|
|
|
519
|
|
|
while j < 0: |
520
|
|
|
|
521
|
|
|
if dist <= dist_: |
522
|
|
|
ret_dists[i, j] = dist |
523
|
|
|
ret_inds[i, j] = p |
524
|
|
|
dp += 1 |
525
|
|
|
dist = dists[i, dp] |
526
|
|
|
p = inds[i, dp] |
527
|
|
|
else: |
528
|
|
|
ret_dists[i, j] = dist_ |
529
|
|
|
ret_inds[i, j] = p_ |
530
|
|
|
dp_ += 1 |
531
|
|
|
if dp_ < n_new_points: |
532
|
|
|
dist_ = dists_[i, dp_] |
533
|
|
|
p_ = inds_[i, dp_] |
534
|
|
|
else: |
535
|
|
|
dist_ = np.inf |
536
|
|
|
j += 1 |
537
|
|
|
|
538
|
|
|
return ret_dists, ret_inds |
539
|
|
|
|
540
|
|
|
|
541
|
|
|
def nearest_neighbors_minval( |
542
|
|
|
motif: FloatArray, cell: FloatArray, min_val: float |
543
|
|
|
) -> Tuple[FloatArray, FloatArray, IntArray]: |
544
|
|
|
"""Return the same ``dists``/PDD matrix as ``nearest_neighbors``, |
545
|
|
|
but with enough columns such that all values in the last column are |
546
|
|
|
at least ``min_val``. Unlike ``nearest_neighbors``, does not take a |
547
|
|
|
query array ``x`` but only finds neighbors to motif points, and |
548
|
|
|
does not return the point cloud or indices of the nearest |
549
|
|
|
neighbors. Used in ``PDD_reconstructable``. |
550
|
|
|
|
551
|
|
|
TODO: this function should be updated in line with |
552
|
|
|
nearest_neighbors. |
553
|
|
|
""" |
554
|
|
|
|
555
|
|
|
# Generate initial cloud of points from the periodic set |
556
|
|
|
int_lat_generator = _integer_lattice_generator(cell.shape[0]) |
557
|
|
|
int_lat_generator = iter(int_lat_generator) |
558
|
|
|
cloud = [] |
559
|
|
|
for _ in range(3): |
560
|
|
|
cloud.append(_lattice_to_cloud(motif, next(int_lat_generator) @ cell)) |
561
|
|
|
cloud = np.concatenate(cloud) |
562
|
|
|
|
563
|
|
|
# Find k neighbors in the point cloud for points in motif |
564
|
|
|
dists_, inds = KDTree( |
565
|
|
|
cloud, leafsize=30, compact_nodes=False, balanced_tree=False |
566
|
|
|
).query(motif, k=cloud.shape[0]) |
567
|
|
|
dists = np.zeros_like(dists_, dtype=np.float64) |
568
|
|
|
|
569
|
|
|
# Add layers & find k nearest neighbors until all distances smaller than |
570
|
|
|
# min_val don't change |
571
|
|
|
max_cdist = np.amax(cdist(motif, motif)) |
572
|
|
|
while True: |
573
|
|
|
|
574
|
|
|
if np.all(dists_[:, -1] >= min_val): |
575
|
|
|
col = np.argwhere(np.all(dists_ >= min_val, axis=0))[0][0] + 1 |
576
|
|
|
if np.array_equal(dists[:, :col], dists_[:, :col]): |
577
|
|
|
break |
578
|
|
|
|
579
|
|
|
dists = dists_ |
580
|
|
|
lattice = next(int_lat_generator) @ cell |
581
|
|
|
closest_dist_bound = np.linalg.norm(lattice, axis=-1) - max_cdist |
582
|
|
|
is_close = closest_dist_bound <= np.amax(dists_[:, -1]) |
583
|
|
|
if not np.any(is_close): |
584
|
|
|
break |
585
|
|
|
|
586
|
|
|
cloud = np.vstack((cloud, _lattice_to_cloud(motif, lattice[is_close]))) |
587
|
|
|
dists_, inds = KDTree( |
588
|
|
|
cloud, leafsize=30, compact_nodes=False, balanced_tree=False |
589
|
|
|
).query(motif, k=cloud.shape[0]) |
590
|
|
|
|
591
|
|
|
k = np.argwhere(np.all(dists >= min_val, axis=0))[0][0] |
592
|
|
|
return dists_[:, 1:k+1], cloud, inds |
593
|
|
|
|
594
|
|
|
|
595
|
|
|
def generate_concentric_cloud( |
596
|
|
|
motif: FloatArray, cell: FloatArray |
597
|
|
|
) -> Iterable[FloatArray]: |
598
|
|
|
"""Generates batches of points from a periodic set given by (motif, |
599
|
|
|
cell) which get successively further away from the origin. |
600
|
|
|
|
601
|
|
|
Each yield gives all points (that have not already been yielded) |
602
|
|
|
which lie in a unit cell whose corner lattice point was generated by |
603
|
|
|
``generate_integer_lattice(motif.shape[1])``. |
604
|
|
|
|
605
|
|
|
Parameters |
606
|
|
|
---------- |
607
|
|
|
motif : :class:`numpy.ndarray` |
608
|
|
|
Cartesian representation of the motif, shape (no points, dims). |
609
|
|
|
cell : :class:`numpy.ndarray` |
610
|
|
|
Cartesian representation of the unit cell, shape (dims, dims). |
611
|
|
|
|
612
|
|
|
Yields |
613
|
|
|
------- |
614
|
|
|
:class:`numpy.ndarray` |
615
|
|
|
Yields arrays of points from the periodic set. |
616
|
|
|
""" |
617
|
|
|
|
618
|
|
|
int_lat_generator = _integer_lattice_generator(cell.shape[0]) |
619
|
|
|
for layer in int_lat_generator: |
620
|
|
|
yield _lattice_to_cloud(motif, layer @ cell) |
621
|
|
|
|