wrap_pos

mdtools.box.wrap_pos(pos, box, out=None, dtype=None)[source]

Shift all particle positions into the primary unit cell.

Parameters:
  • pos (array_like) – The position array which to move into the primary unit cell. Must be either of shape (3,), (n, 3) or (k, n, 3), where n is the number of particles and k is the number of frames.

  • box (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]. box can also be an array of boxes of shape (k, 6) (one box for each frame). If box has shape (k, 6) and pos has shape (n, 3), the latter will be broadcast to (k, n, 3). If box has shape (k, 6) and pos 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. out must not point to any input array.

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

Returns:

pos_wrapped (numpy.ndarray) – The input position array with all positions shifted into the primary unit cell as defined by box.

Shapes

pos

box

pos_wrapped

(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.apply_PBC()

Moves coordinates into the primary unit cell

mdtools.box.wrap()

Shift compounds of an MDAnalysis AtomGroup back into the primary unit cell

Notes

This function has the same functionality as MDAnalysis.lib.distances.apply_PBC(). But in contrast to MDAnalysis.lib.distances.apply_PBC(), this function also accepts position arrays of shape (k, n, 3) and box arrays of shape (k, 6).

Mathematically, the wrapping algorithm can be described by the following formula:

\[\Delta\mathbf{r}^w = \Delta\mathbf{r}^u - \mathbf{L} \biggl\lfloor \mathbf{L}^{-1} \Delta\mathbf{r}^u \biggr\rfloor\]

Here, \(\Delta\mathbf{r}^w\) and \(\Delta\mathbf{r}^u\) are the wrapped and unwrapped coordinates, respectively, and \(\mathbf{L} = \left( \mathbf{L}_x | \mathbf{L}_y | \mathbf{L}_z \right)\) is the matrix of the (triclinic) box vectors.

Examples

pos has shape (3,):

>>> pos = np.array([0, 3, 6])
>>> box = np.array([4, 5, 6, 90, 90, 90])
>>> mdt.box.wrap_pos(pos, box)
array([0., 3., 0.])
>>> box_mat = np.array([[4, 0, 0],
...                     [0, 5, 0],
...                     [0, 0, 6]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([0., 3., 0.])
>>> box = np.array([[4, 5, 6, 90, 90, 90],
...                 [6, 5, 4, 90, 90, 90]])
>>> mdt.box.wrap_pos(pos, box)
array([[[0., 3., 0.]],

       [[0., 3., 2.]]])
>>> box_mat = np.array([[[4, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 6]],
...
...                     [[6, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 4]]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([[[0., 3., 0.]],

       [[0., 3., 2.]]])

pos has shape (n, 3):

>>> pos = np.array([[0, 3, 6],
...                 [1, 5, 7]])
>>> box = np.array([4, 5, 6, 90, 90, 90])
>>> mdt.box.wrap_pos(pos, box)
array([[0., 3., 0.],
       [1., 0., 1.]])
>>> box_mat = np.array([[4, 0, 0],
...                     [0, 5, 0],
...                     [0, 0, 6]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([[0., 3., 0.],
       [1., 0., 1.]])
>>> box = np.array([[4, 5, 6, 90, 90, 90],
...                 [6, 5, 4, 90, 90, 90]])
>>> mdt.box.wrap_pos(pos, box)
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[0., 3., 2.],
        [1., 0., 3.]]])
>>> box_mat = np.array([[[4, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 6]],
...
...                     [[6, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 4]]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([[[0., 3., 0.],
        [1., 0., 1.]],

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

pos has shape (k, n, 3):

>>> pos = np.array([[[0, 3, 6],
...                  [1, 5, 7]],
...
...                 [[2, 6, 10],
...                  [3, 7, 11]]])
>>> box = np.array([4, 5, 6, 90, 90, 90])
>>> mdt.box.wrap_pos(pos, box)
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[2., 1., 4.],
        [3., 2., 5.]]])
>>> box_mat = np.array([[4, 0, 0],
...                     [0, 5, 0],
...                     [0, 0, 6]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[2., 1., 4.],
        [3., 2., 5.]]])
>>> box = np.array([[4, 5, 6, 90, 90, 90],
...                 [6, 5, 4, 90, 90, 90]])
>>> mdt.box.wrap_pos(pos, box)
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[2., 1., 2.],
        [3., 2., 3.]]])
>>> box_mat = np.array([[[4, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 6]],
...
...                     [[6, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 4]]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[2., 1., 2.],
        [3., 2., 3.]]])

