triclinic_box

mdtools.box.triclinic_box(box_mat, dtype=None)[source]

Convert the matrix representation of a simulation box to the length-angle representation.

Convert the matrix representation [[lx1, lx2, lx3], [[lz1, lz2, lz3]], [[lz1, lz2, lz3]]] (as returned by MDAnalysis.coordinates.base.Timestep.triclinic_dimensions) to the length-angle representation [lx, ly, lz, alpha, beta, gamma] (as returned by MDAnalysis.coordinates.base.Timestep.dimensions).

Parameters:
  • box_mat (array_like) – 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.triclinic_dimensions: [[lx1, lx2, lx3], [[lz1, lz2, lz3]], [[lz1, lz2, lz3]]]. box_mat can also be an array of box matrices with an arbitrary number of dimensions as long as the last two dimensions have shape (3, 3). Note that each box matrix must contain the box vectors as rows!

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

Returns:

box (numpy.ndarray of dtype numpy.float32) – The unit cell dimensions of the system in the same format as returned by MDAnalysis.coordinates.base.Timestep.dimensions: [lx, ly, lz, alpha, beta, gamma]. box will have one dimension less than box_mat and its last dimension will have shape (6,).

Shapes

box_mat

box

(3, 3)

(6,)

(k, 3, 3)

(k, 6)

(j, k, 3, 3)

(j, k, 6)

(i, j, k, 3, 3)

(i, j, k, 6)

See also

mdtools.box.triclinic_vectors()

Inverse function: Convert the length-angle representation of a simulation box to the matrix representation.

MDAnalysis.lib.mdamath.triclinic_box()

Similar function: Convert the matrix representation of a simulation box to the length-angle representation.

Notes

The box vectors are converted to their lengths and angles using the following formulas:

\[l_x = |\mathbf{l_x}| = \sqrt{\sum_{i=1}^3 l_{x,i}^2}\]
\[\alpha = \arccos{\frac{\mathbf{l_y}\mathbf{l_z}}{l_y l_z}}\]
\[\beta = \arccos{\frac{\mathbf{l_z}\mathbf{l_x}}{l_z l_x}}\]
\[\gamma = \arccos{\frac{\mathbf{l_x}\mathbf{l_y}}{l_x l_y}}\]

The angles are returned in degrees and are ensured to lie within the interval \(]0, 180[\).

If one of the resulting length-angle representations is invalid, a warning will be raised and a zero vector will be returned for all boxes.

Examples

>>> box_mat = np.eye(3) * 2
>>> mdt.box.triclinic_box(box_mat)
array([ 2.,  2.,  2., 90., 90., 90.])
>>> box_mat = [box_mat, np.eye(3) * 3]
>>> mdt.box.triclinic_box(box_mat)
array([[ 2.,  2.,  2., 90., 90., 90.],
       [ 3.,  3.,  3., 90., 90., 90.]])
>>> box_mat = [box_mat, [np.eye(3) * 4, np.eye(3) * 5]]
>>> mdt.box.triclinic_box(box_mat)
array([[[ 2.,  2.,  2., 90., 90., 90.],
        [ 3.,  3.,  3., 90., 90., 90.]],

       [[ 4.,  4.,  4., 90., 90., 90.],
        [ 5.,  5.,  5., 90., 90., 90.]]])
>>> import MDAnalysis.lib.mdamath as mdamath
>>> box_mat = np.eye(3) * 2
>>> box1 = mdt.box.triclinic_box(box_mat)
>>> box1
array([ 2.,  2.,  2., 90., 90., 90.])
>>> box2 = mdamath.triclinic_box(*box_mat)
>>> box2
array([ 2.,  2.,  2., 90., 90., 90.], dtype=float32)
>>> np.allclose(box1, box2)
True
>>> box_mat = np.arange(9).reshape((3, 3))
>>> box1 = mdt.box.triclinic_box(box_mat)
>>> box1
array([ 2.23606798,  7.07106781, 12.20655562,  4.88390747, 32.57846892,
       27.69456145])
>>> box2 = mdamath.triclinic_box(*box_mat)
>>> box2
array([ 2.236068 ,  7.071068 , 12.206555 ,  4.8839073, 32.57847  ,
       27.694561 ], dtype=float32)
>>> np.allclose(box1, box2)
True
>>> box_mat = np.array([[1, 0, 0], [2, 3, 0], [4, 5, 6]])
>>> box1 = mdt.box.triclinic_box(box_mat)
>>> box1
array([ 1.        ,  3.60555128,  8.77496439, 43.36781871, 62.88085723,
       56.30993247])
>>> box2 = mdamath.triclinic_box(*box_mat)
>>> box2
array([ 1.       ,  3.6055512,  8.774964 , 43.367817 , 62.880856 ,
       56.309933 ], dtype=float32)
>>> np.allclose(box1, box2)
True