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_frames Number of frames in the trajectory
n_atoms Number of atoms 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
topology Topology of the system, describing the organization of atoms into residues, bonds, etc
top Alias for self.topology, describing the organization of atoms into residues, bonds, etc
xyz Cartesian coordinates of each atom in each simulation frame
unitcell_vectors The vectors 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_angles Angles that define the shape of the unit cell in each frame.
__init__(xyz, topology, time=None, unitcell_lengths=None, unitcell_angles=None)

Methods

__init__(xyz, topology[, time, ...])
center_coordinates([mass_weighted]) Center each trajectory frame at the origin (0,0,0).
join(other[, check_topology, ...]) Join two trajectories together along the time/frame axis.
load(filenames, **kwargs) Load a trajectory from disk
openmm_boxes(frame) OpenMM-compatable box vectors of a single frame.
openmm_positions(frame) OpenMM-compatable positions of a single frame.
restrict_atoms(atom_indices) Retain only a subset of the atoms in a trajectory (inplace)
save(filename, **kwargs) Save trajectory to disk, in a format determined by the filename extension
save_binpos(filename[, force_overwrite]) Save trajectory to AMBER BINPOS format
save_dcd(filename[, force_overwrite]) Save trajectory to CHARMM/NAMD DCD format
save_hdf5(filename[, force_overwrite]) Save trajectory to MDTraj HDF5 format
save_lh5(filename) Save trajectory in deprecated MSMBuilder2 LH5 (lossy HDF5) format.
save_mdcrd(filename[, force_overwrite]) Save trajectory to AMBER mdcrd format
save_netcdf(filename[, force_overwrite]) Save trajectory in AMBER NetCDF format
save_pdb(filename[, force_overwrite]) Save trajectory to RCSB PDB format
save_trr(filename[, force_overwrite]) Save trajectory to Gromacs TRR format
save_xtc(filename[, force_overwrite]) Save trajectory to Gromacs XTC format
slice(key[, copy]) Slice trajectory, by extracting one or more frames into a separate object
stack(other) Stack two trajectories along the atom axis
superpose(reference[, frame, atom_indices, ...]) Superpose each conformation in this trajectory upon a reference

Attributes

n_atoms Number of atoms 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
xyz Cartesian coordinates of each atom in each simulation frame
n_frames

Number of frames in the trajectory

Returns:

n_frames : int

The number of frames in the trajectory

n_atoms

Number of atoms in the trajectory

Returns:

n_atoms : int

The number of atoms in the trajectory

n_residues

Number of residues (amino acids) in the trajectory

Returns:

n_residues : int

The number of residues in the trajectory’s topology

top

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

Returns:

topology : md.Topology

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

timestep

Timestep between frames, in picoseconds

Returns:

timestep : float

The timestep between frames, in picoseconds.

unitcell_vectors

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

Returns:

vectors : np.ndarray, shape(n_frames, 3, 3)

Vectors definiing 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, :].

superpose(reference, frame=0, atom_indices=None, parallel=True)

Superpose each conformation in this trajectory upon a reference

Parameters:

reference : md.Trajectory

For each conformation in this trajectory, aligned to a particular reference conformation in another trajectory object.

frame : int

The index of the conformation in reference to align to.

atom_indices : array_like, or None

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

parallel : bool

Use OpenMP to run the superposition in parallel over multiple cores

Returns:

self

join(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:

other : Trajectory or list of Trajectory

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

check_topology : bool

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_frames : bool, optional

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

See also

stack
join two trajectories along the atom axis
stack(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:

other : Trajectory

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
slice(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.

copy : bool, 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.

topology

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

Returns:

topology : md.Topology

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

xyz

Cartesian coordinates of each atom in each simulation frame

Returns:

xyz : np.ndarray, shape=(n_frames, n_atoms, 3)

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

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.

unitcell_angles

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

Returns:

lengths : np.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.

time

The simulation time corresponding to each frame, in picoseconds

Returns:

time : np.ndarray, shape=(n_frames,)

The simulation time corresponding to each frame, in picoseconds

openmm_positions(frame)

OpenMM-compatable positions of a single frame.

Parameters:

frame : int

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

Returns:

positions : list

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

Examples

>>> t = md.load('trajectory.h5')
>>> context.setPositions(t.openmm_positions(0))
openmm_boxes(frame)

OpenMM-compatable box vectors of a single frame.

Parameters:

frame : int

Return box for this single frame.

Returns:

box : tuple

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

Examples

>>> t = md.load('trajectory.h5')
>>> context.setPeriodicBoxVectors(t.openmm_positions(0))
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

save(filename, **kwargs)

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

Parameters:

filename : str

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

Other Parameters:
 

lossy : bool

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

no_models: bool

For .pdb. TODO: Document this?

force_overwrite : bool

For .binpos, .xtc, .dcd. If filename already exists, overwrite it.

save_hdf5(filename, force_overwrite=True)

Save trajectory to MDTraj HDF5 format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filename, if its already there

save_pdb(filename, force_overwrite=True)

Save trajectory to RCSB PDB format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filename, if its already there

save_xtc(filename, force_overwrite=True)

Save trajectory to Gromacs XTC format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filename, if its already there

save_trr(filename, force_overwrite=True)

Save trajectory to Gromacs TRR format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, 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_dcd(filename, force_overwrite=True)

Save trajectory to CHARMM/NAMD DCD format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filenames, if its already there

save_binpos(filename, force_overwrite=True)

Save trajectory to AMBER BINPOS format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filename, if its already there

save_mdcrd(filename, force_overwrite=True)

Save trajectory to AMBER mdcrd format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filename, if its already there

save_netcdf(filename, force_overwrite=True)

Save trajectory in AMBER NetCDF format

Parameters:

filename : str

filesystem path in which to save the trajectory

force_overwrite : bool, default=True

Overwrite anything that exists at filename, if its already there

save_lh5(filename)

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

Parameters:

filename : str

filesystem path in which to save the trajectory

center_coordinates(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_weighted : bool, optional (default = False)

If True, weight atoms by mass when removing COM.

Returns:

self

restrict_atoms(atom_indices)

Retain only a subset of the atoms in a trajectory (inplace)

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

Parameters:

atom_indices : list([int])

List of atom indices to keep.

Returns:

self

Versions