Triclinic box:

>>> pos = np.array([[0, 3, 6],
...                 [1, 5, 7]])
>>> box = np.array([4, 5, 6, 80, 90, 100])
>>> np.round(mdt.box.wrap_pos(pos, box), 3)
array([[0.   , 1.942, 0.094],
       [1.   , 3.942, 1.094]])
>>> box_mat = np.array([[1, 0, 0],
...                     [2, 3, 0],
...                     [4, 5, 6]])
>>> mdt.box.wrap_pos(pos, box_mat)
array([[1., 1., 0.],
       [3., 3., 1.]])

out argument:

>>> pos = np.array([0, 3, 6])
>>> box = np.array([4, 5, 6, 90, 90, 90])
>>> out = np.full_like(box[:3], np.nan, dtype=np.float64)
>>> pos_wrapped = mdt.box.wrap_pos(pos, box, out=out)
>>> pos_wrapped
array([0., 3., 0.])
>>> out
array([0., 3., 0.])
>>> pos_wrapped is out
True
>>> box_mat = np.array([[4, 0, 0],
...                     [0, 5, 0],
...                     [0, 0, 6]])
>>> out = np.full(box_mat.shape[:-1], np.nan, dtype=np.float64)
>>> pos_wrapped = mdt.box.wrap_pos(pos, box_mat, out=out)
>>> pos_wrapped
array([0., 3., 0.])
>>> out
array([0., 3., 0.])
>>> pos_wrapped is out
True
>>> box = np.array([[4, 5, 6, 90, 90, 90],
...                 [6, 5, 4, 90, 90, 90]])
>>> out = np.full_like(box[:,:3], np.nan, dtype=np.float64)
>>> pos_wrapped = mdt.box.wrap_pos(pos, box, out=out)
>>> pos_wrapped
array([[[0., 3., 0.]],

       [[0., 3., 2.]]])
>>> out
array([[[0., 3., 0.]],

       [[0., 3., 2.]]])
>>> pos_wrapped is out
True
>>> box_mat = np.array([[[4, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 6]],
...
...                     [[6, 0, 0],
...                      [0, 5, 0],
...                      [0, 0, 4]]])
>>> out = np.full(box_mat.shape[:-1], np.nan, dtype=np.float64)
>>> pos_wrapped = mdt.box.wrap_pos(pos, box_mat, out=out)
>>> pos_wrapped
array([[[0., 3., 0.]],

       [[0., 3., 2.]]])
>>> out
array([[[0., 3., 0.]],

       [[0., 3., 2.]]])
>>> pos_wrapped is out
True
>>> import MDAnalysis.lib.distances as mdadist
>>> pos = np.array([[0, 3, 6],
...                 [1, 5, 7]])
>>> box = np.array([4, 5, 6, 90, 90, 90])
>>> pos_wrapped1 = mdt.box.wrap_pos(pos, box)
>>> pos_wrapped2 = mdadist.apply_PBC(pos, box)
>>> np.allclose(pos_wrapped1, pos_wrapped2, rtol=0, atol=1e-6)
True
>>> box = np.array([4, 5, 6, 80, 90, 100])
>>> pos_wrapped1 = mdt.box.wrap_pos(pos, box)
>>> pos_wrapped2 = mdadist.apply_PBC(pos, box)
>>> np.allclose(pos_wrapped1, pos_wrapped2, rtol=0, atol=1e-5)
True