vdist

mdtools.box.vdist(pos1, pos2, box=None, out=None, out_tmp=None, dtype=np.float64)[source]

Calculate the distance vectors between two position arrays.

Parameters:
  • pos1, pos2 (array_like) – The two position arrays between which to calculate the distance vectors. Position arrays must be of shape (3,), (n, 3) or (k, n, 3) where n is the number of particles and k is the number of frames. If pos1 and pos2 do not have the same shape, they must be broadcastable to a common shape. The user is responsible to provide inputs that result in a physically meaningful broadcasting!

  • box (None or array_like, optional) – The unit cell dimensions of the system, which can be orthogonal or triclinic and must be provided in the same format as returned by MDAnalysis.coordinates.base.Timestep.dimensions: [lx, ly, lz, alpha, beta, gamma]. If supplied, the minimum image convention is taken into account. box can also be an array of boxes of shape (k, 6) (one box for each frame). If box has shape (k, 6) and the result of pos1 - pos2 has shape (n, 3), the latter will be broadcast to (k, n, 3). If box has shape (k, 6) and the result of pos1 - pos2 has shape (3,), the latter will be broadcast to (k, 1, 3).

    Alternatively, box can be provided in the same format as returned by MDAnalysis.coordinates.base.Timestep.triclinic_dimensions: [[lx1, lx2, lx3], [[lz1, lz2, lz3]], [[lz1, lz2, lz3]]]. box can also be an array of boxes of shape (k, 3, 3) (one box for each frame). Equivalent broadcasting rules as above apply. Providing box already in matrix representation avoids the conversion step from the length-angle to the matrix representation.

  • out (None or numpy.ndarray, optional) – Preallocated array of the given dtype and appropriate shape into which the result is stored. If box is None, out is ignored the final result is stored in out_tmp.

  • out_tmp (None or numpy.ndarray, optional) – Preallocated array of the given dtype into which temporary results are stored. If provided, it must have the same shape as the result of pos1 - pos2. Providing out_tmp can speed up the calculated if this function is called many times.

  • dtype (type, optional) – The data type of the output array. If None, the data type is inferred from the input arrays.

Returns:

dist_vecs (numpy.ndarray) – Distance array containing the result of pos1 - pos2 with the minimum image convention taken into account if box was supplied.

Shapes

pos1 - pos2

box

dist_vecs

(3,)

(6,) or (3, 3)

(3,)

(3,)

(k, 6) or (k, 3, 3)

(k, 1, 3)

(n, 3)

(6,) or (3, 3)

(n, 3)

(n, 3)

(k, 6) or (k, 3, 3)

(k, n, 3)

(k, n, 3)

(6,) or (3, 3)

(k, n, 3)

(k, n, 3)

(k, 6) or (k, 3, 3)

(k, n, 3)

See also

MDAnalysis.lib.distances.minimize_vectors()

Apply the minimum image convention to an array of vectors

MDAnalysis.lib.distances.distance_array()

Calculate all possible distances between a reference set and another configuration

MDAnalysis.lib.distances.self_distance_array()

Calculate all possible distances within a configuration

MDAnalysis.lib.distances.calc_bonds()

Calculate the bond lengths between pairs of atom positions from the two coordinate arrays

MDAnalysis.analysis.distances.dist()

Calculate the distance between atoms in two MDAnalysis AtomGroups

mdtools.numpy_helper_functions.subtract_mic()

Subtract two arrays element-wise respecting the minium image convention

mdtools.numpy_helper_functions.diff_mic()

Calculate the difference between array elements respecting the minimum image convention.

Notes

If your are only interested in the distances itself (i.e. the norm of the distance vectors), consider using MDAnalysis.analysis.distances.dist() or MDAnalysis.lib.distances.calc_bonds() instead.

The function MDAnalysis.lib.distances.minimize_vectors() (new in MDAnalysis version 2.1) has similar functionality as this function but only accepts distance vectors of shape (n, 3) and boxes of shape (6,).

