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
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. | 
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 | 
Number of frames in the trajectory
| Returns: | n_frames : int 
 | 
|---|
Number of atoms in the trajectory
| Returns: | n_atoms : int 
 | 
|---|
Number of residues (amino acids) in the trajectory
| Returns: | n_residues : int 
 | 
|---|
Alias for self.topology, describing the organization of atoms into residues, bonds, etc
| Returns: | topology : md.Topology 
 | 
|---|
Timestep between frames, in picoseconds
| Returns: | timestep : float 
 | 
|---|
The vectors that define the shape of the unit cell in each frame
| Returns: | vectors : np.ndarray, shape(n_frames, 3, 3) 
 | 
|---|
Superpose each conformation in this trajectory upon a reference
| Parameters: | reference : md.Trajectory 
 frame : int 
 atom_indices : array_like, or None 
 parallel : bool 
 | 
|---|---|
| Returns: | self | 
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 
 check_topology : bool 
 discard_overlapping_frames : bool, optional 
 | 
|---|
See also
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 
 | 
|---|
See also
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 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} 
 copy : bool, default=True 
 | 
|---|
Topology of the system, describing the organization of atoms into residues, bonds, etc
| Returns: | topology : md.Topology 
 | 
|---|
Cartesian coordinates of each atom in each simulation frame
| Returns: | xyz : np.ndarray, shape=(n_frames, n_atoms, 3) 
 | 
|---|
Lengths that define the shape of the unit cell in each frame.
| Returns: | lengths : {np.ndarray, shape=(n_frames, 3), None} 
 | 
|---|
Angles that define the shape of the unit cell in each frame.
| Returns: | lengths : np.ndarray, shape=(n_frames, 3) 
 | 
|---|
The simulation time corresponding to each frame, in picoseconds
| Returns: | time : np.ndarray, shape=(n_frames,) 
 | 
|---|
OpenMM-compatable positions of a single frame.
| Parameters: | frame : int 
 | 
|---|---|
| Returns: | positions : list 
 | 
Examples
>>> t = md.load('trajectory.h5')
>>> context.setPositions(t.openmm_positions(0))
OpenMM-compatable box vectors of a single frame.
| Parameters: | frame : int 
 | 
|---|---|
| Returns: | box : tuple 
 | 
Examples
>>> t = md.load('trajectory.h5')
>>> context.setPeriodicBoxVectors(t.openmm_positions(0))
Load a trajectory from disk
| Parameters: | filenames : {str, [str]} 
 | 
|---|---|
| Other Parameters: | |
| As requested by the various load functions – it depends on the extension | |
Save trajectory to disk, in a format determined by the filename extension
| Parameters: | filename : str 
 | 
|---|---|
| Other Parameters: | |
| lossy : bool 
 no_models: bool 
 force_overwrite : bool 
 | |
Save trajectory to MDTraj HDF5 format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory to RCSB PDB format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory to Gromacs XTC format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory to Gromacs TRR format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Notes
Only the xyz coordinates and the time are saved, the velocities and forces in the trr will be zeros
Save trajectory to CHARMM/NAMD DCD format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory to AMBER BINPOS format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory to AMBER mdcrd format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory in AMBER NetCDF format
| Parameters: | filename : str 
 force_overwrite : bool, default=True 
 | 
|---|
Save trajectory in deprecated MSMBuilder2 LH5 (lossy HDF5) format.
| Parameters: | filename : str 
 | 
|---|
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) 
 | 
|---|---|
| Returns: | self | 
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]) 
 | 
|---|---|
| Returns: | self |