-
Notifications
You must be signed in to change notification settings - Fork 637
/
distances.py
1616 lines (1399 loc) · 66.7 KB
/
distances.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# -*- Mode: python; tab-width: 4; indent-tabs-mode:nil; -*-
# vim: tabstop=4 expandtab shiftwidth=4 softtabstop=4
#
# MDAnalysis --- https://www.mdanalysis.org
# Copyright (c) 2006-2017 The MDAnalysis Development Team and contributors
# (see the file AUTHORS for the full list of names)
#
# Released under the GNU Public Licence, v2 or any higher version
#
# Please cite your use of MDAnalysis in published work:
#
# R. J. Gowers, M. Linke, J. Barnoud, T. J. E. Reddy, M. N. Melo, S. L. Seyler,
# D. L. Dotson, J. Domanski, S. Buchoux, I. M. Kenney, and O. Beckstein.
# MDAnalysis: A Python package for the rapid analysis of molecular dynamics
# simulations. In S. Benthall and S. Rostrup editors, Proceedings of the 15th
# Python in Science Conference, pages 102-109, Austin, TX, 2016. SciPy.
# doi: 10.25080/majora-629e541a-00e
#
# N. Michaud-Agrawal, E. J. Denning, T. B. Woolf, and O. Beckstein.
# MDAnalysis: A Toolkit for the Analysis of Molecular Dynamics Simulations.
# J. Comput. Chem. 32 (2011), 2319--2327, doi:10.1002/jcc.21787
#
#
"""Fast distance array computation --- :mod:`MDAnalysis.lib.distances`
===================================================================
Fast C-routines to calculate arrays of distances or angles from coordinate
arrays. Many of the functions also exist in parallel versions, which typically
provide higher performance than the serial code.
The boolean attribute `MDAnalysis.lib.distances.USED_OPENMP` can be checked to
see if OpenMP was used in the compilation of MDAnalysis.
Selection of acceleration ("backend")
-------------------------------------
All functions take the optional keyword `backend`, which determines the type of
acceleration. Currently, the following choices are implemented (`backend` is
case-insensitive):
.. Table:: Available *backends* for accelerated distance functions.
========== ========================= ======================================
*backend* module description
========== ========================= ======================================
"serial" :mod:`c_distances` serial implementation in C/Cython
"OpenMP" :mod:`c_distances_openmp` parallel implementation in C/Cython
with OpenMP
========== ========================= ======================================
.. versionadded:: 0.13.0
Functions
---------
.. autofunction:: distance_array
.. autofunction:: self_distance_array
.. autofunction:: capped_distance
.. autofunction:: self_capped_distance
.. autofunction:: calc_bonds
.. autofunction:: calc_angles
.. autofunction:: calc_dihedrals
.. autofunction:: apply_PBC
.. autofunction:: apply_compact_PBC
.. autofunction:: transform_RtoS
.. autofunction:: transform_StoR
.. autofunction:: augment_coordinates(coordinates, box, r)
.. autofunction:: undo_augment(results, translation, nreal)
"""
import numpy as np
from numpy.lib.utils import deprecate
import itertools
from .util import check_coords, check_box
from .mdamath import triclinic_vectors
from ._augment import augment_coordinates, undo_augment
from .nsgrid import FastNS
# hack to select backend with backend=<backend> kwarg. Note that
# the cython parallel code (prange) in parallel.distances is
# independent from the OpenMP code
import importlib
_distances = {}
_distances['serial'] = importlib.import_module(".c_distances",
package="MDAnalysis.lib")
try:
_distances['openmp'] = importlib.import_module(".c_distances_openmp",
package="MDAnalysis.lib")
except ImportError:
pass
del importlib
def _run(funcname, args=None, kwargs=None, backend="serial"):
"""Helper function to select a backend function `funcname`."""
args = args if args is not None else tuple()
kwargs = kwargs if kwargs is not None else dict()
backend = backend.lower()
try:
func = getattr(_distances[backend], funcname)
except KeyError:
errmsg = (f"Function {funcname} not available with backend {backend} "
f"try one of: {_distances.keys()}")
raise ValueError(errmsg) from None
return func(*args, **kwargs)
# serial versions are always available (and are typically used within
# the core and topology modules)
from .c_distances import (calc_distance_array,
calc_distance_array_ortho,
calc_distance_array_triclinic,
calc_self_distance_array,
calc_self_distance_array_ortho,
calc_self_distance_array_triclinic,
coord_transform,
calc_bond_distance,
calc_bond_distance_ortho,
calc_bond_distance_triclinic,
calc_angle,
calc_angle_ortho,
calc_angle_triclinic,
calc_dihedral,
calc_dihedral_ortho,
calc_dihedral_triclinic,
ortho_pbc,
triclinic_pbc)
from .c_distances_openmp import OPENMP_ENABLED as USED_OPENMP
def _check_result_array(result, shape):
"""Check if the result array is ok to use.
The `result` array must meet the following requirements:
* Must have a shape equal to `shape`.
* Its dtype must be ``numpy.float64``.
Paramaters
----------
result : numpy.ndarray or None
The result array to check. If `result` is `None``, a newly created
array of correct shape and dtype ``numpy.float64`` will be returned.
shape : tuple
The shape expected for the `result` array.
Returns
-------
result : numpy.ndarray (``dtype=numpy.float64``, ``shape=shape``)
The input array or a newly created array if the input was ``None``.
Raises
------
ValueError
If `result` is of incorrect shape.
TypeError
If the dtype of `result` is not ``numpy.float64``.
"""
if result is None:
return np.zeros(shape, dtype=np.float64)
if result.shape != shape:
raise ValueError("Result array has incorrect shape, should be {0}, got "
"{1}.".format(shape, result.shape))
if result.dtype != np.float64:
raise TypeError("Result array must be of type numpy.float64, got {}."
"".format(result.dtype))
# The following two lines would break a lot of tests. WHY?!
# if not coords.flags['C_CONTIGUOUS']:
# raise ValueError("{0} is not C-contiguous.".format(desc))
return result
@check_coords('reference', 'configuration', reduce_result_if_single=False,
check_lengths_match=False)
def distance_array(reference, configuration, box=None, result=None,
backend="serial"):
"""Calculate all possible distances between a reference set and another
configuration.
If there are ``n`` positions in `reference` and ``m`` positions in
`configuration`, a distance array of shape ``(n, m)`` will be computed.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
If a 2D numpy array of dtype ``numpy.float64`` with the shape ``(n, m)``
is provided in `result`, then this preallocated array is filled. This can
speed up calculations.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array of shape ``(3,)`` or ``(n, 3)`` (dtype is
arbitrary, will be converted to ``numpy.float32`` internally).
configuration : numpy.ndarray
Configuration coordinate array of shape ``(3,)`` or ``(m, 3)`` (dtype is
arbitrary, will be converted to ``numpy.float32`` internally).
box : array_like, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
result : numpy.ndarray, optional
Preallocated result array which must have the shape ``(n, m)`` and dtype
``numpy.float64``.
Avoids creating the array which saves time when the function
is called repeatedly.
backend : {'serial', 'OpenMP'}, optional
Keyword selecting the type of acceleration.
Returns
-------
d : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n, m)``)
Array containing the distances ``d[i,j]`` between reference coordinates
``i`` and configuration coordinates ``j``.
.. versionchanged:: 0.13.0
Added *backend* keyword.
.. versionchanged:: 0.19.0
Internal dtype conversion of input coordinates to ``numpy.float32``.
Now also accepts single coordinates as input.
"""
confnum = configuration.shape[0]
refnum = reference.shape[0]
distances = _check_result_array(result, (refnum, confnum))
if len(distances) == 0:
return distances
if box is not None:
boxtype, box = check_box(box)
if boxtype == 'ortho':
_run("calc_distance_array_ortho",
args=(reference, configuration, box, distances),
backend=backend)
else:
_run("calc_distance_array_triclinic",
args=(reference, configuration, box, distances),
backend=backend)
else:
_run("calc_distance_array",
args=(reference, configuration, distances),
backend=backend)
return distances
@check_coords('reference', reduce_result_if_single=False)
def self_distance_array(reference, box=None, result=None, backend="serial"):
"""Calculate all possible distances within a configuration `reference`.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
If a 1D numpy array of dtype ``numpy.float64`` with the shape
``(n*(n-1)/2,)`` is provided in `result`, then this preallocated array is
filled. This can speed up calculations.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array of shape ``(3,)`` or ``(n, 3)`` (dtype is
arbitrary, will be converted to ``numpy.float32`` internally).
box : array_like, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
result : numpy.ndarray, optional
Preallocated result array which must have the shape ``(n*(n-1)/2,)`` and
dtype ``numpy.float64``. Avoids creating the array which saves time when
the function is called repeatedly.
backend : {'serial', 'OpenMP'}, optional
Keyword selecting the type of acceleration.
Returns
-------
d : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n*(n-1)/2,)``)
Array containing the distances ``dist[i,j]`` between reference
coordinates ``i`` and ``j`` at position ``d[k]``. Loop through ``d``:
.. code-block:: python
for i in range(n):
for j in range(i + 1, n):
k += 1
dist[i, j] = d[k]
.. versionchanged:: 0.13.0
Added *backend* keyword.
.. versionchanged:: 0.19.0
Internal dtype conversion of input coordinates to ``numpy.float32``.
"""
refnum = reference.shape[0]
distnum = refnum * (refnum - 1) // 2
distances = _check_result_array(result, (distnum,))
if len(distances) == 0:
return distances
if box is not None:
boxtype, box = check_box(box)
if boxtype == 'ortho':
_run("calc_self_distance_array_ortho",
args=(reference, box, distances),
backend=backend)
else:
_run("calc_self_distance_array_triclinic",
args=(reference, box, distances),
backend=backend)
else:
_run("calc_self_distance_array",
args=(reference, distances),
backend=backend)
return distances
def capped_distance(reference, configuration, max_cutoff, min_cutoff=None,
box=None, method=None, return_distances=True):
"""Calculates pairs of indices corresponding to entries in the `reference`
and `configuration` arrays which are separated by a distance lying within
the specified cutoff(s). Optionally, these distances can be returned as
well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
An automatic guessing of the optimal method to calculate the distances is
included in the function. An optional keyword for the method is also
provided. Users can enforce a particular method with this functionality.
Currently brute force, grid search, and periodic KDtree methods are
implemented.
Parameters
-----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)``.
configuration : numpy.ndarray
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)``.
max_cutoff : float
Maximum cutoff distance between the reference and configuration.
min_cutoff : float, optional
Minimum cutoff distance between reference and configuration.
box : array_like, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
method : {'bruteforce', 'nsgrid', 'pkdtree'}, optional
Keyword to override the automatic guessing of the employed search
method.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` and
`configuration` arrays such that the distance between them lies within
the interval (`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th coordinate in `reference` and the ``j``-th coordinate in
`configuration`.
distances : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``), optional
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between the
coordinates ``reference[pairs[k, 0]]`` and
``configuration[pairs[k, 1]]``.
.. code-block:: python
pairs, distances = capped_distances(reference, configuration,
max_cutoff, return_distances=True)
for k, [i, j] in enumerate(pairs):
coord1 = reference[i]
coord2 = configuration[j]
distance = distances[k]
Note
-----
Currently supports brute force, grid-based, and periodic KDtree search
methods.
See Also
--------
distance_array
MDAnalysis.lib.pkdtree.PeriodicKDTree.search
MDAnalysis.lib.nsgrid.FastNS.search
"""
if box is not None:
box = np.asarray(box, dtype=np.float32)
if box.shape[0] != 6:
raise ValueError("Box Argument is of incompatible type. The "
"dimension should be either None or of the form "
"[lx, ly, lz, alpha, beta, gamma]")
method = _determine_method(reference, configuration, max_cutoff,
min_cutoff=min_cutoff, box=box, method=method)
return method(reference, configuration, max_cutoff, min_cutoff=min_cutoff,
box=box, return_distances=return_distances)
def _determine_method(reference, configuration, max_cutoff, min_cutoff=None,
box=None, method=None):
"""Guesses the fastest method for capped distance calculations based on the
size of the coordinate sets and the relative size of the target volume.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)``.
configuration : numpy.ndarray
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)``.
max_cutoff : float
Maximum cutoff distance between `reference` and `configuration`
coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` and `configuration`
coordinates.
box : numpy.ndarray
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
method : {'bruteforce', 'nsgrid', 'pkdtree'}, optional
Keyword to override the automatic guessing of the employed search
method.
Returns
-------
function : callable
The function implementing the guessed (or deliberatly chosen) method.
"""
methods = {'bruteforce': _bruteforce_capped,
'pkdtree': _pkdtree_capped,
'nsgrid': _nsgrid_capped}
if method is not None:
return methods[method.lower()]
if len(reference) < 10 or len(configuration) < 10:
return methods['bruteforce']
elif len(reference) * len(configuration) >= 1e8:
# CAUTION : for large datasets, shouldnt go into 'bruteforce'
# in any case. Arbitrary number, but can be characterized
return methods['nsgrid']
else:
if box is None:
min_dim = np.array([reference.min(axis=0),
configuration.min(axis=0)])
max_dim = np.array([reference.max(axis=0),
configuration.max(axis=0)])
size = max_dim.max(axis=0) - min_dim.min(axis=0)
elif np.all(box[3:] == 90.0):
size = box[:3]
else:
tribox = triclinic_vectors(box)
size = tribox.max(axis=0) - tribox.min(axis=0)
if np.any(max_cutoff > 0.3*size):
return methods['bruteforce']
else:
return methods['nsgrid']
@check_coords('reference', 'configuration', enforce_copy=False,
reduce_result_if_single=False, check_lengths_match=False)
def _bruteforce_capped(reference, configuration, max_cutoff, min_cutoff=None,
box=None, return_distances=True):
"""Capped distance evaluations using a brute force method.
Computes and returns an array containing pairs of indices corresponding to
entries in the `reference` and `configuration` arrays which are separated by
a distance lying within the specified cutoff(s). Employs naive distance
computations (brute force) to find relevant distances.
Optionally, these distances can be returned as well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally).
configuration : array
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)`` (dtype
will be converted to ``numpy.float32`` internally).
max_cutoff : float
Maximum cutoff distance between `reference` and `configuration`
coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` and `configuration`
coordinates.
box : numpy.ndarray, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` and
`configuration` arrays such that the distance between them lies within
the interval (`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th coordinate in `reference` and the ``j``-th coordinate in
`configuration`.
distances : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``), optional
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between the
coordinates ``reference[pairs[k, 0]]`` and
``configuration[pairs[k, 1]]``.
"""
# Default return values (will be overwritten only if pairs are found):
pairs = np.empty((0, 2), dtype=np.intp)
distances = np.empty((0,), dtype=np.float64)
if len(reference) > 0 and len(configuration) > 0:
_distances = distance_array(reference, configuration, box=box)
if min_cutoff is not None:
mask = np.where((_distances <= max_cutoff) & \
(_distances > min_cutoff))
else:
mask = np.where((_distances <= max_cutoff))
if mask[0].size > 0:
pairs = np.c_[mask[0], mask[1]]
if return_distances:
distances = _distances[mask]
if return_distances:
return pairs, distances
else:
return pairs
@check_coords('reference', 'configuration', enforce_copy=False,
reduce_result_if_single=False, check_lengths_match=False)
def _pkdtree_capped(reference, configuration, max_cutoff, min_cutoff=None,
box=None, return_distances=True):
"""Capped distance evaluations using a KDtree method.
Computes and returns an array containing pairs of indices corresponding to
entries in the `reference` and `configuration` arrays which are separated by
a distance lying within the specified cutoff(s). Employs a (periodic) KDtree
algorithm to find relevant distances.
Optionally, these distances can be returned as well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally).
configuration : array
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)`` (dtype
will be converted to ``numpy.float32`` internally).
max_cutoff : float
Maximum cutoff distance between `reference` and `configuration`
coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` and `configuration`
coordinates.
box : numpy.ndarray, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` and
`configuration` arrays such that the distance between them lies within
the interval (`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th coordinate in `reference` and the ``j``-th coordinate in
`configuration`.
distances : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``), optional
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between the
coordinates ``reference[pairs[k, 0]]`` and
``configuration[pairs[k, 1]]``.
"""
from .pkdtree import PeriodicKDTree # must be here to avoid circular import
# Default return values (will be overwritten only if pairs are found):
pairs = np.empty((0, 2), dtype=np.intp)
distances = np.empty((0,), dtype=np.float64)
if len(reference) > 0 and len(configuration) > 0:
kdtree = PeriodicKDTree(box=box)
cut = max_cutoff if box is not None else None
kdtree.set_coords(configuration, cutoff=cut)
_pairs = kdtree.search_tree(reference, max_cutoff)
if _pairs.size > 0:
pairs = _pairs
if (return_distances or (min_cutoff is not None)):
refA, refB = pairs[:, 0], pairs[:, 1]
distances = calc_bonds(reference[refA], configuration[refB],
box=box)
if min_cutoff is not None:
mask = np.where(distances > min_cutoff)
pairs, distances = pairs[mask], distances[mask]
if return_distances:
return pairs, distances
else:
return pairs
@check_coords('reference', 'configuration', enforce_copy=False,
reduce_result_if_single=False, check_lengths_match=False)
def _nsgrid_capped(reference, configuration, max_cutoff, min_cutoff=None,
box=None, return_distances=True):
"""Capped distance evaluations using a grid-based search method.
Computes and returns an array containing pairs of indices corresponding to
entries in the `reference` and `configuration` arrays which are separated by
a distance lying within the specified cutoff(s). Employs a grid-based search
algorithm to find relevant distances.
Optionally, these distances can be returned as well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally).
configuration : array
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)`` (dtype
will be converted to ``numpy.float32`` internally).
max_cutoff : float
Maximum cutoff distance between `reference` and `configuration`
coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` and `configuration`
coordinates.
box : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``), optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` and
`configuration` arrays such that the distance between them lies within
the interval (`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th coordinate in `reference` and the ``j``-th coordinate in
`configuration`.
distances : numpy.ndarray, optional
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between the
coordinates ``reference[pairs[k, 0]]`` and
``configuration[pairs[k, 1]]``.
"""
# Default return values (will be overwritten only if pairs are found):
pairs = np.empty((0, 2), dtype=np.intp)
distances = np.empty((0,), dtype=np.float64)
if len(reference) > 0 and len(configuration) > 0:
if box is None:
# create a pseudobox
# define the max range
# and supply the pseudobox
# along with only one set of coordinates
pseudobox = np.zeros(6, dtype=np.float32)
all_coords = np.concatenate([reference, configuration])
lmax = all_coords.max(axis=0)
lmin = all_coords.min(axis=0)
# Using maximum dimension as the box size
boxsize = (lmax-lmin).max()
# to avoid failures for very close particles but with
# larger cutoff
boxsize = np.maximum(boxsize, 2 * max_cutoff)
pseudobox[:3] = boxsize + 2.2*max_cutoff
pseudobox[3:] = 90.
shiftref, shiftconf = reference.copy(), configuration.copy()
# Extra padding near the origin
shiftref -= lmin - 0.1*max_cutoff
shiftconf -= lmin - 0.1*max_cutoff
gridsearch = FastNS(max_cutoff, shiftconf, box=pseudobox, pbc=False)
results = gridsearch.search(shiftref)
else:
gridsearch = FastNS(max_cutoff, configuration, box=box)
results = gridsearch.search(reference)
pairs = results.get_pairs()
if return_distances or (min_cutoff is not None):
distances = results.get_pair_distances()
if min_cutoff is not None:
idx = distances > min_cutoff
pairs, distances = pairs[idx], distances[idx]
if return_distances:
return pairs, distances
else:
return pairs
def self_capped_distance(reference, max_cutoff, min_cutoff=None, box=None,
method=None, return_distances=True):
"""Calculates pairs of indices corresponding to entries in the `reference`
array which are separated by a distance lying within the specified
cutoff(s). Optionally, these distances can be returned as well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
An automatic guessing of the optimal method to calculate the distances is
included in the function. An optional keyword for the method is also
provided. Users can enforce a particular method with this functionality.
Currently brute force, grid search, and periodic KDtree methods are
implemented.
Parameters
-----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)``.
max_cutoff : float
Maximum cutoff distance between `reference` coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` coordinates.
box : array_like, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
method : {'bruteforce', 'nsgrid', 'pkdtree'}, optional
Keyword to override the automatic guessing of the employed search
method.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` array
such that the distance between them lies within the interval
(`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th and the ``j``-th coordinate in `reference`.
distances : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``)
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between the
coordinates ``reference[pairs[k, 0]]`` and ``reference[pairs[k, 1]]``.
.. code-block:: python
pairs, distances = self_capped_distances(reference, max_cutoff,
return_distances=True)
for k, [i, j] in enumerate(pairs):
coord1 = reference[i]
coord2 = reference[j]
distance = distances[k]
Note
-----
Currently supports brute force, grid-based, and periodic KDtree search
methods.
See Also
--------
self_distance_array
MDAnalysis.lib.pkdtree.PeriodicKDTree.search
MDAnalysis.lib.nsgrid.FastNS.self_search
.. versionchanged:: 0.20.0
Added `return_distances` keyword.
"""
if box is not None:
box = np.asarray(box, dtype=np.float32)
if box.shape[0] != 6:
raise ValueError("Box Argument is of incompatible type. The "
"dimension should be either None or of the form "
"[lx, ly, lz, alpha, beta, gamma]")
method = _determine_method_self(reference, max_cutoff,
min_cutoff=min_cutoff,
box=box, method=method)
return method(reference, max_cutoff, min_cutoff=min_cutoff, box=box,
return_distances=return_distances)
def _determine_method_self(reference, max_cutoff, min_cutoff=None, box=None,
method=None):
"""Guesses the fastest method for capped distance calculations based on the
size of the `reference` coordinate set and the relative size of the target
volume.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)``.
max_cutoff : float
Maximum cutoff distance between `reference` coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` coordinates.
box : numpy.ndarray
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
method : {'bruteforce', 'nsgrid', 'pkdtree'}, optional
Keyword to override the automatic guessing of the employed search
method.
Returns
-------
function : callable
The function implementing the guessed (or deliberatly chosen) method.
"""
methods = {'bruteforce': _bruteforce_capped_self,
'pkdtree': _pkdtree_capped_self,
'nsgrid': _nsgrid_capped_self}
if method is not None:
return methods[method.lower()]
if len(reference) < 100:
return methods['bruteforce']
if box is None:
min_dim = np.array([reference.min(axis=0)])
max_dim = np.array([reference.max(axis=0)])
size = max_dim.max(axis=0) - min_dim.min(axis=0)
elif np.all(box[3:] == 90.0):
size = box[:3]
else:
tribox = triclinic_vectors(box)
size = tribox.max(axis=0) - tribox.min(axis=0)
if max_cutoff < 0.03*size.min():
return methods['pkdtree']
else:
return methods['nsgrid']
@check_coords('reference', enforce_copy=False, reduce_result_if_single=False)
def _bruteforce_capped_self(reference, max_cutoff, min_cutoff=None, box=None,
return_distances=True):
"""Capped distance evaluations using a brute force method.
Computes and returns an array containing pairs of indices corresponding to
entries in the `reference` array which are separated by a distance lying
within the specified cutoff(s). Employs naive distance computations (brute
force) to find relevant distances. Optionally, these distances can be
returned as well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally).
max_cutoff : float
Maximum cutoff distance between `reference` coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` coordinates.
box : numpy.ndarray, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` array
such that the distance between them lies within the interval
(`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th and the ``j``-th coordinate in `reference`.
distances : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``), optional
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between the
coordinates ``reference[pairs[k, 0]]`` and
``configuration[pairs[k, 1]]``.
.. versionchanged:: 0.20.0
Added `return_distances` keyword.
"""
# Default return values (will be overwritten only if pairs are found):
pairs = np.empty((0, 2), dtype=np.intp)
distances = np.empty((0,), dtype=np.float64)
N = len(reference)
# We're searching within a single coordinate set, so we need at least two
# coordinates to find distances between them.
if N > 1:
distvec = self_distance_array(reference, box=box)
dist = np.full((N, N), np.finfo(np.float64).max, dtype=np.float64)
dist[np.triu_indices(N, 1)] = distvec
if min_cutoff is not None:
mask = np.where((dist <= max_cutoff) & (dist > min_cutoff))
else:
mask = np.where((dist <= max_cutoff))
if mask[0].size > 0:
pairs = np.c_[mask[0], mask[1]]
distances = dist[mask]
if return_distances:
return pairs, distances
return pairs
@check_coords('reference', enforce_copy=False, reduce_result_if_single=False)
def _pkdtree_capped_self(reference, max_cutoff, min_cutoff=None, box=None,
return_distances=True):
"""Capped distance evaluations using a KDtree method.
Computes and returns an array containing pairs of indices corresponding to
entries in the `reference` array which are separated by a distance lying
within the specified cutoff(s). Employs a (periodic) KDtree algorithm to
find relevant distances. Optionally, these distances can be returned as
well.
If the optional argument `box` is supplied, the minimum image convention is
applied when calculating distances. Either orthogonal or triclinic boxes are
supported.
Parameters
----------
reference : numpy.ndarray
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally).
max_cutoff : float
Maximum cutoff distance between `reference` coordinates.
min_cutoff : float, optional
Minimum cutoff distance between `reference` coordinates.
box : numpy.ndarray, optional
The unitcell dimensions of the system, which can be orthogonal or
triclinic and must be provided in the same format as returned by
:attr:`MDAnalysis.coordinates.base.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
Returns
-------
pairs : numpy.ndarray (``dtype=numpy.int64``, ``shape=(n_pairs, 2)``)
Pairs of indices, corresponding to coordinates in the `reference` array
such that the distance between them lies within the interval
(`min_cutoff`, `max_cutoff`].
Each row in `pairs` is an index pair ``[i, j]`` corresponding to the
``i``-th and the ``j``-th coordinate in `reference`.
distances : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n_pairs,)``)
Distances corresponding to each pair of indices. Only returned if
`return_distances` is ``True``. ``distances[k]`` corresponds to the
``k``-th pair returned in `pairs` and gives the distance between
the coordinates ``reference[pairs[k, 0]]`` and
``reference[pairs[k, 1]]``.
.. versionchanged:: 0.20.0
Added `return_distances` keyword.
"""
from .pkdtree import PeriodicKDTree # must be here to avoid circular import
# Default return values (will be overwritten only if pairs are found):
pairs = np.empty((0, 2), dtype=np.intp)
distances = np.empty((0,), dtype=np.float64)
# We're searching within a single coordinate set, so we need at least two
# coordinates to find distances between them.
if len(reference) > 1: