Table Of Contents

Previous topic

5.16. pdb.extensions – Extensions to Bio.PDB

Next topic

5.18. TRZ trajectory I/O — MDAnalysis.coordinates.TRZ

This Page

5.17. The Gromacs xtc/trr library libxdrfile2

libxdrfile2 provides an interface to some high-level functions in the libxdrfile2 library, a derivative of the Gromacs libxdrfile library. Only functions required for reading and processing whole trajectories are exposed at the moment; low-level routines to read individual numbers are not provided. In addition, libxdrfile2 exposes functions to allow fast frame indexing and XDR file seeking.

The functions querying the numbers of atoms in a trajectory frame (read_xtc_natoms() and read_trr_natoms()) open a file themselves and only require the file name.

All other functions operate on a XDRFILE object, which is a special file handle for xdr files. Any xdr-based trajectory file (xtc or trr format) always has to be opened with xdrfile_open(). When done, close the trajectory with xdrfile_close().

The functions fill or read existing arrays of coordinates; they never allocate these arrays themselves. Hence they need to be setup outside libxdrfile2 as numpy arrays. The exception to these are the indexing ones functions that take care of array allocation and transference to a garbage-collectable memory object.

Changed in version 0.8.0.

libxdrfile2 is now used instead of libxdrfile. libxdrfile2 is based on libxdrfile but interfaces with libxdrfile2, which has xdr seeking and indexing capabilities. Unlike libxdrfile and libxdrfile before them, libxdrfile2 and libxdrfile2 are distributed under the GNU GENERAL PUBLIC LICENSE, version 2 (or higher).

5.17.1. Example: Reading from a xtc

In the example we read coordinate frames from an existing xtc trajectory:

import numpy as np
from libxdrfile2 import xdrfile_open, xdrfile_close, read_xtc_natoms, read_xtc, DIM, exdrOK
xtc = 'md.xtc'

# get number of atoms
natoms = read_xtc_natoms(xtc)

# allocate coordinate array of the right size and type
# (the type float32 is crucial to match the underlying C-code!!)
x = np.zeros((natoms, DIM), dtype=np.float32)
# allocate unit cell box
box = np.zeros((DIM, DIM), dtype=np.float32)

# open file
XTC = xdrfile_open(xtc, 'r')

# loop through file until return status signifies end or a problem
# (it should become exdrENDOFFILE on the last iteration)
status = exdrOK
while status == exdrOK:
   status,step,time,prec = read_xtc(XTC, box, x)
   # do something with x
   centre = x.mean(axis=0)
   print 'Centre of geometry at %(time)g ps: %(centre)r' % vars()

# finally close file
xdrfile_close(XTC)

Note that only the contents of the coordinate and unitcell arrays x and box change.

5.17.2. Functions and constants

The module defines a number of constants such as DIM or the Status symbols.

MDAnalysis.coordinates.xdrfile.libxdrfile2.DIM

The number of cartesian dimensions for which the underlying C-code was compiled; this is most certainly 3.

5.17.2.1. Status symbols

A number of symbols are exported; they all start with the letters exdr. Important ones are listed here:

MDAnalysis.coordinates.xdrfile.libxdrfile2.exdrOK

Success of xdr file read/write operation.

MDAnalysis.coordinates.xdrfile.libxdrfile2.exdrCLOSE

xdr file is closed

MDAnalysis.coordinates.xdrfile.libxdrfile2.exdrENDOFFILE

end of file was reached (response of read_xtc() and read_trr() after the last read frame)

MDAnalysis.coordinates.xdrfile.libxdrfile2.exdrFILENOTFOUND

xdrfile_open() cannot find the requested file

5.17.2.2. Opening and closing of XDR files

Two low-level functions are used to obtain a XDRFILE object (a file handle) to access xdr files such as xtc or trr trajectories.

MDAnalysis.coordinates.xdrfile.libxdrfile2.xdrfile_open(path, mode) → XDRFILE[source]

Open path and returns a XDRFILE handle that is required by other functions.

Arguments:
path

file name

mode

‘r’ for reading and ‘w’ for writing

Returns:

XDRFILE handle

MDAnalysis.coordinates.xdrfile.libxdrfile2.xdrfile_close(XDRFILE) → status[source]

Close the xdrfile pointed to by XDRFILE.

Warning

Closing an already closed file will lead to a crash with a double-free pointer error.

5.17.2.3. XTC functions

The XTC trajectory format is a lossy compression format that only stores coordinates. Compression level is determined by the precision argument to the write_xtc() function. Coordinates (Gromacs uses nm natively) are multiplied by precision and truncated to the integer part. A typical value is 1000.0, which gives an accuracy of 1/100 of an Angstroem.

The advantage of XTC over TRR is its significantly reduced size.

MDAnalysis.coordinates.xdrfile.libxdrfile2.read_xtc_natoms(fn) → natoms[source]

Read the number of atoms natoms from a xtc file fn.

Arguments:
fn

file name of an xtc file

Raises:

IOError if the supplied filed is not a XTC or if it is not readable.

MDAnalysis.coordinates.xdrfile.libxdrfile2.read_xtc_numframes(fn) -> (numframes, offsets)[source]

