Restore gmxapi._gmxapi.add_mdmodule() Python functionality.
[alexxy/gromacs.git] / python_packaging / sample_restraint / tests / test_plugin.py
1 # The myplugin module must be locatable by Python.
2 # If you configured CMake in the build directory ``/path/to/repo/build`` then,
3 # assuming you are in ``/path/to/repo``, run the tests with something like
4 #     PYTHONPATH=./cmake-build-debug/src/pythonmodule mpiexec -n 2 python -m mpi4py -m pytest tests/
5
6 # This test is not currently run automatically in any way. Build the module, point your PYTHONPATH at it,
7 # and run pytest in the tests directory.
8
9 import logging
10 import os
11
12 try:
13     import mpi4py.MPI as _MPI
14 except (ImportError, ModuleNotFoundError):
15     _MPI = None
16
17 import gmxapi as gmx
18 from gmxapi.simulation.context import Context
19 from gmxapi.simulation.workflow import WorkElement, from_tpr
20 from gmxapi import version as gmx_version
21 import pytest
22
23 # create console handler
24 ch = logging.StreamHandler()
25 ch.setLevel(logging.DEBUG)
26 # create formatter and add it to the handler
27 formatter = logging.Formatter('%(asctime)s:%(name)s:%(levelname)s: %(message)s')
28 ch.setFormatter(formatter)
29 # add the handlers to the logger
30 logging.getLogger().addHandler(ch)
31
32 logger = logging.getLogger()
33
34
35 def test_import():
36     # Suppress inspection warning outside of testing context.
37     # noinspection PyUnresolvedReferences
38     import myplugin
39     assert myplugin
40
41
42 @pytest.mark.usefixtures("cleandir")
43 def test_binding_protocol(spc_water_box, mdrun_kwargs):
44     """Test that gmxapi successfully attaches MD plugins."""
45     import myplugin
46
47     if _MPI is not None:
48         _size = _MPI.COMM_WORLD.Get_size()
49         _rank = _MPI.COMM_WORLD.Get_rank()
50     else:
51         _size = 1
52         _rank = 0
53
54     tpr_filename = spc_water_box
55     logger.info("Testing plugin potential with input file {}".format(os.path.abspath(tpr_filename)))
56
57     assert gmx.version.api_is_at_least(0, 2, 1)
58     md = from_tpr([tpr_filename] * _size, append_output=False, **mdrun_kwargs)
59
60     potential = WorkElement(namespace="myplugin",
61                             operation="null_restraint",
62                             params={'sites': [1, 4]})
63     potential.name = "null restraint"
64     md.add_dependency(potential)
65
66     context = Context(md)
67
68     with context as session:
69         session.run()
70
71     # See also #3038, #3145, #4079
72     assert isinstance(context.potentials, list)
73     assert len(context.potentials) > 0
74     for restraint in context.potentials:
75         if isinstance(restraint, myplugin.NullRestraint):
76             assert restraint.count() > 1
77
78
79 @pytest.mark.usefixtures("cleandir")
80 def test_ensemble_potential_nompi(spc_water_box, mdrun_kwargs):
81     """Test ensemble potential without an ensemble.
82     """
83     tpr_filename = spc_water_box
84     logger.info("Testing plugin potential with input file {}".format(os.path.abspath(tpr_filename)))
85
86     assert gmx.version.api_is_at_least(0, 0, 5)
87     md = from_tpr([tpr_filename], append_output=False, **mdrun_kwargs)
88
89     # Create a WorkElement for the potential
90     params = {'sites': [1, 4],
91               'nbins': 10,
92               'binWidth': 0.1,
93               'min_dist': 0.,
94               'max_dist': 10.,
95               'experimental': [1.] * 10,
96               'nsamples': 1,
97               'sample_period': 0.001,
98               'nwindows': 4,
99               'k': 10000.,
100               'sigma': 1.}
101     potential = WorkElement(namespace="myplugin",
102                             operation="ensemble_restraint",
103                             params=params)
104     # Note that we could flexibly capture accessor methods as workflow elements, too. Maybe we can
105     # hide the extra Python bindings by letting myplugin.HarmonicRestraint automatically convert
106     # to a WorkElement when add_dependency is called on it.
107     potential.name = "ensemble_restraint"
108     md.add_dependency(potential)
109
110     context = Context(md)
111
112     with context as session:
113         session.run()
114
115
116 @pytest.mark.withmpi_only
117 @pytest.mark.usefixtures("cleandir")
118 def test_ensemble_potential_withmpi(spc_water_box, mdrun_kwargs):
119     tpr_filename = spc_water_box
120
121     logger.info("Testing plugin potential with input file {}".format(os.path.abspath(tpr_filename)))
122
123     assert gmx_version.api_is_at_least(0, 0, 5)
124     md = from_tpr([tpr_filename, tpr_filename], append_output=False, **mdrun_kwargs)
125
126     # Create a WorkElement for the potential
127     params = {'sites': [1, 4],
128               'nbins': 10,
129               'binWidth': 0.1,
130               'min_dist': 0.,
131               'max_dist': 10.,
132               'experimental': [0.5] * 10,
133               'nsamples': 1,
134               'sample_period': 0.001,
135               'nwindows': 4,
136               'k': 10000.,
137               'sigma': 1.}
138
139     potential = WorkElement(namespace="myplugin",
140                             operation="ensemble_restraint",
141                             params=params)
142     # Note that we could flexibly capture accessor methods as workflow elements, too. Maybe we can
143     # hide the extra Python bindings by letting myplugin.HarmonicRestraint automatically convert
144     # to a WorkElement when add_dependency is called on it.
145     potential.name = "ensemble_restraint"
146     md.add_dependency(potential)
147
148     context = Context(md)
149     with context as session:
150         session.run()