mdtraj.Trajectory

class mdtraj.Trajectory(xyz, topology, time=None, unitcell_lengths=None, unitcell_angles=None)

Container object for a molecular dynamics trajectory

A Trajectory represents a collection of one or more molecular structures, generally (but not necessarily) from a molecular dynamics trajectory. The Trajectory stores a number of fields describing the system through time, including the cartesian coordinates of each atoms (xyz), the topology of the molecular system (topology), and information about the unitcell if appropriate (unitcell_vectors, unitcell_length, unitcell_angles).

A Trajectory should generally be constructed by loading a file from disk. Trajectories can be loaded from (and saved to) the PDB, XTC, TRR, DCD, binpos, NetCDF or MDTraj HDF5 formats.

Trajectory supports fancy indexing, so you can extract one or more frames from a Trajectory as a separate trajectory. For example, to form a trajectory with every other frame, you can slice with traj[::2].

Trajectory uses the nanometer, degree & picosecond unit system.

See also

mdtraj.load

High-level function that loads files and returns an md.Trajectory

Examples

>>> # loading a trajectory
>>> import mdtraj as md
>>> md.load('trajectory.xtc', top='native.pdb')
<mdtraj.Trajectory with 1000 frames, 22 atoms at 0x1058a73d0>
>>> # slicing a trajectory
>>> t = md.load('trajectory.h5')
>>> print(t)
<mdtraj.Trajectory with 100 frames, 22 atoms>
>>> print(t[::2])
<mdtraj.Trajectory with 50 frames, 22 atoms>
>>> # calculating the average distance between two atoms
>>> import mdtraj as md
>>> import numpy as np
>>> t = md.load('trajectory.h5')
>>> np.mean(np.sqrt(np.sum((t.xyz[:, 0, :] - t.xyz[:, 21, :])**2, axis=1)))
Attributes
n_framesint

Number of frames in the trajectory

n_atomsint

Number of atoms in the trajectory

n_residuesint

Number of residues (amino acids) in the trajectory

timenp.ndarray, shape=(n_frames,)

The simulation time corresponding to each frame, in picoseconds

timestepfloat

Timestep between frames, in picoseconds

topologymd.Topology

Topology of the system, describing the organization of atoms into residues, bonds, etc

topmd.Topology

Alias for self.topology, describing the organization of atoms into residues, bonds, etc

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

Cartesian coordinates of each atom in each simulation frame

unitcell_vectors{np.ndarray, shape=(n_frames, 3, 3), None}

The vectors that define the shape of the unit cell in each frame

unitcell_lengths{np.ndarray, shape=(n_frames, 3), None}

Lengths that define the shape of the unit cell in each frame.

unitcell_angles{np.ndarray, shape=(n_frames, 3), None}

Angles that define the shape of the unit cell in each frame.

Methods

atom_slice(self, atom_indices[, inplace])

Create a new trajectory from a subset of atoms

center_coordinates(self[, mass_weighted])

Center each trajectory frame at the origin (0,0,0).

image_molecules(self[, inplace, …])

Recenter and apply periodic boundary conditions to the molecules in each frame of the trajectory.

join(self, other[, check_topology, …])

Join two trajectories together along the time/frame axis.

load(filenames, \*\*kwargs)

Load a trajectory from disk

make_molecules_whole(self[, inplace, …])

Only make molecules whole

openmm_boxes(self, frame)

OpenMM-compatable box vectors of a single frame.

openmm_positions(self, frame)

OpenMM-compatable positions of a single frame.

remove_solvent(self[, exclude, inplace])

Create a new trajectory without solvent atoms

restrict_atoms(\*args, \*\*kwargs)

DEPRECATED: restrict_atoms was replaced by atom_slice and will be removed in 2.0

save(self, filename, \*\*kwargs)

Save trajectory to disk, in a format determined by the filename extension

save_amberrst7(self, filename[, force_overwrite])

Save trajectory in AMBER ASCII restart format

save_binpos(self, filename[, force_overwrite])