If box is provided, the minimum image convention is taken into account using algorithm C4 from Deiters[1]. This algorithm not only works for wrapped coordinates but for any image coordinates. Mathematically, it can be described by the following formula:

\[\Delta\mathbf{r}^{mic} = \Delta\mathbf{r} - \mathbf{L} \biggl\lfloor \mathbf{L}^{-1} \Delta\mathbf{r} + \mathbf{\frac{1}{2}} \biggr\rfloor\]

Here, \(\Delta\mathbf{r} = \mathbf{r}_1 - \mathbf{r}_2\) is the difference of the position vectors and \(\mathbf{L} = \left( \mathbf{L}_x | \mathbf{L}_y | \mathbf{L}_z \right)\) is the matrix of the (triclinic) box vectors.

References

Examples

Shape of pos1 and pos2 is (3,):

>>> pos1 = np.array([0, 2, 4])
>>> pos2 = np.array([5, 3, 1])
>>> mdt.box.vdist(pos1, pos2)
array([-5., -1.,  3.])
>>> # Shape of box is (6,), shape of output is (3,).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([ 1., -1., -1.])
>>> # Shape of box is (3, 3), shape of output is (3,).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([ 1., -1., -1.])
>>> # Shape of box is (k, 6), shape of output is (k, 1, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.]],

       [[-1., -1., -1.]]])
>>> # Shape of box is (k, 3, 3), shape of output is (k, 1, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.]],

       [[-1., -1., -1.]]])

Shape of pos1 and pos2 is (n, 3):

>>> pos1 = np.array([[0, 2, 4],
...                  [5, 3, 1]])
>>> pos2 = np.array([[5, 3, 1],
...                  [0, 2, 4]])
>>> mdt.box.vdist(pos1, pos2)
array([[-5., -1.,  3.],
       [ 5.,  1., -3.]])
>>> # Shape of box is (6,), shape of output is (n, 3).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[ 1., -1., -1.],
       [-1., -1., -1.]])
>>> # Shape of box is (3, 3), shape of output is (n, 3).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[ 1., -1., -1.],
       [-1., -1., -1.]])
>>> # Shape of box is (k, 6), shape of output is (k, n, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1., -1., -1.],
        [-1.,  1.,  1.]]])
>>> # Shape of box is (k, 3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1., -1., -1.],
        [-1.,  1.,  1.]]])

Shape of pos1 and pos2 is (k, n, 3):

>>> pos1 = np.array([[[0, 2, 4],
...                   [5, 3, 1]],
...
...                  [[4, 0, 2],
...                   [5, 3, 1]]])
>>> pos2 = np.array([[[5, 3, 1],
...                   [0, 2, 4]],
...
...                  [[5, 3, 1],
...                   [4, 0, 2]]])
>>> mdt.box.vdist(pos1, pos2)
array([[[-5., -1.,  3.],
        [ 5.,  1., -3.]],

       [[-1., -3.,  1.],
        [ 1.,  3., -1.]]])
>>> # Shape of box is (6,), shape of output is (k, n, 3).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1., -1., -1.],
        [ 1., -1., -1.]]])
>>> # Shape of box is (3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1., -1., -1.],
        [ 1., -1., -1.]]])
>>> # Shape of box is (k, 6), shape of output is (k, n, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1.,  0.,  1.],
        [-1.,  0., -1.]]])
>>> # Shape of box is (k, 3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1.,  0.,  1.],
        [-1.,  0., -1.]]])

Shape of pos1 is (3,) and shape of pos2 is (n, 3):

>>> pos1 = np.array([0, 2, 4])
>>> pos2 = np.array([[5, 3, 1],
...                  [0, 2, 4]])
>>> mdt.box.vdist(pos1, pos2)
array([[-5., -1.,  3.],
       [ 0.,  0.,  0.]])
>>> # Shape of box is (6,), shape of output is (n, 3).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[ 1., -1., -1.],
       [ 0.,  0.,  0.]])
