Clean up some gmxapi sources and tests.
[alexxy/gromacs.git] / python_packaging / src / test / test_fileio_low_level.py
1 #
2 # This file is part of the GROMACS molecular simulation package.
3 #
4 # Copyright (c) 2021, by the GROMACS development team, led by
5 # Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 # and including many others, as listed in the AUTHORS file in the
7 # top-level source directory and at http://www.gromacs.org.
8 #
9 # GROMACS is free software; you can redistribute it and/or
10 # modify it under the terms of the GNU Lesser General Public License
11 # as published by the Free Software Foundation; either version 2.1
12 # of the License, or (at your option) any later version.
13 #
14 # GROMACS is distributed in the hope that it will be useful,
15 # but WITHOUT ANY WARRANTY; without even the implied warranty of
16 # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17 # Lesser General Public License for more details.
18 #
19 # You should have received a copy of the GNU Lesser General Public
20 # License along with GROMACS; if not, see
21 # http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 # Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23 #
24 # If you want to redistribute modifications to GROMACS, please
25 # consider that scientific software is very special. Version
26 # control is crucial - bugs must be traceable. We will be happy to
27 # consider code for inclusion in the official distribution, but
28 # derived work must not be called official GROMACS. Details are found
29 # in the README & COPYING files - if they are missing, get the
30 # official version at http://www.gromacs.org.
31 #
32 # To help us fund GROMACS development, we humbly ask that you cite
33 # the research papers on the package. Check out http://www.gromacs.org.
34
35 """Test gmxapi._gmxapi file I/O routines."""
36 import os
37 import tempfile
38
39 import pytest
40 import gmxapi
41 import gmxapi.simulation.fileio
42 from gmxapi import _gmxapi
43
44
45 @pytest.mark.usefixtures('cleandir')
46 def test_core_read_tpr(spc_water_box):
47     tpr_filehandle: _gmxapi.TprFile = _gmxapi.read_tprfile(os.fsencode(spc_water_box))
48     parameters: _gmxapi.SimulationParameters = tpr_filehandle.params()
49     assert 'nsteps' in parameters.extract()
50     assert 'foo' not in parameters.extract()
51     assert parameters.extract()['nsteps'] == 2
52
53
54 @pytest.mark.usefixtures('cleandir')
55 def test_core_rewrite_tprfile(spc_water_box):
56     """Test _gmxapi.rewrite_tprfile for update of end_time.
57
58     Set a new end time that is 5000 steps later than the original. Read dt
59     from file to avoid floating point round-off errors.
60
61     Transitively test gmx.fileio.read_tpr()
62     """
63     tpr_filename = spc_water_box
64     additional_steps = 5000
65     tpr_filehandle: _gmxapi.TprFile = _gmxapi.read_tprfile(os.fsencode(spc_water_box))
66     params = tpr_filehandle.params().extract()
67     dt = params['dt']
68     nsteps = params['nsteps']
69     init_step = params['init-step']
70     initial_endtime = (init_step + nsteps) * dt
71     new_endtime = initial_endtime + additional_steps * dt
72     _, temp_filename = tempfile.mkstemp(suffix='.tpr')
73     _gmxapi.rewrite_tprfile(source=tpr_filename, destination=temp_filename, end_time=new_endtime)
74     tprfile = gmxapi.simulation.fileio.TprFile(temp_filename, 'r')
75     with tprfile as fh:
76         params = gmxapi.simulation.fileio.read_tpr(fh).parameters.extract()
77         dt = params['dt']
78         nsteps = params['nsteps']
79         init_step = params['init-step']
80         assert (init_step + nsteps) * dt == new_endtime
81
82     os.unlink(temp_filename)
83
84
85 @pytest.mark.usefixtures('cleandir')
86 def test_core_read_and_write_tpr_file(spc_water_box):
87     """Test gmx.fileio.write_tpr_file() using gmx.core API.
88     """
89     additional_steps = 5000
90
91     tpr_filehandle: _gmxapi.TprFile = _gmxapi.read_tprfile(os.fsencode(spc_water_box))
92     sim_input: _gmxapi.SimulationParameters = tpr_filehandle.params()
93     params: dict = sim_input.extract()
94     nsteps = params['nsteps']
95     init_step = params['init-step']
96
97     # Choose a new nsteps to check integer parameter setting.
98     new_nsteps = init_step + additional_steps
99     # Choose a new dt to check floating point parameter setting
100     new_dt = params['dt'] * 2.
101
102     sim_input.set('nsteps', new_nsteps)
103     sim_input.set('dt', new_dt)
104
105     _, temp_filename = tempfile.mkstemp(suffix='.tpr')
106     _gmxapi.write_tprfile(temp_filename, sim_input)
107
108     tprfile = gmxapi.simulation.fileio.TprFile(temp_filename, 'r')
109     with tprfile as fh:
110         params = gmxapi.simulation.fileio.read_tpr(fh).parameters.extract()
111         # Note that we have chosen an exactly representable dt for spc_water_box.
112         # Otherwise, we would have to use pytest.approx with a suitable tolerance.
113         assert params['dt'] == new_dt
114         assert params['nsteps'] != nsteps
115         assert params['nsteps'] == new_nsteps
116
117     os.unlink(temp_filename)