Save trajectory to AMBER BINPOS format

save_dcd(self, filename[, force_overwrite])

Save trajectory to CHARMM/NAMD DCD format

save_dtr(self, filename[, force_overwrite])

Save trajectory to DESMOND DTR format

save_gro(self, filename[, force_overwrite, …])

Save trajectory in Gromacs .gro format

save_hdf5(self, filename[, force_overwrite])

Save trajectory to MDTraj HDF5 format

save_lammpstrj(self, filename[, force_overwrite])

Save trajectory to LAMMPS custom dump format

save_lh5(self, filename[, force_overwrite])

Save trajectory in deprecated MSMBuilder2 LH5 (lossy HDF5) format.

save_mdcrd(self, filename[, force_overwrite])

Save trajectory to AMBER mdcrd format

save_netcdf(self, filename[, force_overwrite])

Save trajectory in AMBER NetCDF format

save_netcdfrst(self, filename[, force_overwrite])

Save trajectory in AMBER NetCDF restart format

save_pdb(self, filename[, force_overwrite, …])

Save trajectory to RCSB PDB format

save_tng(self, filename[, force_overwrite])

Save trajectory to Gromacs TNG format

save_trr(self, filename[, force_overwrite])

Save trajectory to Gromacs TRR format

save_xtc(self, filename[, force_overwrite])

Save trajectory to Gromacs XTC format

save_xyz(self, filename[, force_overwrite])

Save trajectory to .xyz format.

slice(self, key[, copy])

Slice trajectory, by extracting one or more frames into a separate object

smooth(self, width[, order, atom_indices, …])

Smoothen a trajectory using a zero-delay Buttersworth filter.

stack(self, other)

Stack two trajectories along the atom axis

superpose(self, reference[, frame, …])

Superpose each conformation in this trajectory upon a reference

__init__(self, xyz, topology, time=None, unitcell_lengths=None, unitcell_angles=None)

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

Methods

__init__(self, xyz, topology[, time, …])

Initialize self.

atom_slice(self, atom_indices[, inplace])

Create a new trajectory from a subset of atoms

center_coordinates(self[, mass_weighted])

Center each trajectory frame at the origin (0,0,0).

image_molecules(self[, inplace, …])

Recenter and apply periodic boundary conditions to the molecules in each frame of the trajectory.

join(self, other[, check_topology, …])

Join two trajectories together along the time/frame axis.

load(filenames, \*\*kwargs)

Load a trajectory from disk

make_molecules_whole(self[, inplace, …])

Only make molecules whole

openmm_boxes(self, frame)

OpenMM-compatable box vectors of a single frame.

openmm_positions(self, frame)

OpenMM-compatable positions of a single frame.

remove_solvent(self[, exclude, inplace])

Create a new trajectory without solvent atoms

restrict_atoms(\*args, \*\*kwargs)

DEPRECATED: restrict_atoms was replaced by atom_slice and will be removed in 2.0

save(self, filename, \*\*kwargs)

Save trajectory to disk, in a format determined by the filename extension

save_amberrst7(self, filename[, force_overwrite])

Save trajectory in AMBER ASCII restart format

save_binpos(self, filename[, force_overwrite])

Save trajectory to AMBER BINPOS format

save_dcd(self, filename[, force_overwrite])

Save trajectory to CHARMM/NAMD DCD format

save_dtr(self, filename[, force_overwrite])

Save trajectory to DESMOND DTR format

save_gro(self, filename[, force_overwrite, …])

Save trajectory in Gromacs .gro format

save_hdf5(self, filename[, force_overwrite])

Save trajectory to MDTraj HDF5 format

save_lammpstrj(self, filename[, force_overwrite])

Save trajectory to LAMMPS custom dump format

save_lh5(self, filename[, force_overwrite])

Save trajectory in deprecated MSMBuilder2 LH5 (lossy HDF5) format.

save_mdcrd(self, filename[, force_overwrite])

Save trajectory to AMBER mdcrd format