>>> # Shape of box is (3, 3), shape of output is (n, 3).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[ 1., -1., -1.],
       [ 0.,  0.,  0.]])
>>> # Shape of box is (k, 6), shape of output is (k, n, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [ 0.,  0.,  0.]],

       [[-1., -1., -1.],
        [ 0.,  0.,  0.]]])
>>> # Shape of box is (k, 3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [ 0.,  0.,  0.]],

       [[-1., -1., -1.],
        [ 0.,  0.,  0.]]])

Shape of pos1 is (3,) and shape of pos2 is (k, n, 3):

>>> pos1 = np.array([0, 2, 4])
>>> pos2 = np.array([[[5, 3, 1],
...                   [0, 2, 4]],
...
...                  [[5, 3, 1],
...                   [4, 0, 2]]])
>>> mdt.box.vdist(pos1, pos2)
array([[[-5., -1.,  3.],
        [ 0.,  0.,  0.]],

       [[-5., -1.,  3.],
        [-4.,  2.,  2.]]])
>>> # Shape of box is (6,), shape of output is (k, n, 3).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [ 0.,  0.,  0.]],

       [[ 1., -1., -1.],
        [-1.,  0.,  0.]]])
>>> # Shape of box is (3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [ 0.,  0.,  0.]],

       [[ 1., -1., -1.],
        [-1.,  0.,  0.]]])
>>> # Shape of box is (k, 6), shape of output is (k, n, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [ 0.,  0.,  0.]],

       [[-1., -1., -1.],
        [ 0., -1., -2.]]])
>>> # Shape of box is (k, 3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [ 0.,  0.,  0.]],

       [[-1., -1., -1.],
        [ 0., -1., -2.]]])

Shape of pos1 is (n, 3) and shape of pos2 is (k, n, 3):

>>> pos1 = np.array([[0, 2, 4],
...                  [5, 3, 1]])
>>> pos2 = np.array([[[5, 3, 1],
...                   [0, 2, 4]],
...
...                  [[5, 3, 1],
...                   [4, 0, 2]]])
>>> mdt.box.vdist(pos1, pos2)
array([[[-5., -1.,  3.],
        [ 5.,  1., -3.]],

       [[-5., -1.,  3.],
        [ 1.,  3., -1.]]])
>>> # Shape of box is (6,), shape of output is (k, n, 3).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[ 1., -1., -1.],
        [ 1., -1., -1.]]])
>>> # Shape of box is (3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[ 1., -1., -1.],
        [ 1., -1., -1.]]])
>>> # Shape of box is (k, 6), shape of output is (k, n, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> mdt.box.vdist(pos1, pos2, box=box)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1., -1., -1.],
        [-1.,  0., -1.]]])
>>> # Shape of box is (k, 3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[[ 1., -1., -1.],
        [-1., -1., -1.]],

       [[-1., -1., -1.],
        [-1.,  0., -1.]]])

Triclinic boxes:

>>> pos1 = np.array([[0, 2, 4],
...                  [5, 3, 1]])
>>> pos2 = np.array([[5, 3, 1],
...                  [0, 2, 4]])
>>> mdt.box.vdist(pos1, pos2)
array([[-5., -1.,  3.],
       [ 5.,  1., -3.]])
>>> box = np.array([3, 2, 2, 80, 90, 100])
>>> np.round(mdt.box.vdist(pos1, pos2, box=box), 3)
array([[ 0.653,  0.264, -0.937],
       [-0.653, -0.264,  0.937]])
>>> box_mat = np.array([[1, 0, 0],
...                     [2, 3, 0],
...                     [4, 5, 6]])
>>> mdt.box.vdist(pos1, pos2, box=box_mat)
array([[-2., -3., -3.],
       [-2., -2., -3.]])

out and out_tmp arguments:

