ba96040a8ff3b2752330a4c16e044308555be2c5
[alexxy/gromacs.git] / python_packaging / sample_restraint / tests / test_binding.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 import gmxapi as gmx
13 from gmxapi.simulation.context import ParallelArrayContext
14 from gmxapi.simulation.workflow import WorkElement, from_tpr
15 from gmxapi import version as gmx_version
16 import pytest
17 from tests.conftest import withmpi_only
18
19 logging.getLogger().setLevel(logging.DEBUG)
20 # create console handler
21 ch = logging.StreamHandler()
22 ch.setLevel(logging.DEBUG)
23 # create formatter and add it to the handler
24 formatter = logging.Formatter('%(asctime)s:%(name)s:%(levelname)s: %(message)s')
25 ch.setFormatter(formatter)
26 # add the handlers to the logger
27 logging.getLogger().addHandler(ch)
28
29 logger = logging.getLogger()
30
31
32 def test_import():
33     import myplugin
34     assert myplugin
35
36
37 @pytest.mark.usefixtures("cleandir")
38 def test_ensemble_potential_nompi(tpr_filename):
39     """Test ensemble potential without an ensemble.
40     """
41     print("Testing plugin potential with input file {}".format(os.path.abspath(tpr_filename)))
42
43     assert gmx.version.api_is_at_least(0, 0, 5)
44     md = from_tpr([tpr_filename], append_output=False)
45
46     # Create a WorkElement for the potential
47     params = {'sites': [1, 4],
48               'nbins': 10,
49               'binWidth': 0.1,
50               'min_dist': 0.,
51               'max_dist': 10.,
52               'experimental': [1.] * 10,
53               'nsamples': 1,
54               'sample_period': 0.001,
55               'nwindows': 4,
56               'k': 10000.,
57               'sigma': 1.}
58     potential = WorkElement(namespace="myplugin",
59                             operation="ensemble_restraint",
60                             params=params)
61     # Note that we could flexibly capture accessor methods as workflow elements, too. Maybe we can
62     # hide the extra Python bindings by letting myplugin.HarmonicRestraint automatically convert
63     # to a WorkElement when add_dependency is called on it.
64     potential.name = "ensemble_restraint"
65     md.add_dependency(potential)
66
67     context = ParallelArrayContext(md)
68
69     with context as session:
70         session.run()
71
72
73 @withmpi_only
74 @pytest.mark.usefixtures("cleandir")
75 def test_ensemble_potential_withmpi(tpr_filename):
76     import os
77
78     from mpi4py import MPI
79     rank = MPI.COMM_WORLD.Get_rank()
80
81     rank_dir = os.path.join(os.getcwd(), str(rank))
82     os.mkdir(rank_dir)
83
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, tpr_filename], append_output=False)
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': [0.5] * 10,
96               'nsamples': 1,
97               'sample_period': 0.001,
98               'nwindows': 4,
99               'k': 10000.,
100               'sigma': 1.}
101
102     potential = WorkElement(namespace="myplugin",
103                             operation="ensemble_restraint",
104                             params=params)
105     # Note that we could flexibly capture accessor methods as workflow elements, too. Maybe we can
106     # hide the extra Python bindings by letting myplugin.HarmonicRestraint automatically convert
107     # to a WorkElement when add_dependency is called on it.
108     potential.name = "ensemble_restraint"
109     md.add_dependency(potential)
110
111     context = ParallelArrayContext(md)
112     with context as session:
113         session.run()