save_netcdf(self, filename[, force_overwrite])

Save trajectory in AMBER NetCDF format

save_netcdfrst(self, filename[, force_overwrite])

Save trajectory in AMBER NetCDF restart format

save_pdb(self, filename[, force_overwrite, …])

Save trajectory to RCSB PDB format

save_tng(self, filename[, force_overwrite])

Save trajectory to Gromacs TNG format

save_trr(self, filename[, force_overwrite])

Save trajectory to Gromacs TRR format

save_xtc(self, filename[, force_overwrite])

Save trajectory to Gromacs XTC format

save_xyz(self, filename[, force_overwrite])

Save trajectory to .xyz format.

slice(self, key[, copy])

Slice trajectory, by extracting one or more frames into a separate object

smooth(self, width[, order, atom_indices, …])

Smoothen a trajectory using a zero-delay Buttersworth filter.

stack(self, other)

Stack two trajectories along the atom axis

superpose(self, reference[, frame, …])

Superpose each conformation in this trajectory upon a reference

Attributes

n_atoms

Number of atoms in the trajectory

n_chains

Number of chains in the trajectory

n_frames

Number of frames in the trajectory

n_residues

Number of residues (amino acids) in the trajectory

time

The simulation time corresponding to each frame, in picoseconds

timestep

Timestep between frames, in picoseconds

top

Alias for self.topology, describing the organization of atoms into residues, bonds, etc

topology

Topology of the system, describing the organization of atoms into residues, bonds, etc

unitcell_angles

Angles that define the shape of the unit cell in each frame.

unitcell_lengths

Lengths that define the shape of the unit cell in each frame.

unitcell_vectors

The vectors that define the shape of the unit cell in each frame

unitcell_volumes

Volumes of unit cell for each frame.

xyz

Cartesian coordinates of each atom in each simulation frame

atom_slice(self, atom_indices, inplace=False)

Create a new trajectory from a subset of atoms

Parameters
atom_indicesarray-like, dtype=int, shape=(n_atoms)

List of indices of atoms to retain in the new trajectory.

inplacebool, default=False

If True, the operation is done inplace, modifying self. Otherwise, a copy is returned with the sliced atoms, and self is not modified.

Returns
trajmd.Trajectory

The return value is either self, or the new trajectory, depending on the value of inplace.

See also

stack

stack multiple trajectories along the atom axis

center_coordinates(self, mass_weighted=False)

Center each trajectory frame at the origin (0,0,0).

This method acts inplace on the trajectory. The centering can be either uniformly weighted (mass_weighted=False) or weighted by the mass of each atom (mass_weighted=True).

Parameters
mass_weightedbool, optional (default = False)

If True, weight atoms by mass when removing COM.

Returns
self
image_molecules(self, inplace=False, anchor_molecules=None, other_molecules=None, sorted_bonds=None, make_whole=True)

Recenter and apply periodic boundary conditions to the molecules in each frame of the trajectory.

This method is useful for visualizing a trajectory in which molecules were not wrapped to the periodic unit cell, or in which the macromolecules are not centered with respect to the solvent. It tries to be intelligent in deciding what molecules to center, so you can simply call it and trust that it will “do the right thing”.

Parameters
inplacebool, default=False

If False, a new Trajectory is created and returned. If True, this Trajectory is modified directly.

anchor_moleculeslist of atom sets, optional, default=None

Molecule that should be treated as an “anchor”. These molecules will be centered in the box and put near each other. If not specified, anchor molecules are guessed using a heuristic.

other_moleculeslist of atom sets, optional, default=None

Molecules that are not anchors. If not specified, these will be molecules other than the anchor molecules

sorted_bondsarray of shape (n_bonds, 2)

Pairs of atom indices that define bonds, in sorted order. If not specified, these will be determined from the trajectory’s topology. Only relevant if make_whole is True.

make_wholebool

Whether to make molecules whole.

Returns
trajmd.Trajectory

The return value is either self or the new trajectory, depending on the value of inplace.