>>> # box is None.
>>> pos1 = np.array([0, 2, 4])
>>> pos2 = np.array([5, 3, 1])
>>> out_tmp = np.full_like(pos1, np.nan, dtype=np.float64)
>>> dist_vecs = mdt.box.vdist(pos1, pos2, out_tmp=out_tmp)
>>> dist_vecs
array([-5., -1.,  3.])
>>> out_tmp
array([-5., -1.,  3.])
>>> dist_vecs is out_tmp
True
>>> # Shape of box is (6,), shape of output is (3,).
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> out = np.full_like(box[:3], np.nan, dtype=np.float64)
>>> dist_vecs = mdt.box.vdist(pos1, pos2, box=box, out=out)
>>> dist_vecs
array([ 1., -1., -1.])
>>> out
array([ 1., -1., -1.])
>>> dist_vecs is out
True
>>> # Shape of box is (3, 3), shape of output is (3,).
>>> box_mat = np.array([[3, 0, 0],
...                     [0, 2, 0],
...                     [0, 0, 2]])
>>> out = np.full(box_mat.shape[:-1], np.nan, dtype=np.float64)
>>> dist_vecs = mdt.box.vdist(pos1, pos2, box=box_mat, out=out)
>>> dist_vecs
array([ 1., -1., -1.])
>>> out
array([ 1., -1., -1.])
>>> dist_vecs is out
True
>>> # Shape of box is (k, 6), shape of output is (k, 1, 3).
>>> box = np.array([[3, 2, 2, 90, 90, 90],
...                 [2, 3, 4, 90, 90, 90]])
>>> out = np.full_like(box[:,:3], np.nan, dtype=np.float64)
>>> out_tmp = np.full_like(pos1, np.nan, dtype=np.float64)
>>> dist_vecs = mdt.box.vdist(
...     pos1, pos2, box=box, out=out, out_tmp=out_tmp
... )
>>> dist_vecs
array([[[ 1., -1., -1.]],

       [[-1., -1., -1.]]])
>>> out
array([[[ 1., -1., -1.]],

       [[-1., -1., -1.]]])
>>> dist_vecs is out
True
>>> # Shape of box is (k, 3, 3), shape of output is (k, n, 3).
>>> box_mat = np.array([[[3, 0, 0],
...                      [0, 2, 0],
...                      [0, 0, 2]],
...
...                     [[2, 0, 0],
...                      [0, 3, 0],
...                      [0, 0, 4]]])
>>> out = np.full(box_mat.shape[:-1], np.nan, dtype=np.float64)
>>> out_tmp = np.full_like(pos1, np.nan, dtype=np.float64)
>>> dist_vecs = mdt.box.vdist(
...     pos1, pos2, box=box_mat, out=out, out_tmp=out_tmp
... )
>>> dist_vecs
array([[[ 1., -1., -1.]],

       [[-1., -1., -1.]]])
>>> out
array([[[ 1., -1., -1.]],

       [[-1., -1., -1.]]])
>>> dist_vecs is out
True
>>> import MDAnalysis.lib.distances as mdadist
>>> import numpy as np
>>> pos1 = np.array([[0, 2, 4],
...                  [5, 3, 1]])
>>> pos2 = np.array([[5, 3, 1],
...                  [0, 2, 4]])
>>> box = np.array([3, 2, 2, 90, 90, 90])
>>> dist_vecs1 = mdt.box.vdist(pos1, pos2, box)
>>> dists1 = np.linalg.norm(dist_vecs1, axis=-1)
>>> dists2 = mdadist.calc_bonds(pos1, pos2, box)
>>> np.allclose(dists1, dists2, rtol=0, atol=1e-6)
True
>>> box = np.array([3, 2, 2, 80, 90, 100])
>>> dist_vecs1 = mdt.box.vdist(pos1, pos2, box)
>>> dists1 = np.linalg.norm(dist_vecs1, axis=-1)
>>> dists2 = mdadist.calc_bonds(pos1, pos2, box)
>>> np.allclose(dists1, dists2, rtol=0, atol=1e-5)
True