Restrict gmxapi pytest tests to `-nt 2`
[alexxy/gromacs.git] / python_packaging / src / test / test_mdrun.py
1 #
2 # This file is part of the GROMACS molecular simulation package.
3 #
4 # Copyright (c) 2019,2020, 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 gromacs.mdrun operation.
36
37 Factory produces deferred execution operation.
38
39 TODO: Factory expects input for conformation, topology, simulation parameters, simulation state.
40
41 TODO: Factory accepts additional keyword input to indicate binding
42  to the "potential" interface.
43 """
44
45 import logging
46 import os
47 import pytest
48
49 import gmxapi as gmx
50
51 # TODO: (#3573) Normalize the handling of run-time arguments.
52 from gmxapi.simulation.mdrun import ResourceManager as _ResourceManager
53 # Note that *threads* argument causes errors for MPI-enabled GROMACS.
54 # Ref #3563 and #3573
55 _ResourceManager.mdrun_kwargs = {'threads': 2}
56
57 # Configure the `logging` module before proceeding any further.
58 gmx.logger.setLevel(logging.WARNING)
59
60 try:
61     from mpi4py import MPI
62     rank_number = MPI.COMM_WORLD.Get_rank()
63 except ImportError:
64     rank_number = 0
65     rank_tag = ''
66     MPI = None
67 else:
68     rank_tag = 'rank{}:'.format(rank_number)
69
70 # Use this formatter to improve the caplog log records.
71 formatter = logging.Formatter(rank_tag + '%(name)s:%(levelname)s: %(message)s')
72
73 # For additional console logging, create and attach a stream handler.
74 # For example:
75 #    ch = logging.StreamHandler()
76 #    ch.setFormatter(formatter)
77 #    logging.getLogger().addHandler(ch)
78
79
80 @pytest.mark.usefixtures('cleandir')
81 def test_run_from_tpr(spc_water_box):
82     assert os.path.exists(spc_water_box)
83
84     md = gmx.mdrun(spc_water_box)
85     md.run()
86     # TODO: better handling of output on unused MPI ranks.
87
88
89 @pytest.mark.withmpi_only
90 @pytest.mark.usefixtures('cleandir')
91 def test_run_trivial_ensemble(spc_water_box, caplog):
92     from mpi4py import MPI
93     current_rank = MPI.COMM_WORLD.Get_rank()
94     with caplog.at_level(logging.DEBUG):
95         caplog.handler.setFormatter(formatter)
96         with caplog.at_level(logging.WARNING, 'gmxapi'), \
97                 caplog.at_level(logging.DEBUG, 'gmxapi.mdrun'), \
98                 caplog.at_level(logging.DEBUG, 'gmxapi.modify_input'), \
99                 caplog.at_level(logging.DEBUG, 'gmxapi.read_tpr'), \
100                 caplog.at_level(logging.DEBUG, 'gmxapi.simulation'):
101
102             tpr_filename = spc_water_box
103             ensemble_width = 2
104             simulation_input = gmx.read_tpr([tpr_filename] * ensemble_width)
105             assert simulation_input.output.ensemble_width == ensemble_width
106             assert len(simulation_input.output._simulation_input.result()) == ensemble_width
107             md = gmx.mdrun(simulation_input)
108             assert md.output.ensemble_width == ensemble_width
109             md.run()
110
111             output_directory = md.output._work_dir.result()
112             logging.info('output_directory result: {}'.format(str(output_directory)))
113             assert len(output_directory) == 2
114
115             # Note that the 'cleandir' test fixture will clean up the output directory on
116             # other ranks, so only check the current rank. Generally, our behavior
117             # is undefined if the client removes the working directory while the job
118             # is in progress. We can consider adding some sort of synchronization at
119             # the end of the job if running in temporary directories becomes an
120             # important use case outside of testing.
121             assert output_directory[0] != output_directory[1]
122             assert os.path.exists(output_directory[current_rank])
123             assert os.path.exists(md.output.trajectory.result()[current_rank])
124
125
126 @pytest.mark.usefixtures('cleandir')
127 def test_run_from_read_tpr_op(spc_water_box, caplog):
128     with caplog.at_level(logging.DEBUG):
129         caplog.handler.setFormatter(formatter)
130         with caplog.at_level(logging.DEBUG, 'gmxapi'):
131             simulation_input = gmx.read_tpr(spc_water_box)
132             md = gmx.mdrun(input=simulation_input)
133
134             md.run()
135             if rank_number == 0:
136                 assert os.path.exists(md.output.trajectory.result())
137
138
139 @pytest.mark.usefixtures('cleandir')
140 def test_run_from_modify_input_op(spc_water_box, caplog):
141     with caplog.at_level(logging.DEBUG):
142
143         simulation_input = gmx.read_tpr(spc_water_box)
144         modified_input = gmx.modify_input(input=simulation_input, parameters={'nsteps': 4})
145         md = gmx.mdrun(input=modified_input)
146
147         md.run()