join(self, other, check_topology=True, discard_overlapping_frames=False)

Join two trajectories together along the time/frame axis.

This method joins trajectories along the time axis, giving a new trajectory of length equal to the sum of the lengths of self and other. It can also be called by using self + other

Parameters
otherTrajectory or list of Trajectory

One or more trajectories to join with this one. These trajectories are appended to the end of this trajectory.

check_topologybool

Ensure that the topology of self and other are identical before joining them. If false, the resulting trajectory will have the topology of self.

discard_overlapping_framesbool, optional

If True, compare coordinates at trajectory edges to discard overlapping frames. Default: False.

See also

stack

join two trajectories along the atom axis

static load(filenames, **kwargs)

Load a trajectory from disk

Parameters
filenames{str, [str]}

Either a string or list of strings

Other Parameters
As requested by the various load functions – it depends on the extension
make_molecules_whole(self, inplace=False, sorted_bonds=None)

Only make molecules whole

Parameters
inplacebool

If False, a new Trajectory is created and returned. If True, this Trajectory is modified directly.

sorted_bondsarray of shape (n_bonds, 2)

Pairs of atom indices that define bonds, in sorted order. If not specified, these will be determined from the trajectory’s topology.

See also

image_molecules
property n_atoms

Number of atoms in the trajectory

Returns
n_atomsint

The number of atoms in the trajectory

property n_chains

Number of chains in the trajectory

Returns
n_chainsint

The number of chains in the trajectory’s topology

property n_frames

Number of frames in the trajectory

Returns
n_framesint

The number of frames in the trajectory

property n_residues

Number of residues (amino acids) in the trajectory

Returns
n_residuesint

The number of residues in the trajectory’s topology

openmm_boxes(self, frame)

OpenMM-compatable box vectors of a single frame.

Parameters
frameint

Return box for this single frame.

Returns
boxtuple

The periodic box vectors for this frame, formatted for input to OpenMM.

Examples

>>> t = md.load('trajectory.h5')
>>> context.setPeriodicBoxVectors(t.openmm_positions(0))
openmm_positions(self, frame)

OpenMM-compatable positions of a single frame.

Parameters
frameint

The index of frame of the trajectory that you wish to extract

Returns
positionslist

The cartesian coordinates of specific trajectory frame, formatted for input to OpenMM

Examples

>>> t = md.load('trajectory.h5')
>>> context.setPositions(t.openmm_positions(0))
remove_solvent(self, exclude=None, inplace=False)

Create a new trajectory without solvent atoms

Parameters
excludearray-like, dtype=str, shape=(n_solvent_types)

List of solvent residue names to retain in the new trajectory.

inplacebool, default=False

The return value is either self, or the new trajectory, depending on the value of inplace.

Returns
trajmd.Trajectory

The return value is either self, or the new trajectory, depending on the value of inplace.

restrict_atoms(*args, **kwargs)

DEPRECATED: restrict_atoms was replaced by atom_slice and will be removed in 2.0

Retain only a subset of the atoms in a trajectory

Deletes atoms not in atom_indices, and re-indexes those that remain

Parameters
atom_indicesarray-like, dtype=int, shape=(n_atoms)

List of atom indices to keep.

inplacebool, default=True

If True, the operation is done inplace, modifying self. Otherwise, a copy is returned with the restricted atoms, and self is not modified.

Returns
trajmd.Trajectory

The return value is either self, or the new trajectory, depending on the value of inplace.

save(self, filename, **kwargs)

Save trajectory to disk, in a format determined by the filename extension

Parameters
filenamestr

filesystem path in which to save the trajectory. The extension will be parsed and will control the format.

Other Parameters
lossybool

For .h5 or .lh5, whether or not to use compression.

no_models: bool

For .pdb. TODO: Document this?

force_overwritebool

If filename already exists, overwrite it.

save_amberrst7(self, filename, force_overwrite=True)

Save trajectory in AMBER ASCII restart format