Read through the whole trajectory headers to obtain the total number of frames. The process is speeded up by reading frame headers for the amount of data in the frame, and then skipping directly to the next header. An array of frame offsets is also returned, which can later be used to seek direcly to arbitrary frames in the trajectory.

Arguments:
fn

file name of an xtc file

Returns:

The function returns a tuple containing numframes

an int with the total frame count in the trajectory

offsets

a numpy array of int64 recording the starting byte offset of each frame

Raises:

IOError if the supplied filed is not a XTC or if it is not readable.

MDAnalysis.coordinates.xdrfile.libxdrfile2.read_xtc(XDRFILE, box, x) -> (status, step, time, precision)[source]

Read the next frame from the opened xtc trajectory into x.

Arguments:
XDRFILE

open XDRFILE object

box

pre-allocated numpy array((DIM,DIM),dtype=numpy.float32) which is filled with the unit cell box vectors

x

pre-allocated numpy array((natoms, DIM),dtype=numpy.float32) which is updated with the coordinates from the frame

Returns:

The function returns a tuple containing status

integer status (0 = exdrOK), see Status symbols for other values)

step

simulation step

time

simulation time in ps

precision

precision of the lossy xtc format (typically 1000.0)

MDAnalysis.coordinates.xdrfile.libxdrfile2.write_xtc(XDRFILE, step, time, box, x, prec) → status[source]

Write the next frame x to the opened xtc trajectory.

Arguments:
XDRFILE

open XDRFILE object (writable)

step

simulation step

time

time step in ps

box

numpy array((DIM,DIM),dtype=numpy.float32) which contains the unit cell box vectors

x

numpy array((natoms, DIM),dtype=nump.float32) which contains the coordinates from the frame

precision

precision of the lossy xtc format (typically 1000.0)

Returns:

status, integer status (0 = OK), see the libxdrfile2.exdr* constants under Status symbols for other values)

5.17.2.4. TRR functions

TRR is the Gromacs native full-feature trajectory storage format. It can contain position coordinates, velocities and forces, and the lambda value for free energy perturbation calculations. Velocities and forces are optional in the sense that they can be all zero.

MDAnalysis.coordinates.xdrfile.libxdrfile2.read_trr_natoms(fn) → natoms[source]

Read the number of atoms natoms from a trr file fn.

Arguments:
fn

file name of a trr file

Raises:

IOError if the supplied filed is not a TRR or if it is not readable.

MDAnalysis.coordinates.xdrfile.libxdrfile2.read_trr_numframes(fn) -> (numframes, offsets)[source]

Read through the whole trajectory headers to obtain the total number of frames. The process is speeded up by reading frame headers for the amount of data in the frame, and then skipping directly to the next header. An array of frame offsets is also returned, which can later be used to seek direcly to arbitrary frames in the trajectory.

Arguments:
fn

file name of an xtc file

Returns:

The function returns a tuple containing numframes

an int with the total frame count in the trajectory

offsets

a numpy array of int64 recording the starting byte offset of each frame

Raises:

IOError if the supplied filed is not a TRR or if it is not readable.

MDAnalysis.coordinates.xdrfile.libxdrfile2.read_trr(XDRFILE, box, x, v, f) -> (status, step, time, lambda)[source]

Read the next frame from the opened trr trajectory into x, v, and f.

Arguments:
XDRFILE

open XDRFILE object

box

pre-allocated numpy array((DIM,DIM),dtype=numpy.float32) which is filled with the unit cell box vectors

x

pre-allocated numpy array((natoms, DIM),dtype=nump.float32) which is updated with the coordinates from the frame

v

pre-allocated numpy array((natoms, DIM),dtype=nump.float32) which is updated with the velocities from the frame

f

pre-allocated numpy array((natoms, DIM),dtype=nump.float32) which is updated with the forces from the frame

Returns:

The function returns a tuple containing status

integer status (0 = exdrOK), see the libxdrfile2.exdr* constants under Status symbols for other values)

step

simulation step

time

simulation time in ps

lambda

current lambda value (only interesting for free energy perturbation)

has_x

boolean indicating whether coordinates were read from the TRR

has_v

boolean indicating whether velocities were read from the TRR

has_f

boolean indicating whether forces were read from the TRR

MDAnalysis.coordinates.xdrfile.libxdrfile2.write_trr(XDRFILE, step, time, lambda, box, x, v, f) → status[source]

Write the next frame to the opened trr trajectory.

Arguments:
XDRFILE

open XDRFILE object (writable)

step

simulation step

time

time step in ps

lambda

free energy lambda value (typically 0.0)

box

numpy array((DIM,DIM),dtype=numpy.float32) which contains the unit cell box vectors

x

numpy array((natoms, DIM),dtype=nump.float32) which contains the coordinates from the frame

v

numpy array((natoms, DIM),dtype=nump.float32) which contains the velocities from the frame

f

numpy array((natoms, DIM),dtype=nump.float32) which contains the forces from the frame

Changed in version 0.8.0: either one of x, v, or f can now be set as a natom,0-DIM numpy array((natom, 0),dtype=nump.float32). This will cause the corresponding property to be skipped when writing to file.

Returns:status, integer status (0 = OK), see the libxdrfile2.exdr* constants under Status symbols for other values)