Classes to read and write DCD binary trajectories, the format used by CHARMM, NAMD but also LAMMPS. Trajectories can be read regardless of system-endianness as this is auto-detected.
Generally, DCD trajectories produced by any code can be read (with the DCDReader) although there can be issues with the unitcell (simulation box) representation (see Timestep.dimensions). DCDs can also be written but the DCDWriter follows recent NAMD/VMD convention for the unitcell but still writes AKMA time. Reading and writing these trajectories within MDAnalysis will work seamlessly but if you process those trajectories with other tools you might need to watch out that time and unitcell dimensions are correctly interpreted.
Note
The DCD file format is not well defined. In particular, NAMD and CHARMM use it differently. Currently, MDAnalysis tries to guess the correct format for the unitcell representation but it can be wrong. Check the unitcell dimensions, especially for triclinic unitcells (see Issue 187 and Timestep.dimensions). A second potential issue are the units of time which are AKMA for the DCDReader (following CHARMM) but ps for NAMD. As a workaround one can employ the configurable MDAnalysis.coordinates.LAMMPS.DCDReader for NAMD trajectories.
See also
The MDAnalysis.coordinates.LAMMPS module provides a more flexible DCD reader/writer.
The classes in this module are the reference implementations for the Trajectory API.
Make a new Timestep containing a subset of the original Timestep.
ts.copy_slice(slice(start, stop, skip)) ts.copy_slice([list of indices])
Returns: | A Timestep object of the same type containing all header information and all atom information relevent to the selection. |
---|
Note
The selection must be a 0 based slice or array of the atom indices in this Timestep
New in version 0.8.
unitcell dimensions (A, B, C, alpha, beta, gamma)
lengths A, B, C are in the MDAnalysis length unit (Å), and angles are in degrees.
dimensions is read-only because it transforms the actual format of the unitcell (which differs between different trajectory formats) to the representation described here, which is used everywhere in MDAnalysis.
The ordering of the angles in the unitcell is the same as in recent versions of VMD’s DCDplugin (2013), namely the X-PLOR DCD format: The original unitcell is read as [A, gamma, B, beta, alpha, C] from the DCD file (actually, the direction cosines are stored instead of the angles but the underlying C code already does this conversion); if any of these values are < 0 or if any of the angles are > 180 degrees then it is assumed it is a new-style CHARMM unitcell (at least since c36b2) in which box vectors were recorded.
Warning
The DCD format is not well defined. Check your unit cell dimensions carefully, especially when using triclinic boxes. Different software packages implement different conventions and MDAnalysis is currently implementing the newer NAMD/VMD convention and tries to guess the new CHARMM one. Old CHARMM trajectories might give wrong unitcell values. For more details see Issue 187.
Changed in version 0.9.0: Unitcell is now interpreted in the newer NAMD DCD format as [A, gamma, B, beta, alpha, C] instead of the old MDAnalysis/CHARMM ordering [A, alpha, B, beta, gamma, C]. We attempt to detect the new CHARMM DCD unitcell format (see Issue 187 for a discussion).
volume of the unitcell
Reads from a DCD file
Data: |
|
---|---|
Methods: |
|
Note
The DCD file format is not well defined. In particular, NAMD and CHARMM use it differently. Currently, MDAnalysis tries to guess the correct format for the unitcell representation but it can be wrong. Check the unitcell dimensions, especially for triclinic unitcells (see Issue 187 and Timestep.dimensions). A second potential issue are the units of time (TODO).
Changed in version 0.9.0: The underlying DCD reader (written in C and derived from the catdcd/molfile plugin code of VMD) is now reading the unitcell in NAMD ordering: [A, B, C, sin(gamma), sin(beta), sin(alpha)]. See Issue 187 for further details.
Returns a writer appropriate for filename.
Sets the default keywords start, step and delta (if available). numatoms is always set from Reader.numatoms.
See also
Reader.Writer() and MDAnalysis.Writer()
Returns a DCDWriter for filename with the same parameters as this DCD.
All values can be changed through keyword arguments.
Arguments: |
|
---|---|
Keywords: |
|
Returns: |
Note
The keyword arguments set the low-level attributes of the DCD according to the CHARMM format. The time between two frames would be delta * step !
See also
DCDWriter has detailed argument description
Specific implementation of trajectory closing.
In-place conversion of forces array force from native units to base units.
By default, the input force is modified in place and also returned.
New in version 0.7.7.
In-place conversion of force array force from base units to native units.
By default, the input force is modified in place and also returned.
New in version 0.7.7.
In-place conversion of coordinate array x from native units to base units.
By default, the input x is modified in place and also returned.
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
Conversion of coordinate array x from base units to native units.
By default, the input x is modified in place and also returned.
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
Convert time t from native units to base units.
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
Convert time t from base units to native units.
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
In-place conversion of velocities array v from native units to base units.
By default, the input v is modified in place and also returned.
New in version 0.7.5.
In-place conversion of coordinate array v from base units to native units.
By default, the input v is modified in place and also returned.
New in version 0.7.5.
Populate a TimeseriesCollection object with timeseries computed from the trajectory
Arguments: |
|
---|
Time between two trajectory frames in picoseconds.
Frame number of the current time step.
This is a simple short cut to Timestep.frame.
Forward one step to next frame.
Position at beginning of trajectory
Time of the current frame in MDAnalysis time units (typically ps).
time = Timestep.frame * Reader.dt
Return a subset of coordinate data for an AtomGroup
Arguments: |
|
---|
Total length of the trajectory numframes * dt.
Writes to a DCD file
Typical usage:
with DCDWriter("new.dcd", u.atoms.numberOfAtoms()) as w:
for ts in u.trajectory
w.write_next_timestep(ts)
Keywords are available to set some of the low-level attributes of the DCD.
Methods: | d = DCDWriter(dcdfilename, numatoms, start, step, delta, remarks) |
---|
Note
The Writer will write the unit cell information to the DCD in a format compatible with NAMD and older CHARMM versions, namely the unit cell lengths in Angstrom and the angle cosines (see Timestep). Newer versions of CHARMM (at least c36b2) store the matrix of the box vectors. Writing this matrix to a DCD is currently not supported (although reading is supported with the DCDReader); instead the angle cosines are written, which might make the DCD file unusable in CHARMM itself. See Issue 187 for further information.
The writing behavior of the DCDWriter is identical to that of the DCD molfile plugin of VMD with the exception that by default it will use AKMA time units.
Create a new DCDWriter
Arguments:
- filename
name of output file
- numatoms
number of atoms in dcd file
- start
starting timestep
- step
skip between subsequent timesteps (indicate that step MD integrator steps (!) make up one trajectory frame); default is 1.
- delta
timestep (MD integrator time step (!), in AKMA units); default is 20.45482949774598 (corresponding to 1 ps).
- remarks
comments to annotate dcd file
- dt
Override step and delta so that the DCD records that dt ps lie between two frames. (It sets step = 1 and delta = AKMA(dt).) The default is None, in which case step and delta are used.
- convert_units
units are converted to the MDAnalysis base format; None selects the value of MDAnalysis.core.flags [‘convert_lengths’]. (see Flags)
Note
The keyword arguments set the low-level attributes of the DCD according to the CHARMM format. The time between two frames would be delta * step ! For convenience, one can alternatively supply the dt keyword (see above) to just tell the writer that it should record “There are dt ps between each frame”.
Specific implementation of trajectory closing.
Read dimensions from timestep ts and return appropriate native unitcell.
See also
In-place conversion of forces array force from native units to base units.
By default, the input force is modified in place and also returned.
New in version 0.7.7.
In-place conversion of force array force from base units to native units.
By default, the input force is modified in place and also returned.
New in version 0.7.7.
In-place conversion of coordinate array x from native units to base units.
By default, the input x is modified in place and also returned.
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
Conversion of coordinate array x from base units to native units.
By default, the input x is modified in place and also returned.
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
Convert time t from native units to base units.
By default, the input t is modified in place and also returned (although note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
Convert time t from base units to native units.
By default, the input t is modified in place and also returned. (Also note that scalar values t are passed by value in Python and hence an in-place modification has no effect on the caller.)
Changed in version 0.7.5: Keyword inplace can be set to False so that a modified copy is returned unless no conversion takes place, in which case the reference to the unmodified x is returned.
In-place conversion of velocities array v from native units to base units.
By default, the input v is modified in place and also returned.
New in version 0.7.5.
In-place conversion of coordinate array v from base units to native units.
By default, the input v is modified in place and also returned.
New in version 0.7.5.
Returns True if all values are within limit values of their formats.
Due to rounding, the test is asymmetric (and min is supposed to be negative):
min < x <= max
Arguments: |
|
---|---|
Returns: | boolean |
Write current timestep, using the supplied obj.
The argument should be a AtomGroup or a Universe or a Timestep instance.
Note
The size of the obj must be the same as the number of atom provided when setting up the trajectory.
write a new timestep to the dcd file
ts - timestep object containing coordinates to be written to dcd file
Changed in version 0.7.5: Raises ValueError instead of generic Exception if wrong number of atoms supplied and NoDataError if no coordinates to be written.