Parameters
filenamestr

filesystem path in which to save the restart

force_overwritebool, default=True

Overwrite anything that exists at filename, if it’s already there

Notes

Amber restart files can only store a single frame. If only one frame exists, “filename” will be written. Otherwise, “filename.#” will be written, where # is a zero-padded number from 1 to the total number of frames in the trajectory

save_binpos(self, filename, force_overwrite=True)

Save trajectory to AMBER BINPOS format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

save_dcd(self, filename, force_overwrite=True)

Save trajectory to CHARMM/NAMD DCD format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filenames, if its already there

save_dtr(self, filename, force_overwrite=True)

Save trajectory to DESMOND DTR format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filenames, if its already there

save_gro(self, filename, force_overwrite=True, precision=3)

Save trajectory in Gromacs .gro format

Parameters
filenamestr

Path to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at that filename if it exists

precisionint, default=3

The number of decimal places to use for coordinates in GRO file

save_hdf5(self, filename, force_overwrite=True)

Save trajectory to MDTraj HDF5 format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

save_lammpstrj(self, filename, force_overwrite=True)

Save trajectory to LAMMPS custom dump format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

save_lh5(self, filename, force_overwrite=True)

Save trajectory in deprecated MSMBuilder2 LH5 (lossy HDF5) format.

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if it’s already there

save_mdcrd(self, filename, force_overwrite=True)

Save trajectory to AMBER mdcrd format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

save_netcdf(self, filename, force_overwrite=True)

Save trajectory in AMBER NetCDF format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if it’s already there

save_netcdfrst(self, filename, force_overwrite=True)

Save trajectory in AMBER NetCDF restart format

Parameters
filenamestr

filesystem path in which to save the restart

force_overwritebool, default=True

Overwrite anything that exists at filename, if it’s already there

Notes

NetCDF restart files can only store a single frame. If only one frame exists, “filename” will be written. Otherwise, “filename.#” will be written, where # is a zero-padded number from 1 to the total number of frames in the trajectory

save_pdb(self, filename, force_overwrite=True, bfactors=None)

Save trajectory to RCSB PDB format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

bfactorsarray_like, default=None, shape=(n_frames, n_atoms) or (n_atoms,)

Save bfactors with pdb file. If the array is two dimensional it should contain a bfactor for each atom in each frame of the trajectory. Otherwise, the same bfactor will be saved in each frame.

save_tng(self, filename, force_overwrite=True)

Save trajectory to Gromacs TNG format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

save_trr(self, filename, force_overwrite=True)

Save trajectory to Gromacs TRR format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

Notes

Only the xyz coordinates and the time are saved, the velocities and forces in the trr will be zeros

save_xtc(self, filename, force_overwrite=True)

Save trajectory to Gromacs XTC format

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

save_xyz(self, filename, force_overwrite=True)

Save trajectory to .xyz format.

Parameters
filenamestr

filesystem path in which to save the trajectory

force_overwritebool, default=True

Overwrite anything that exists at filename, if its already there

slice(self, key, copy=True)

Slice trajectory, by extracting one or more frames into a separate object

This method can also be called using index bracket notation, i.e traj[1] == traj.slice(1)

Parameters
key{int, np.ndarray, slice}

The slice to take. Can be either an int, a list of ints, or a slice object.

copybool, default=True

Copy the arrays after slicing. If you set this to false, then if you modify a slice, you’ll modify the original array since they point to the same data.

smooth(self, width, order=3, atom_indices=None, inplace=False)

Smoothen a trajectory using a zero-delay Buttersworth filter. Please note that for optimal results the trajectory should be properly aligned prior to smoothing (see md.Trajectory.superpose).

Parameters
widthint

This acts very similar to the window size in a moving average smoother. In this implementation, the frequency of the low-pass filter is taken to be two over this width, so it’s like “half the period” of the sinusiod where the filter starts to kick in. Must be an integer greater than one.

orderint, optional, default=3

The order of the filter. A small odd number is recommended. Higher order filters cutoff more quickly, but have worse numerical properties.

