mdtraj.formats.GroTrajectoryFile

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

Interface for reading and writing to GROMACS GRO files.

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?

See also

load_gro

High-level wrapper that returns a md.Trajectory

Attributes
n_atomsint

The number of atoms in the file

topologymd.Topology

The topology. TODO(rmcgibbo) note about chain

Methods

close(self)

Close the file

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

Read data from a molecular dynamics trajectory in the GROMACS GRO format.

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

Read a trajectory from a gro file

seek(self, offset[, whence])

Move to a new file position

tell(self)

Current file position

write(self, coordinates, topology[, time, …])

Write one or more frames of a molecular dynamics trajectory to disk in the GROMACS GRO format.

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

Initialize self. See help(type(self)) for accurate signature.

Methods

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

Initialize self.

close(self)

Close the file

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

Read data from a molecular dynamics trajectory in the GROMACS GRO format.

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

Read a trajectory from a gro file

seek(self, offset[, whence])

Move to a new file position

tell(self)

Current file position

write(self, coordinates, topology[, time, …])

Write one or more frames of a molecular dynamics trajectory to disk in the GROMACS GRO format.

Attributes

distance_unit

close(self)

Close the file

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

Read data from a molecular dynamics trajectory in the GROMACS GRO format.

Parameters
n_framesint, optional

If n_frames is not None, the next n_frames of data from the file will be read. Otherwise, all of the frames in the file will be read.

strideint, optional

If stride is not None, read only every stride-th frame from disk.

atom_indicesnp.ndarray, dtype=int, optional

The specific indices of the atoms you’d like to retrieve. If not supplied, all of the atoms will be retrieved.

Returns
coordinatesnp.ndarray, shape=(n_frames, n_atoms, 3)

The cartesian coordinates of the atoms, in units of nanometers.

timenp.ndarray, None

The time corresponding to each frame, in units of picoseconds, or None if no time information is present in the trajectory.

unitcell_vectorsnp.ndarray, shape=(n_frames, 3, 3)

The box vectors in each frame, in units of nanometers

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

Read a trajectory from a gro file

Parameters
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.

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, coordinates, topology, time=None, unitcell_vectors=None, precision=3)

Write one or more frames of a molecular dynamics trajectory to disk in the GROMACS GRO format.

Parameters
coordinatesnp.ndarray, dtype=np.float32, shape=(n_frames, n_atoms, 3)

The cartesian coordinates of each atom, in units of nanometers.

topologymdtraj.Topology

The Topology defining the model to write.

timenp.ndarray, dtype=float32, shape=(n_frames), optional

The simulation time corresponding to each frame, in picoseconds. If not supplied, the numbers 0..n_frames will be written.

unitcell_vectorsnp.ndarray, dtype=float32, shape=(n_frames, 3, 3), optional

The periodic box vectors of the simulation in each frame, in nanometers.

precisionint, optional

The number of decimal places to print for coordinates. Default is 3