Skip to content

Commit d67bccb

Browse files
committed
ENH: add TPR position read support
* Add support for reading positions from GMX `.tpr` files at `tpx` version `134`.
1 parent 5656fe4 commit d67bccb

File tree

4 files changed

+150
-0
lines changed

4 files changed

+150
-0
lines changed

package/MDAnalysis/coordinates/TPR.py

Lines changed: 73 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,73 @@
1+
# -*- Mode: python; tab-width: 4; indent-tabs-mode:nil; coding:utf-8 -*-
2+
# vim: tabstop=4 expandtab shiftwidth=4 softtabstop=4
3+
4+
from . import base
5+
from ..lib import util
6+
from .timestep import Timestep
7+
import MDAnalysis.topology.tpr.utils as tpr_utils
8+
import MDAnalysis.topology.tpr.setting as S
9+
10+
11+
import numpy as np
12+
13+
14+
class TPRReader(base.SingleFrameReaderBase):
15+
# TODO: reduce duplication with `TPRparser`;
16+
# we could share some state for the position
17+
# in the binary file to avoid re-reading topology
18+
# or perhaps combine the topology and coordinate reading
19+
# with some inheritance shenanigans?
20+
format = "TPR"
21+
units = {"length": "nm"}
22+
_Timestep = Timestep
23+
24+
def _read_first_frame(self):
25+
# Read header/move over topology
26+
# TODO: reduce duplication with TPRparser perhaps...
27+
with util.openany(self.filename, mode='rb') as infile:
28+
tprf = infile.read()
29+
data = tpr_utils.TPXUnpacker(tprf)
30+
try:
31+
th = tpr_utils.read_tpxheader(data) # tpxheader
32+
except (EOFError, ValueError):
33+
msg = f"{self.filename}: Invalid tpr file or cannot be recognized"
34+
logger.critical(msg)
35+
raise IOError(msg)
36+
37+
self.ts = ts = self._Timestep(th.natoms, **self._ts_kwargs)
38+
self.n_atoms = th.natoms
39+
40+
# Starting with gromacs 2020 (tpx version 119), the body of the file
41+
# is encoded differently. We change the unpacker accordingly.
42+
if th.fver >= S.tpxv_AddSizeField and th.fgen >= 27:
43+
actual_body_size = len(data.get_buffer()) - data.get_position()
44+
if actual_body_size == 4 * th.sizeOfTprBody:
45+
# See issue #2428.
46+
msg = (
47+
"TPR files produced with beta versions of gromacs 2020 "
48+
"are not supported."
49+
)
50+
logger.critical(msg)
51+
raise IOError(msg)
52+
data = tpr_utils.TPXUnpacker2020.from_unpacker(data)
53+
54+
state_ngtc = th.ngtc # done init_state() in src/gmxlib/tpxio.c
55+
if th.bBox:
56+
tpr_utils.extract_box_info(data, th.fver)
57+
58+
if state_ngtc > 0:
59+
if th.fver < 69: # redundancy due to different versions
60+
tpr_utils.ndo_real(data, state_ngtc)
61+
tpr_utils.ndo_real(data, state_ngtc) # relevant to Berendsen tcoupl_lambda
62+
63+
if th.bTop:
64+
tpr_top = tpr_utils.do_mtop(data, th.fver,
65+
tpr_resid_from_one=True)
66+
else:
67+
msg = f"{self.filename}: No topology found in tpr file"
68+
logger.critical(msg)
69+
raise IOError(msg)
70+
71+
if th.bX:
72+
self.ts._pos = np.asarray(tpr_utils.ndo_rvec(data, th.natoms),
73+
dtype=np.float32)

package/MDAnalysis/coordinates/__init__.py

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -776,6 +776,7 @@ class can choose an appropriate reader automatically.
776776
from . import PDB
777777
from . import PDBQT
778778
from . import PQR
779+
from . import TPR
779780
from . import TRC
780781
from . import TRJ
781782
from . import TRR

package/MDAnalysis/topology/tpr/utils.py

Lines changed: 42 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -423,6 +423,48 @@ def do_mtop(data, fver, tpr_resid_from_one=False):
423423
elements = Elements(np.array(elements, dtype=object))
424424
top.add_TopologyAttr(elements)
425425

426+
# NOTE: the tpr striding code below serves the
427+
# purpose of placing us in a suitable "seek" position
428+
# in the binary file such that we are well placed
429+
# for reading coords and velocities if they are present
430+
# the order of operations is based on an analysis of
431+
# the C++ code in do_mtop() function in the GROMACS
432+
# source at:
433+
# src/gromacs/fileio/tpxio.cpp
434+
# TODO: expand tpx version support for striding to
435+
# the coordinates
436+
if fver == 134:
437+
# TODO: the following value is important, and not sure
438+
# how to access programmatically yet...
439+
# from GMX source code:
440+
# api/legacy/include/gromacs/topology/topology_enums.h
441+
# worst case scenario we hard code it based on
442+
# tpx/GMX version?
443+
SimulationAtomGroupType_size = 10
444+
n_atoms = data.unpack_int()
445+
interm = data.unpack_uchar()
446+
ngrid = data.unpack_int()
447+
grid_spacing = data.unpack_int()
448+
n_elements = grid_spacing ** 2
449+
for i in range(ngrid):
450+
for j in range(n_elements):
451+
ndo_real(data, 4)
452+
for i in range(SimulationAtomGroupType_size):
453+
group_size = data.unpack_int()
454+
ndo_int(data, group_size)
455+
n_group_names = data.unpack_int()
456+
for i in range(n_group_names):
457+
data.unpack_int()
458+
for i in range(SimulationAtomGroupType_size):
459+
n_grp_numbers = data.unpack_int()
460+
if n_grp_numbers != 0:
461+
for i in range(n_grp_numbers):
462+
data.unpack_uchar()
463+
im_excl_grp_size = data.unpack_int()
464+
ndo_int(data, im_excl_grp_size)
465+
# TODO: why is this needed?
466+
data.unpack_int()
467+
426468
return top
427469

428470

Lines changed: 34 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,34 @@
1+
from MDAnalysisTests.datafiles import (TPR2024_4_bonded,
2+
TPR_EXTRA_2024_4,
3+
TPR2024_4)
4+
import MDAnalysis as mda
5+
6+
7+
import pytest
8+
from numpy.testing import assert_allclose, assert_equal
9+
10+
11+
@pytest.mark.parametrize("tpr_file, exp_first_atom, exp_last_atom, exp_shape", [
12+
(TPR2024_4_bonded,
13+
[4.446, 4.659, 2.384],
14+
[4.446, 4.659, 2.384],
15+
(14, 3),
16+
),
17+
# same coordinates, different shape
18+
(TPR_EXTRA_2024_4,
19+
[4.446, 4.659, 2.384],
20+
[4.446, 4.659, 2.384],
21+
(18, 3),
22+
),
23+
# different coordinates and different shape
24+
(TPR2024_4,
25+
[3.25000e-01, 1.00400e+00, 1.03800e+00],
26+
[-2.56000e-01, 1.37300e+00, 3.59800e+00],
27+
(2263, 3),
28+
),
29+
])
30+
def test_basic_read_tpr(tpr_file, exp_first_atom, exp_last_atom, exp_shape):
31+
u = mda.Universe(tpr_file)
32+
assert_allclose(u.atoms.positions[0, ...], exp_first_atom)
33+
assert_allclose(u.atoms.positions[-1, ...], exp_last_atom)
34+
assert_equal(u.atoms.positions.shape, exp_shape)

0 commit comments

Comments
 (0)