atom_indicesarray-like, dtype=int, shape=(n_atoms), default=None

List of indices of atoms to retain in the new trajectory. Default is set to None, which applies smoothing to all atoms.

inplacebool, default=False

The return value is either self, or the new trajectory, depending on the value of inplace.

Returns
trajmd.Trajectory

The return value is either self, or the new smoothed trajectory, depending on the value of inplace.

References

1

“FiltFilt”. Scipy Cookbook. SciPy. <http://www.scipy.org/Cookbook/FiltFilt>.

stack(self, other)

Stack two trajectories along the atom axis

This method joins trajectories along the atom axis, giving a new trajectory with a number of atoms equal to the sum of the number of atoms in self and other.

Parameters
otherTrajectory

The other trajectory to join

See also

join

join two trajectories along the time/frame axis.

Notes

The resulting trajectory will have the unitcell and time information the left operand.

Examples

>>> t1 = md.load('traj1.h5')
>>> t2 = md.load('traj2.h5')
>>> # even when t2 contains no unitcell information
>>> t2.unitcell_vectors = None
>>> stacked = t1.stack(t2)
>>> # the stacked trajectory inherits the unitcell information
>>> # from the first trajectory
>>> np.all(stacked.unitcell_vectors == t1.unitcell_vectors)
True
superpose(self, reference, frame=0, atom_indices=None, ref_atom_indices=None, parallel=True)

Superpose each conformation in this trajectory upon a reference

Parameters
referencemd.Trajectory

Align self to a particular frame in reference

frameint

The index of the conformation in reference to align to.

atom_indicesarray_like, or None

The indices of the atoms to superpose. If not supplied, all atoms will be used.

ref_atom_indicesarray_like, or None

Use these atoms on the reference structure. If not supplied, the same atom indices will be used for this trajectory and the reference one.

parallelbool

Use OpenMP to run the superposition in parallel over multiple cores

Returns
self
property time

The simulation time corresponding to each frame, in picoseconds

Returns
timenp.ndarray, shape=(n_frames,)

The simulation time corresponding to each frame, in picoseconds

property timestep

Timestep between frames, in picoseconds

Returns
timestepfloat

The timestep between frames, in picoseconds.

property top

Alias for self.topology, describing the organization of atoms into residues, bonds, etc

Returns
topologymd.Topology

The topology object, describing the organization of atoms into residues, bonds, etc

property topology

Topology of the system, describing the organization of atoms into residues, bonds, etc

Returns
topologymd.Topology

The topology object, describing the organization of atoms into residues, bonds, etc

property unitcell_angles

Angles that define the shape of the unit cell in each frame.

Returns
lengthsnp.ndarray, shape=(n_frames, 3)

The angles between the three unitcell vectors in each frame, alpha, beta, and gamma. alpha' gives the angle between vectors ``b and c, beta gives the angle between vectors c and a, and gamma gives the angle between vectors a and b. The angles are in degrees.

property unitcell_lengths

Lengths that define the shape of the unit cell in each frame.

Returns
lengths{np.ndarray, shape=(n_frames, 3), None}

Lengths of the unit cell in each frame, in nanometers, or None if the Trajectory contains no unitcell information.

property unitcell_vectors

The vectors that define the shape of the unit cell in each frame

Returns
vectorsnp.ndarray, shape(n_frames, 3, 3)

Vectors defining the shape of the unit cell in each frame. The semantics of this array are that the shape of the unit cell in frame i are given by the three vectors, value[i, 0, :], value[i, 1, :], and value[i, 2, :].

property unitcell_volumes

Volumes of unit cell for each frame.

Returns
volumes{np.ndarray, shape=(n_frames), None}

Volumes of the unit cell in each frame, in nanometers^3, or None if the Trajectory contains no unitcell information.

property xyz

Cartesian coordinates of each atom in each simulation frame

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

A three dimensional numpy array, with the cartesian coordinates of each atoms in each frame.