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 1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717 1718 1719 1720 1721 1722 1723 1724 1725 1726 1727 1728 1729 1730 1731 1732 1733 1734 1735 1736 1737 1738 1739 1740 1741 1742 1743 1744 1745 1746 1747 1748 1749 1750 1751 1752 1753 1754 1755 1756 1757 1758 1759 1760 1761 1762 1763 1764 1765 1766 1767 1768 1769 1770 1771 1772 1773 1774 1775 1776 1777 1778 1779 1780 1781 1782 1783 1784 1785 1786 1787 1788 1789 1790 1791 1792 1793 1794 1795 1796 1797 1798 1799 1800 1801 1802 1803 1804 1805 1806 1807 1808 1809 1810 1811 1812 1813 1814 1815 1816 1817 1818 1819 1820 1821 1822 1823 1824 1825 1826 1827 1828 1829 1830 1831 1832 1833 1834 1835 1836 1837 1838 1839 1840 1841 1842 1843 1844 1845 1846 1847 1848 1849 1850 1851 1852 1853 1854 1855 1856 1857 1858 1859 1860 1861 1862 1863 1864 1865 1866 1867 1868 1869 1870 1871 1872 1873 1874 1875 1876 1877 1878 1879 1880 1881 1882 1883 1884 1885 1886 1887 1888 1889 1890 1891 1892 1893 1894 1895 1896 1897 1898 1899 1900 1901 1902 1903 1904 1905 1906 1907 1908 1909 1910 1911 1912 1913 1914 1915 1916 1917 1918 1919 1920 1921 1922 1923 1924 1925 1926 1927 1928 1929 1930 1931 1932 1933 1934 1935 1936 1937 1938 1939 1940 1941 1942 1943 1944 1945 1946 1947 1948 1949 1950 1951 1952 1953 1954 1955 1956 1957 1958 1959 1960 1961 1962 1963 1964 1965 1966 1967 1968 1969 1970 1971 1972 1973 1974 1975 1976 1977 1978 1979 1980 1981 1982 1983 1984 1985 1986 1987 1988 1989 1990 1991 1992 1993 1994 1995 1996 1997 1998 1999 2000 2001 2002 2003 2004 2005 2006 2007 2008 2009 2010 2011 2012 2013 2014 2015 2016 2017 2018 2019 2020 2021 2022 2023 2024 2025 2026 2027 2028 2029 2030 2031 2032 2033 2034 2035 2036 2037 2038 2039 2040 2041 2042 2043 2044 2045 2046 2047 2048 2049 2050 2051 2052 2053 2054 2055 2056 2057 2058 2059 2060 2061 2062 2063 2064 2065 2066 2067 2068 2069 2070 2071 2072 2073 2074 2075 2076 2077 2078 2079 2080 2081 2082 2083 2084 2085 2086 2087 2088 2089 2090 2091 2092 2093 2094 2095 2096 2097 2098 2099 2100 2101 2102 2103 2104 2105 2106 2107 2108 2109 2110 2111 2112 2113 2114 2115 2116 2117 2118 2119 2120 2121 2122 2123 2124 2125 2126 2127 2128 2129 2130 2131 2132 2133 2134 2135 2136 2137 2138 2139 2140 2141 2142 2143 2144 2145
|
# -*- 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 Lesser GNU Public Licence, v2.1 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. Distance functions can accept a NumPy :class:`np.ndarray` or an
:class:`~MDAnalysis.core.groups.AtomGroup`. 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
"distopia" `_distopia` SIMD-accelerated implementation
with the `distopia`_ library
========== ========================= ======================================
Use of the distopia library
---------------------------
MDAnalysis has developed a standalone library, `distopia`_ for accelerating
the distance functions in this module using explicit SIMD vectorisation.
This can provide many-fold speedups in calculating distances. Distopia is
under active development and as such only a selection of functions in this
module are covered. Consult the following table to see if the function
you wish to use is covered by distopia. For more information see the
`distopia documentation`_.
.. table:: Functions available using the `distopia`_ backend.
:align: center
+-------------------------------------------------------+
| Functions |
+=======================================================+
| :func:`MDAnalysis.lib.distances.calc_bonds` |
+-------------------------------------------------------+
| :func:`MDAnalysis.lib.distances.calc_angles` |
+-------------------------------------------------------+
| :func:`MDAnalysis.lib.distances.calc_dihedrals` |
+-------------------------------------------------------+
| :func:`MDAnalysis.lib.distances.distance_array` |
+-------------------------------------------------------+
| :func:`MDAnalysis.lib.distances.self_distance_array` |
+-------------------------------------------------------+
If `distopia`_ is installed, the functions in this table will accept the key
'distopia' for the `backend` keyword argument. The variable
:data:`HAS_DISTOPIA` is set to ``True`` if distopia is available.
If the distopia backend is selected the `distopia` library will be used to
calculate the distances. Note that for functions listed in this table
**distopia is not the default backend and must be explicitly selected.**
.. Note::
Due to the use of Instruction Set Architecture (`ISA`_) specific SIMD
intrinsics in distopia via `HWY`_, the precision of your results may
depend on the ISA available on your machine. However, in all tested cases
distopia satisfied the accuracy thresholds used to the functions in this
module. Please document any issues you encounter with distopia's accuracy
in the `relevant distopia issue`_ on the MDAnalysis GitHub repository.
.. _distopia: https://github.com/MDAnalysis/distopia
.. _distopia documentation: https://www.mdanalysis.org/distopia
.. _ISA: https://en.wikipedia.org/wiki/Instruction_set_architecture
.. _HWY: https://github.com/google/highway
.. _relevant distopia issue: https://github.com/MDAnalysis/mdanalysis/issues/3915
.. versionadded:: 0.13.0
.. versionchanged:: 2.3.0
Distance functions can now accept an
:class:`~MDAnalysis.core.groups.AtomGroup` or an :class:`np.ndarray`
.. versionchanged:: 2.5.0
Interface to the `distopia`_ package added.
.. versionchanged:: 2.9.0
Distopia support greatly expanded (with distopia ≥ 0.4.0).
Constants
---------
.. data:: HAS_DISTOPIA
This variable is ``True`` if the :mod:`distopia` package has been
installed and is available as a `backend`. Otherwise it is
``False``.
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:: transform_RtoS
.. autofunction:: transform_StoR
.. autofunction:: augment_coordinates(coordinates, box, r)
.. autofunction:: undo_augment(results, translation, nreal)
.. autofunction:: minimize_vectors(vectors, box)
"""
import numpy as np
import numpy.typing as npt
from typing import Union, Optional, Callable
from typing import TYPE_CHECKING
if TYPE_CHECKING: # pragma: no cover
from ..core.groups import AtomGroup
from .util import check_coords, check_box
from .mdamath import triclinic_vectors
from ._augment import augment_coordinates, undo_augment
from .nsgrid import FastNS
from .c_distances import _minimize_vectors_ortho, _minimize_vectors_triclinic
from ._distopia import HAS_DISTOPIA
# 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
if HAS_DISTOPIA:
_distances["distopia"] = importlib.import_module(
"._distopia", package="MDAnalysis.lib"
)
del importlib
def _run(
funcname: str,
args: Optional[tuple] = None,
kwargs: Optional[dict] = None,
backend: str = "serial",
) -> Callable:
"""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 (
_UINT64_MAX,
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: Optional[npt.NDArray], shape: tuple
) -> npt.NDArray:
"""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,
allow_atomgroup=True,
)
def distance_array(
reference: Union[npt.NDArray, "AtomGroup"],
configuration: Union[npt.NDArray, "AtomGroup"],
box: Optional[npt.NDArray] = None,
result: Optional[npt.NDArray] = None,
backend: str = "serial",
) -> npt.NDArray:
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array of shape ``(3,)`` or ``(n, 3)`` (dtype is
arbitrary, will be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
configuration : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Configuration coordinate array of shape ``(3,)`` or ``(m, 3)`` (dtype is
arbitrary, will be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.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', 'distopia'}, 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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.9.0
Added support for the `distopia` backend.
"""
confnum = configuration.shape[0]
refnum = reference.shape[0]
# check resulting array will not overflow UINT64_MAX
if refnum * confnum > _UINT64_MAX:
raise ValueError(
f"Size of resulting array {refnum * confnum} elements"
" larger than size of maximum integer"
)
distances = _check_result_array(result, (refnum, confnum))
if len(distances) == 0:
return distances
if backend == "distopia":
# distopia requires that all the input arrays are the same type,
# while MDAnalysis allows for mixed types, this should be changed
# pre 3.0.0 release see issue #3707
distances = distances.astype(np.float32)
box = np.asarray(box).astype(np.float32) if box is not None else None
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,
)
if backend == "distopia":
# mda expects the result to be in float64, so we need to convert it back
# to float64, change for 3.0, see #3707
distances = distances.astype(np.float64)
if result is not None:
result[:] = distances
return distances
@check_coords("reference", reduce_result_if_single=False, allow_atomgroup=True)
def self_distance_array(
reference: Union[npt.NDArray, "AtomGroup"],
box: Optional[npt.NDArray] = None,
result: Optional[npt.NDArray] = None,
backend: str = "serial",
) -> npt.NDArray:
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array of shape ``(3,)`` or ``(n, 3)`` (dtype is
arbitrary, will be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.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', 'distopia'}, 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):
dist[i, j] = dist[j, i] = d[k]
k += 1
.. versionchanged:: 0.13.0
Added *backend* keyword.
.. versionchanged:: 0.19.0
Internal dtype conversion of input coordinates to ``numpy.float32``.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.9.0
Added support for the `distopia` backend.
"""
refnum = reference.shape[0]
distnum = refnum * (refnum - 1) // 2
# check resulting array will not overflow UINT64_MAX
if distnum > _UINT64_MAX:
raise ValueError(
f"Size of resulting array {distnum} elements larger"
" than size of maximum integer"
)
distances = _check_result_array(result, (distnum,))
if len(distances) == 0:
return distances
if backend == "distopia":
# distopia requires that all the input arrays are the same type,
# while MDAnalysis allows for mixed types, this should be changed
# pre 3.0.0 release see issue #3707
distances = distances.astype(np.float32)
box = np.asarray(box).astype(np.float32) if box is not None else None
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,
)
if backend == "distopia":
# mda expects the result to be in float64, so we need to convert it back
# to float64, change for 3.0, see #3707
distances = distances.astype(np.float64)
if result is not None:
result[:] = distances
return distances
@check_coords(
"reference",
"configuration",
enforce_copy=False,
reduce_result_if_single=False,
check_lengths_match=False,
allow_atomgroup=True,
)
def capped_distance(
reference: Union[npt.NDArray, "AtomGroup"],
configuration: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
method: Optional[str] = None,
return_distances: Optional[bool] = True,
backend: Optional[str] = "serial",
):
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)``. Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
configuration : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)``. Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.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.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
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]
See Also
--------
distance_array
MDAnalysis.lib.pkdtree.PeriodicKDTree.search
MDAnalysis.lib.nsgrid.FastNS.search
.. versionchanged:: 1.0.1
nsgrid was temporarily removed and replaced with pkdtree due to issues
relating to its reliability and accuracy (Issues #2919, #2229, #2345,
#2670, #2930)
.. versionchanged:: 1.0.2
nsgrid enabled again
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.10.0
Added the "backend" argument to select the type of acceleration of
the distance calculations.
"""
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]"
)
# The check_coords decorator made sure that reference and configuration
# are arrays of positions. Mypy does not know about that so we have to
# tell it.
reference_positions: npt.NDArray = reference # type: ignore
configuration_positions: npt.NDArray = configuration # type: ignore
function = _determine_method(
reference_positions,
configuration_positions,
max_cutoff,
min_cutoff=min_cutoff,
box=box,
method=method,
)
if function.__name__ == "_nsgrid_capped":
return function(
reference,
configuration,
max_cutoff,
min_cutoff=min_cutoff,
box=box,
return_distances=return_distances,
)
else:
return function(
reference,
configuration,
max_cutoff,
min_cutoff=min_cutoff,
box=box,
return_distances=return_distances,
backend=backend,
)
def _determine_method(
reference: npt.NDArray,
configuration: npt.NDArray,
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
method: Optional[str] = None,
) -> Callable:
"""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.timestep.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.
.. versionchanged:: 1.0.1
nsgrid was temporarily removed and replaced with pkdtree due to issues
relating to its reliability and accuracy (Issues #2919, #2229, #2345,
#2670, #2930)
.. versionchanged:: 1.1.0
enabled nsgrid again
"""
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,
allow_atomgroup=True,
)
def _bruteforce_capped(
reference: Union[npt.NDArray, "AtomGroup"],
configuration: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
return_distances: Optional[bool] = True,
backend: Optional[str] = "serial",
):
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
configuration : array or :class:`~MDAnalysis.core.groups.AtomGroup`
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)`` (dtype
will be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
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]]``.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.10.0
Added the "backend" argument to select the type of acceleration of
the distance calculations.
"""
# 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, backend=backend
)
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,
allow_atomgroup=True,
)
def _pkdtree_capped(
reference: Union[npt.NDArray, "AtomGroup"],
configuration: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
return_distances: Optional[bool] = True,
backend: Optional[str] = "serial",
):
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
configuration : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)`` (dtype
will be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
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]]``.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.10.0
Added the "backend" argument to select the type of acceleration of
the distance calculations.
"""
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,
backend=backend,
)
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,
allow_atomgroup=True,
)
def _nsgrid_capped(
reference: Union[npt.NDArray, "AtomGroup"],
configuration: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
return_distances: Optional[bool] = 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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
configuration : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Configuration coordinate array with shape ``(3,)`` or ``(m, 3)`` (dtype
will be converted to ``numpy.float32`` internally). Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.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]]``.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
"""
# 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.0
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
@check_coords(
"reference",
enforce_copy=False,
reduce_result_if_single=False,
check_lengths_match=False,
allow_atomgroup=True,
)
def self_capped_distance(
reference: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
method: Optional[str] = None,
return_distances: Optional[bool] = True,
backend: Optional[str] = "serial",
):
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)``. Also
accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.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.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
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.
.. versionchanged:: 1.0.1
nsgrid was temporarily removed and replaced with pkdtree due to issues
relating to its reliability and accuracy (Issues #2919, #2229, #2345,
#2670, #2930)
.. versionchanged:: 1.0.2
enabled nsgrid again
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.10.0
Added the "backend" argument to select the type of acceleration of
the distance calculations.
"""
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]"
)
# The check_coords decorator made sure that reference is an
# array of positions. Mypy does not know about that so we have to
# tell it.
reference_positions: npt.NDArray = reference # type: ignore
function = _determine_method_self(
reference_positions,
max_cutoff,
min_cutoff=min_cutoff,
box=box,
method=method,
)
if function.__name__ == "_nsgrid_capped_self":
return function(
reference,
max_cutoff,
min_cutoff=min_cutoff,
box=box,
return_distances=return_distances,
)
else:
return function(
reference,
max_cutoff,
min_cutoff=min_cutoff,
box=box,
return_distances=return_distances,
backend=backend,
)
def _determine_method_self(
reference: npt.NDArray,
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
method: Optional[str] = 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.timestep.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.
.. versionchanged:: 1.0.1
nsgrid was temporarily removed and replaced with pkdtree due to issues
relating to its reliability and accuracy (Issues #2919, #2229, #2345,
#2670, #2930)
.. versionchanged:: 1.0.2
enabled nsgrid again
"""
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,
allow_atomgroup=True,
)
def _bruteforce_capped_self(
reference: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
return_distances: Optional[bool] = True,
backend: Optional[str] = "serial",
):
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.10.0
Added the "backend" argument to select the type of acceleration of
the distance calculations.
"""
# 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, backend=backend)
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,
allow_atomgroup=True,
)
def _pkdtree_capped_self(
reference: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
return_distances: Optional[bool] = True,
backend: Optional[str] = "serial",
):
"""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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
return_distances : bool, optional
If set to ``True``, distances will also be returned.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.10.0
Added the "backend" argument to select the type of acceleration of
the distance calculations.
"""
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:
kdtree = PeriodicKDTree(box=box)
cut = max_cutoff if box is not None else None
kdtree.set_coords(reference, cutoff=cut)
_pairs = kdtree.search_pairs(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], reference[refB], box=box, backend=backend
)
if min_cutoff is not None:
idx = distances > min_cutoff
pairs, distances = pairs[idx], distances[idx]
if return_distances:
return pairs, distances
return pairs
@check_coords(
"reference",
enforce_copy=False,
reduce_result_if_single=False,
allow_atomgroup=True,
)
def _nsgrid_capped_self(
reference: Union[npt.NDArray, "AtomGroup"],
max_cutoff: float,
min_cutoff: Optional[float] = None,
box: Optional[npt.NDArray] = None,
return_distances: Optional[bool] = 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` array 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 or :class:`~MDAnalysis.core.groups.AtomGroup`
Reference coordinate array with shape ``(3,)`` or ``(n, 3)`` (dtype will
be converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
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, 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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
"""
# 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:
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)
lmax = reference.max(axis=0)
lmin = reference.min(axis=0)
# Using maximum dimension as the box size
boxsize = (lmax - lmin).max()
# to avoid failures of very close particles
# but with larger cutoff
if boxsize < 2 * max_cutoff:
# just enough box size so that NSGrid doesnot fails
sizefactor = 2.2 * max_cutoff / boxsize
else:
sizefactor = 1.2
pseudobox[:3] = sizefactor * boxsize
pseudobox[3:] = 90.0
shiftref = reference.copy()
# Extra padding near the origin
shiftref -= lmin - 0.1 * boxsize
gridsearch = FastNS(max_cutoff, shiftref, box=pseudobox, pbc=False)
results = gridsearch.self_search()
else:
gridsearch = FastNS(max_cutoff, reference, box=box)
results = gridsearch.self_search()
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
return pairs
@check_coords("coords")
def transform_RtoS(coords, box, backend="serial"):
"""Transform an array of coordinates from real space to S space (a.k.a.
lambda space)
S space represents fractional space within the unit cell for this system.
Reciprocal operation to :meth:`transform_StoR`.
Parameters
----------
coords : numpy.ndarray
A ``(3,)`` or ``(n, 3)`` array of coordinates (dtype is arbitrary, will
be converted to ``numpy.float32`` internally).
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
backend : {'serial', 'OpenMP'}, optional
Keyword selecting the type of acceleration.
Returns
-------
newcoords : numpy.ndarray (``dtype=numpy.float32``, ``shape=coords.shape``)
An array containing fractional coordiantes.
.. versionchanged:: 0.13.0
Added *backend* keyword.
.. versionchanged:: 0.19.0
Internal dtype conversion of input coordinates to ``numpy.float32``.
Now also accepts (and, likewise, returns) a single coordinate.
"""
if len(coords) == 0:
return coords
boxtype, box = check_box(box)
if boxtype == "ortho":
box = np.diag(box)
box = box.astype(np.float64)
# Create inverse matrix of box
# need order C here
inv = np.array(np.linalg.inv(box), order="C")
_run("coord_transform", args=(coords, inv), backend=backend)
return coords
@check_coords("coords")
def transform_StoR(coords, box, backend="serial"):
"""Transform an array of coordinates from S space into real space.
S space represents fractional space within the unit cell for this system.
Reciprocal operation to :meth:`transform_RtoS`
Parameters
----------
coords : numpy.ndarray
A ``(3,)`` or ``(n, 3)`` array of coordinates (dtype is arbitrary, will
be converted to ``numpy.float32`` internally).
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
backend : {'serial', 'OpenMP'}, optional
Keyword selecting the type of acceleration.
Returns
-------
newcoords : numpy.ndarray (``dtype=numpy.float32``, ``shape=coords.shape``)
An array containing real space coordiantes.
.. versionchanged:: 0.13.0
Added *backend* keyword.
.. versionchanged:: 0.19.0
Internal dtype conversion of input coordinates to ``numpy.float32``.
Now also accepts (and, likewise, returns) a single coordinate.
"""
if len(coords) == 0:
return coords
boxtype, box = check_box(box)
if boxtype == "ortho":
box = np.diag(box)
box = box.astype(np.float64)
_run("coord_transform", args=(coords, box), backend=backend)
return coords
@check_coords("coords1", "coords2", allow_atomgroup=True)
def calc_bonds(
coords1: Union[npt.NDArray, "AtomGroup"],
coords2: Union[npt.NDArray, "AtomGroup"],
box: Optional[npt.NDArray] = None,
result: Optional[npt.NDArray] = None,
backend: str = "serial",
) -> npt.NDArray:
"""Calculates the bond lengths between pairs of atom positions from the two
coordinate arrays `coords1` and `coords2`, which must contain the same
number of coordinates. ``coords1[i]`` and ``coords2[i]`` represent the
positions of atoms connected by the ``i``-th bond. If single coordinates are
supplied, a single distance will be returned.
In comparison to :meth:`distance_array` and :meth:`self_distance_array`,
which calculate distances between all possible combinations of coordinates,
:meth:`calc_bonds` only calculates distances between pairs of coordinates,
similar to::
numpy.linalg.norm(a - b) for a, b in zip(coords1, coords2)
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 numpy array of dtype ``numpy.float64`` with shape ``(n,)`` (for ``n``
coordinate pairs) is provided in `result`, then this preallocated array is
filled. This can speed up calculations.
Parameters
----------
coords1 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` for one half of a
single or ``n`` bonds, respectively (dtype is arbitrary, will be
converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
coords2 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` for the other half of
a single or ``n`` bonds, respectively (dtype is arbitrary, will be
converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
result : numpy.ndarray, optional
Preallocated result array of dtype ``numpy.float64`` and shape ``(n,)``
(for ``n`` coordinate pairs). Avoids recreating the array in repeated
function calls.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration. Defaults to 'serial'.
Returns
-------
bondlengths : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n,)``) or
numpy.float64 Array containing the bond lengths between each pair of
coordinates. If two single coordinates were supplied, their distance is
returned as a single number instead of an array.
.. versionadded:: 0.8
.. 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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
.. versionchanged:: 2.5.0
Can now optionally use the fast distance functions from distopia
"""
numatom = coords1.shape[0]
bondlengths = _check_result_array(result, (numatom,))
if backend == "distopia":
# distopia requires that all the input arrays are the same type,
# while MDAnalysis allows for mixed types, this should be changed
# pre 3.0.0 release see issue #3707
bondlengths = bondlengths.astype(np.float32)
box = np.asarray(box).astype(np.float32) if box is not None else None
if numatom > 0:
if box is not None:
boxtype, box = check_box(box)
if boxtype == "ortho":
_run(
"calc_bond_distance_ortho",
args=(coords1, coords2, box, bondlengths),
backend=backend,
)
else:
_run(
"calc_bond_distance_triclinic",
args=(coords1, coords2, box, bondlengths),
backend=backend,
)
else:
_run(
"calc_bond_distance",
args=(coords1, coords2, bondlengths),
backend=backend,
)
if backend == "distopia":
# mda expects the result to be in float64, so we need to convert it back
# to float64, change for 3.0, see #3707
bondlengths = bondlengths.astype(np.float64)
if result is not None:
result[:] = bondlengths
return bondlengths
@check_coords("coords1", "coords2", "coords3", allow_atomgroup=True)
def calc_angles(
coords1: Union[npt.NDArray, "AtomGroup"],
coords2: Union[npt.NDArray, "AtomGroup"],
coords3: Union[npt.NDArray, "AtomGroup"],
box: Optional[npt.NDArray] = None,
result: Optional[npt.NDArray] = None,
backend: str = "serial",
) -> npt.NDArray:
"""Calculates the angles formed between triplets of atom positions from the
three coordinate arrays `coords1`, `coords2`, and `coords3`. All coordinate
arrays must contain the same number of coordinates.
The coordinates in `coords2` represent the apices of the angles::
2---3
/
1
Configurations where the angle is undefined (e.g., when coordinates 1 or 3
of a triplet coincide with coordinate 2) result in a value of **zero** for
that angle.
If the optional argument `box` is supplied, periodic boundaries are taken
into account when constructing the connecting vectors between coordinates,
i.e., the minimum image convention is applied for the vectors forming the
angles. Either orthogonal or triclinic boxes are supported.
If a numpy array of dtype ``numpy.float64`` with shape ``(n,)`` (for ``n``
coordinate triplets) is provided in `result`, then this preallocated array
is filled. This can speed up calculations.
Parameters
----------
coords1 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Array of shape ``(3,)`` or ``(n, 3)`` containing the coordinates of one
side of a single or ``n`` angles, respectively (dtype is arbitrary, will
be converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
coords2 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Array of shape ``(3,)`` or ``(n, 3)`` containing the coordinates of the
apices of a single or ``n`` angles, respectively (dtype is arbitrary,
will be converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
coords3 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Array of shape ``(3,)`` or ``(n, 3)`` containing the coordinates of the
other side of a single or ``n`` angles, respectively (dtype is
arbitrary, will be converted to ``numpy.float32`` internally).
Also accepts an :class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
result : numpy.ndarray, optional
Preallocated result array of dtype ``numpy.float64`` and shape ``(n,)``
(for ``n`` coordinate triplets). Avoids recreating the array in repeated
function calls.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
Returns
-------
angles : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n,)``) or numpy.float64
Array containing the angles between each triplet of coordinates. Values
are returned in radians (rad). If three single coordinates were
supplied, the angle is returned as a single number instead of an array.
.. versionadded:: 0.8
.. versionchanged:: 0.9.0
Added optional box argument to account for periodic boundaries in
calculation
.. 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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
"""
numatom = coords1.shape[0]
angles = _check_result_array(result, (numatom,))
if backend == "distopia":
# distopia requires that all the input arrays are the same type,
# while MDAnalysis allows for mixed types, this should be changed
# pre 3.0.0 release see issue #3707
angles = angles.astype(np.float32)
box = np.asarray(box).astype(np.float32) if box is not None else None
if numatom > 0:
if box is not None:
boxtype, box = check_box(box)
if boxtype == "ortho":
_run(
"calc_angle_ortho",
args=(coords1, coords2, coords3, box, angles),
backend=backend,
)
else:
_run(
"calc_angle_triclinic",
args=(coords1, coords2, coords3, box, angles),
backend=backend,
)
else:
_run(
"calc_angle",
args=(coords1, coords2, coords3, angles),
backend=backend,
)
if backend == "distopia":
# mda expects the result to be in float64, so we need to convert it back
# to float64, change for 3.0, see #3707
angles = angles.astype(np.float64)
if result is not None:
result[:] = angles
return angles
@check_coords("coords1", "coords2", "coords3", "coords4", allow_atomgroup=True)
def calc_dihedrals(
coords1: Union[npt.NDArray, "AtomGroup"],
coords2: Union[npt.NDArray, "AtomGroup"],
coords3: Union[npt.NDArray, "AtomGroup"],
coords4: Union[npt.NDArray, "AtomGroup"],
box: Optional[npt.NDArray] = None,
result: Optional[npt.NDArray] = None,
backend: str = "serial",
) -> npt.NDArray:
r"""Calculates the dihedral angles formed between quadruplets of positions
from the four coordinate arrays `coords1`, `coords2`, `coords3`, and
`coords4`, which must contain the same number of coordinates.
The dihedral angle formed by a quadruplet of positions (1,2,3,4) is
calculated around the axis connecting positions 2 and 3 (i.e., the angle
between the planes spanned by positions (1,2,3) and (2,3,4))::
4
|
2-----3
/
1
If all coordinates lie in the same plane, the cis configuration corresponds
to a dihedral angle of zero, and the trans configuration to :math:`\pi`
radians (180 degrees). Configurations where the dihedral angle is undefined
(e.g., when all coordinates lie on the same straight line) result in a value
of ``nan`` (not a number) for that dihedral.
If the optional argument `box` is supplied, periodic boundaries are taken
into account when constructing the connecting vectors between coordinates,
i.e., the minimum image convention is applied for the vectors forming the
dihedral angles. Either orthogonal or triclinic boxes are supported.
If a numpy array of dtype ``numpy.float64`` with shape ``(n,)`` (for ``n``
coordinate quadruplets) is provided in `result` then this preallocated array
is filled. This can speed up calculations.
Parameters
----------
coords1 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` containing the 1st
positions in dihedrals (dtype is arbitrary, will be converted to
``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
coords2 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` containing the 2nd
positions in dihedrals (dtype is arbitrary, will be converted to
``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
coords3 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` containing the 3rd
positions in dihedrals (dtype is arbitrary, will be converted to
``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
coords4 : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` containing the 4th
positions in dihedrals (dtype is arbitrary, will be converted to
``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
result : numpy.ndarray, optional
Preallocated result array of dtype ``numpy.float64`` and shape ``(n,)``
(for ``n`` coordinate quadruplets). Avoids recreating the array in
repeated function calls.
backend : {'serial', 'OpenMP', 'distopia'}, optional
Keyword selecting the type of acceleration.
Returns
-------
dihedrals : numpy.ndarray (``dtype=numpy.float64``, ``shape=(n,)``) or numpy.float64
Array containing the dihedral angles formed by each quadruplet of
coordinates. Values are returned in radians (rad). If four single
coordinates were supplied, the dihedral angle is returned as a single
number instead of an array. The range of dihedral angle is
:math:`(-\pi, \pi)`.
.. versionadded:: 0.8
.. versionchanged:: 0.9.0
Added optional box argument to account for periodic boundaries in
calculation
.. versionchanged:: 0.11.0
Renamed from calc_torsions to calc_dihedrals
.. 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.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
"""
numatom = coords1.shape[0]
dihedrals = _check_result_array(result, (numatom,))
if backend == "distopia":
# distopia requires that all the input arrays are the same type,
# while MDAnalysis allows for mixed types, this should be changed
# pre 3.0.0 release see issue #3707
dihedrals = dihedrals.astype(np.float32)
box = np.asarray(box).astype(np.float32) if box is not None else None
if numatom > 0:
if box is not None:
boxtype, box = check_box(box)
if boxtype == "ortho":
_run(
"calc_dihedral_ortho",
args=(coords1, coords2, coords3, coords4, box, dihedrals),
backend=backend,
)
else:
_run(
"calc_dihedral_triclinic",
args=(coords1, coords2, coords3, coords4, box, dihedrals),
backend=backend,
)
else:
_run(
"calc_dihedral",
args=(coords1, coords2, coords3, coords4, dihedrals),
backend=backend,
)
if backend == "distopia":
# mda expects the result to be in float64, so we need to convert it back
# to float64, change for 3.0, see #3707
dihedrals = dihedrals.astype(np.float64)
if result is not None:
result[:] = dihedrals
return dihedrals
@check_coords("coords", allow_atomgroup=True)
def apply_PBC(
coords: Union[npt.NDArray, "AtomGroup"],
box: Optional[npt.NDArray] = None,
backend: str = "serial",
) -> npt.NDArray:
"""Moves coordinates into the primary unit cell.
Parameters
----------
coords : numpy.ndarray or :class:`~MDAnalysis.core.groups.AtomGroup`
Coordinate array of shape ``(3,)`` or ``(n, 3)`` (dtype is arbitrary,
will be converted to ``numpy.float32`` internally). Also accepts an
:class:`~MDAnalysis.core.groups.AtomGroup`.
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
backend : {'serial', 'OpenMP'}, optional
Keyword selecting the type of acceleration.
Returns
-------
newcoords : numpy.ndarray (``dtype=numpy.float32``, ``shape=coords.shape``)
Array containing coordinates that all lie within the primary unit cell
as defined by `box`.
.. versionadded:: 0.8
.. versionchanged:: 0.13.0
Added *backend* keyword.
.. versionchanged:: 0.19.0
Internal dtype conversion of input coordinates to ``numpy.float32``.
Now also accepts (and, likewise, returns) single coordinates.
.. versionchanged:: 2.3.0
Can now accept an :class:`~MDAnalysis.core.groups.AtomGroup` as an
argument in any position and checks inputs using type hinting.
"""
# coords is an array, the check_coords decorator made sure of that.
# Mypy, however, is not aware of that so we have to tell it explicitly.
coords_array: npt.NDArray = coords # type: ignore
if len(coords_array) == 0:
return coords_array
boxtype, box = check_box(box)
if boxtype == "ortho":
_run("ortho_pbc", args=(coords_array, box), backend=backend)
else:
_run("triclinic_pbc", args=(coords_array, box), backend=backend)
return coords_array
@check_coords("vectors", enforce_copy=False, enforce_dtype=False)
def minimize_vectors(vectors: npt.NDArray, box: npt.NDArray) -> npt.NDArray:
"""Apply minimum image convention to an array of vectors
This function is required for calculating the correct vectors between two
points. A naive approach of ``ag1.positions - ag2.positions`` will not
provide the minimum vectors between particles, even if all particles are
within the primary unit cell (box).
Parameters
----------
vectors : numpy.ndarray
Vector array of shape ``(n, 3)``, either float32 or float64. These
represent many vectors (such as between two particles).
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.timestep.Timestep.dimensions`:
``[lx, ly, lz, alpha, beta, gamma]``.
Returns
-------
minimized_vectors : numpy.ndarray
Same shape and dtype as input. The vectors from the input, but
minimized according to the size of the box.
.. versionadded:: 2.1.0
"""
boxtype, box = check_box(box)
output = np.empty_like(vectors)
# use box which is same precision as input vectors
box = box.astype(vectors.dtype)
if boxtype == "ortho":
_minimize_vectors_ortho(vectors, box, output)
else:
_minimize_vectors_triclinic(vectors, box.ravel(), output)
return output
|