Skip to content

Commit 7583a2e

Browse files
committed
ENH: PR 4873 revisions
* Expand `TPRReader()` support to include velocity handling, and add tests/functionality for an additional tpx version (`133`). [ci skip] [skip azp]
1 parent 4af047d commit 7583a2e

File tree

3 files changed

+38
-8
lines changed

3 files changed

+38
-8
lines changed

package/MDAnalysis/coordinates/TPR.py

Lines changed: 5 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -18,7 +18,7 @@ class TPRReader(base.SingleFrameReaderBase):
1818
# or perhaps combine the topology and coordinate reading
1919
# with some inheritance shenanigans?
2020
format = "TPR"
21-
units = {"length": "nm"}
21+
units = {"length": "nm", "velocity": "nm/ps"}
2222
_Timestep = Timestep
2323

2424
def _read_first_frame(self):
@@ -71,3 +71,7 @@ def _read_first_frame(self):
7171
if th.bX:
7272
self.ts._pos = np.asarray(tpr_utils.ndo_rvec(data, th.natoms),
7373
dtype=np.float32)
74+
if th.bV:
75+
self.ts._velocities = np.asarray(tpr_utils.ndo_rvec(data, th.natoms),
76+
dtype=np.float32)
77+
self.ts.has_velocities = True

package/MDAnalysis/topology/tpr/utils.py

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -433,7 +433,7 @@ def do_mtop(data, fver, tpr_resid_from_one=False):
433433
# src/gromacs/fileio/tpxio.cpp
434434
# TODO: expand tpx version support for striding to
435435
# the coordinates
436-
if fver == 134:
436+
if fver >= 133:
437437
# TODO: the following value is important, and not sure
438438
# how to access programmatically yet...
439439
# from GMX source code:
Lines changed: 32 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -1,34 +1,60 @@
11
from MDAnalysisTests.datafiles import (TPR2024_4_bonded,
22
TPR_EXTRA_2024_4,
3-
TPR2024_4)
3+
TPR2024_4,
4+
TPR2024)
45
import MDAnalysis as mda
56

67

78
import pytest
9+
import numpy as np
810
from numpy.testing import assert_allclose, assert_equal
911

1012

11-
@pytest.mark.parametrize("tpr_file, exp_first_atom, exp_last_atom, exp_shape", [
12-
(TPR2024_4_bonded,
13+
@pytest.mark.parametrize("tpr_file, exp_first_atom, exp_last_atom, exp_shape, exp_vel_first_atom, exp_vel_last_atom", [
14+
(TPR2024_4_bonded, # tpx 134
1315
[4.446, 4.659, 2.384],
1416
[4.446, 4.659, 2.384],
1517
(14, 3),
18+
np.zeros(3),
19+
np.zeros(3),
1620
),
1721
# same coordinates, different shape
18-
(TPR_EXTRA_2024_4,
22+
(TPR_EXTRA_2024_4, # tpx 134
1923
[4.446, 4.659, 2.384],
2024
[4.446, 4.659, 2.384],
2125
(18, 3),
26+
np.zeros(3),
27+
np.zeros(3),
2228
),
2329
# different coordinates and different shape
24-
(TPR2024_4,
30+
(TPR2024_4, # tpx 134
2531
[3.25000e-01, 1.00400e+00, 1.03800e+00],
2632
[-2.56000e-01, 1.37300e+00, 3.59800e+00],
2733
(2263, 3),
34+
np.zeros(3),
35+
np.zeros(3),
36+
),
37+
(TPR2024, # tpx 133
38+
[3.25000e-01, 1.00400e+00, 1.03800e+00],
39+
[-2.56000e-01, 1.37300e+00, 3.59800e+00],
40+
(2263, 3),
41+
np.zeros(3),
42+
np.zeros(3),
2843
),
2944
])
30-
def test_basic_read_tpr(tpr_file, exp_first_atom, exp_last_atom, exp_shape):
45+
def test_basic_read_tpr(tpr_file,
46+
exp_first_atom,
47+
exp_last_atom,
48+
exp_shape,
49+
exp_vel_first_atom,
50+
exp_vel_last_atom):
51+
# verify basic ability to read positions and
52+
# velocities from GMX .tpr files
53+
# expected values are from gmx dump
3154
u = mda.Universe(tpr_file)
3255
assert_allclose(u.atoms.positions[0, ...], exp_first_atom)
3356
assert_allclose(u.atoms.positions[-1, ...], exp_last_atom)
3457
assert_equal(u.atoms.positions.shape, exp_shape)
58+
assert_allclose(u.atoms.velocities[0, ...], exp_vel_first_atom)
59+
assert_allclose(u.atoms.velocities[-1, ...], exp_vel_last_atom)
60+
assert_equal(u.atoms.velocities.shape, exp_shape)

0 commit comments

Comments
 (0)