mdtraj.formats.LAMMPSTrajectoryFile

class mdtraj.formats.LAMMPSTrajectoryFile(filename, mode='r', force_overwrite=True)

Interface for reading and writing to a LAMMPS lammpstrj files. This is a file-like object, that both reading or writing depending on the mode flag. It implements the context manager protocol, so you can also use it with the python ‘with’ statement.

Parameters
filenamestr

The filename to open. A path to a file on disk.

mode{‘r’, ‘w’}

The mode in which to open the file, either ‘r’ for read or ‘w’ for write.

force_overwritebool

If opened in write mode, and a file by the name of filename already exists on disk, should we overwrite it?

Methods

close(self)

Close the lammpstrj file.

parse_box(self, style)

Extract lengths and angles from a frame.

read(self[, n_frames, stride, atom_indices])

Read data from a lammpstrj file.

read_as_traj(self, topology[, n_frames, …])

Read a trajectory from a lammpstrj file

seek(self, offset[, whence])

Move to a new file position.

tell(self)

Current file position.

write(self, xyz, cell_lengths[, …])

Write one or more frames of data to a lammpstrj file.

write_box(self, lengths, angles, mins)

Write the box lines in the header of a frame.

__init__(self, filename, mode='r', force_overwrite=True)

Open a LAMMPS lammpstrj file for reading/writing.

Methods

__init__(self, filename[, mode, force_overwrite])

Open a LAMMPS lammpstrj file for reading/writing.

close(self)

Close the lammpstrj file.

parse_box(self, style)

Extract lengths and angles from a frame.

read(self[, n_frames, stride, atom_indices])

Read data from a lammpstrj file.

read_as_traj(self, topology[, n_frames, …])

Read a trajectory from a lammpstrj file

seek(self, offset[, whence])

Move to a new file position.

tell(self)

Current file position.

write(self, xyz, cell_lengths[, …])

Write one or more frames of data to a lammpstrj file.

write_box(self, lengths, angles, mins)

Write the box lines in the header of a frame.

Attributes

distance_unit

close(self)

Close the lammpstrj file.

parse_box(self, style)

Extract lengths and angles from a frame.

Parameters
stylestr

Type of box, ‘triclinic’ or ‘orthogonal’.

Returns
lengthsndarray

angles : ndarray

Notes

For more info on how LAMMPS defines boxes: http://lammps.sandia.gov/doc/Section_howto.html#howto_12

read(self, n_frames=None, stride=None, atom_indices=None)

Read data from a lammpstrj file.

Parameters
n_framesint, None

The number of frames you would like to read from the file. If None, all of the remaining frames will be loaded.

stridenp.ndarray, optional

Read only every stride-th frame.

atom_indicesarray_like, optional

If not none, then read only a subset of the atoms coordinates from the file.

Returns
xyznp.ndarray, shape=(n_frames, n_atoms, 3), dtype=np.float32
cell_lengthsnp.ndarray, None

The lengths (a,b,c) of the unit cell for each frame, or None if the information is not present in the file.

cell_anglesnp.ndarray, None

The angles (lpha, eta, gamma) defining the unit cell for each frame, or None if the information is not present in the file.

read_as_traj(self, topology, n_frames=None, stride=None, atom_indices=None)

Read a trajectory from a lammpstrj file

Parameters
topologyTopology

The system topology

n_framesint, optional

If positive, then read only the next n_frames frames. Otherwise read all of the frames in the file.

stridenp.ndarray, optional

Read only every stride-th frame.

atom_indicesarray_like, optional

If not none, then read only a subset of the atoms coordinates from the file. This may be slightly slower than the standard read because it required an extra copy, but will save memory.

Returns
trajectoryTrajectory

A trajectory object containing the loaded portion of the file.

See also

read

Returns the raw data from the file

Notes

If coordinates are specified in more than one style, the first complete trio of x/y/z coordinates will be read in according to the following order:

  1. x,y,z (unscaled coordinates)

  2. xs,ys,zs (scaled atom coordinates)

  3. xu,yu,zu (unwrapped atom coordinates)

  4. xsu,ysu,zsu (scaled unwrapped atom coordinates)

E.g., if the file contains x, y, z, xs, ys, zs then x, y, z will be used.

if the file contains x, y, xs, ys, zs then xs, ys, zs will be used.

seek(self, offset, whence=0)

Move to a new file position.

Parameters
offsetint

A number of frames.

whence{0, 1, 2}

0: offset from start of file, offset should be >=0. 1: move relative to the current position, positive or negative 2: move relative to the end of file, offset should be <= 0. Seeking beyond the end of a file is not supported

tell(self)

Current file position.

Returns
offsetint

The current frame in the file.

write(self, xyz, cell_lengths, cell_angles=None, types=None, unit_set='real')

Write one or more frames of data to a lammpstrj file.

Parameters
xyznp.ndarray, shape=(n_frames, n_atoms, 3)

The cartesian coordinates of the atoms to write. By convention, the lengths should be in units of angstroms.

cell_lengthsnp.ndarray, dtype=np.double, shape=(n_frames, 3)

The lengths (a,b,c) of the unit cell for each frame. By convention, the lengths should be in units of angstroms.

cell_anglesnp.ndarray, dtype=np.double, shape=(n_frames, 3)

The angles (lpha, eta, gamma) defining the unit cell for each frame. (Units of degrees).

typesnp.ndarray, shape(3, ), dtype=int

The numeric type of each particle.

unit_setstr, optional

The LAMMPS unit set that the simulation was performed in. See http://lammps.sandia.gov/doc/units.html for options. Currently supported unit sets: ‘real’.

write_box(self, lengths, angles, mins)

Write the box lines in the header of a frame.

Parameters
lengthsnp.ndarray, dtype=np.double, shape=(3, )

The lengths (a,b,c) of the unit cell for each frame.

anglesnp.ndarray, dtype=np.double, shape=(3, )

The angles (lpha, eta, gamma) defining the unit cell for each frame.

minsnp.ndarray, dtype=np.double, shape=(3, )

The minimum coordinates in the x-, y- and z-directions.