wrap_pos

mdtools.box.wrap_pos(pos, box, mda_backend=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). In this case a 2-dimensional position array is interpreted to contain one position per frame instead of n positions for one frame.

  • mda_backend ({None, 'serial', 'OpenMP'}, optional) – The backend to parse to MDAnalysis.lib.distances.apply_PBC(). If None, it will be set to 'OpenMP' if more than one CPU is available (as determined by mdtools.run_time_info.get_num_CPUs()) or to 'serial' if only one CPU is available.

Returns:

pos_wrapped (numpy.ndarray) – Array of the same shape as pos containing the positions that all lie within the primary unit cell as defined by box.

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 uses MDAnalysis.lib.distances.apply_PBC() to wrap the input positions into the primary unit cell. 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).

Examples

>>> pos = np.array([0, 3, 6])
>>> mdt.box.wrap_pos(pos, [4, 5, 6, 90, 90, 90])
array([0., 3., 0.], dtype=float32)
>>> pos = np.array([[0, 3, 6],
...                 [1, 5, 7]])
>>> mdt.box.wrap_pos(pos, [4, 5, 6, 90, 90, 90])
array([[0., 3., 0.],
       [1., 0., 1.]], dtype=float32)
>>> 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., 3.]], dtype=float32)
>>> pos = np.array([[[0, 3, 6],
...                  [1, 5, 7]],
...
...                 [[2, 6, 10],
...                  [3, 7, 11]]])
>>> mdt.box.wrap_pos(pos, [4, 5, 6, 90, 90, 90])
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[2., 1., 4.],
        [3., 2., 5.]]], dtype=float32)
>>> mdt.box.wrap_pos(pos, box)
array([[[0., 3., 0.],
        [1., 0., 1.]],

       [[2., 1., 2.],
        [3., 2., 3.]]], dtype=float32)