File: test_fileio_low_level.py

package info (click to toggle)
gromacs 2026~rc-1
  • links: PTS, VCS
  • area: main
  • in suites: experimental
  • size: 274,216 kB
  • sloc: xml: 3,831,143; cpp: 686,111; ansic: 75,300; python: 21,171; sh: 3,553; perl: 2,246; yacc: 644; fortran: 397; lisp: 265; makefile: 174; lex: 125; awk: 68; csh: 39
file content (117 lines) | stat: -rw-r--r-- 4,689 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
#
# This file is part of the GROMACS molecular simulation package.
#
# Copyright 2021- The GROMACS Authors
# and the project initiators Erik Lindahl, Berk Hess and David van der Spoel.
# Consult the AUTHORS/COPYING files and https://www.gromacs.org for details.
#
# GROMACS is free software; you can redistribute it and/or
# modify it under the terms of the GNU Lesser General Public License
# as published by the Free Software Foundation; either version 2.1
# of the License, or (at your option) any later version.
#
# GROMACS is distributed in the hope that it will be useful,
# but WITHOUT ANY WARRANTY; without even the implied warranty of
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
# Lesser General Public License for more details.
#
# You should have received a copy of the GNU Lesser General Public
# License along with GROMACS; if not, see
# https://www.gnu.org/licenses, or write to the Free Software Foundation,
# Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
#
# If you want to redistribute modifications to GROMACS, please
# consider that scientific software is very special. Version
# control is crucial - bugs must be traceable. We will be happy to
# consider code for inclusion in the official distribution, but
# derived work must not be called official GROMACS. Details are found
# in the README & COPYING files - if they are missing, get the
# official version at https://www.gromacs.org.
#
# To help us fund GROMACS development, we humbly ask that you cite
# the research papers on the package. Check out https://www.gromacs.org.

"""Test gmxapi._gmxapi file I/O routines."""
import os
import tempfile

import pytest
import gmxapi
import gmxapi.simulation.fileio
from gmxapi import _gmxapi


@pytest.mark.usefixtures("cleandir")
def test_core_read_tpr(spc_water_box):
    tpr_filehandle: _gmxapi.TprFile = _gmxapi.read_tprfile(os.fsencode(spc_water_box))
    parameters: _gmxapi.SimulationParameters = tpr_filehandle.params()
    assert "nsteps" in parameters.extract()
    assert "foo" not in parameters.extract()
    assert parameters.extract()["nsteps"] == 2


@pytest.mark.usefixtures("cleandir")
def test_core_rewrite_tprfile(spc_water_box):
    """Test _gmxapi.rewrite_tprfile for update of end_time.

    Set a new end time that is 5000 steps later than the original. Read dt
    from file to avoid floating point round-off errors.

    Transitively test gmx.fileio.read_tpr()
    """
    tpr_filename = spc_water_box
    additional_steps = 5000
    tpr_filehandle: _gmxapi.TprFile = _gmxapi.read_tprfile(os.fsencode(spc_water_box))
    params = tpr_filehandle.params().extract()
    dt = params["dt"]
    nsteps = params["nsteps"]
    init_step = params["init-step"]
    initial_endtime = (init_step + nsteps) * dt
    new_endtime = initial_endtime + additional_steps * dt
    _, temp_filename = tempfile.mkstemp(suffix=".tpr")
    _gmxapi.rewrite_tprfile(
        source=tpr_filename, destination=temp_filename, end_time=new_endtime
    )
    tprfile = gmxapi.simulation.fileio.TprFile(temp_filename, "r")
    with tprfile as fh:
        params = gmxapi.simulation.fileio.read_tpr(fh).parameters.extract()
        dt = params["dt"]
        nsteps = params["nsteps"]
        init_step = params["init-step"]
        assert (init_step + nsteps) * dt == new_endtime

    os.unlink(temp_filename)


@pytest.mark.usefixtures("cleandir")
def test_core_read_and_write_tpr_file(spc_water_box):
    """Test gmx.fileio.write_tpr_file() using gmx.core API."""
    additional_steps = 5000

    tpr_filehandle: _gmxapi.TprFile = _gmxapi.read_tprfile(os.fsencode(spc_water_box))
    sim_input: _gmxapi.SimulationParameters = tpr_filehandle.params()
    params: dict = sim_input.extract()
    nsteps = params["nsteps"]
    init_step = params["init-step"]

    # Choose a new nsteps to check integer parameter setting.
    new_nsteps = init_step + additional_steps
    # Choose a new dt to check floating point parameter setting
    new_dt = params["dt"] * 2.0

    sim_input.set("nsteps", new_nsteps)
    sim_input.set("dt", new_dt)

    _, temp_filename = tempfile.mkstemp(suffix=".tpr")
    _gmxapi.write_tprfile(temp_filename, sim_input)

    tprfile = gmxapi.simulation.fileio.TprFile(temp_filename, "r")
    with tprfile as fh:
        params = gmxapi.simulation.fileio.read_tpr(fh).parameters.extract()
        # Note that we have chosen an exactly representable dt for spc_water_box.
        # Otherwise, we would have to use pytest.approx with a suitable tolerance.
        assert params["dt"] == new_dt
        assert params["nsteps"] != nsteps
        assert params["nsteps"] == new_nsteps

    os.unlink(temp_filename)