4f7301e2dfa0bc359d921bbd558d5e5077e63df4
[alexxy/gromacs.git] / python_packaging / src / gmxapi / export_tprfile.cpp
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 /*! \file
36  * \brief Exports TPR I/O tools during Python module initialization.
37  *
38  * Provides _gmxapi.SimulationParameters and _gmxapi.TprFile classes, as well
39  * as module functions read_tprfile, write_tprfile, copy_tprfile, and rewrite_tprfile.
40  *
41  * TprFile is a Python object that holds a gmxapicompat::TprReadHandle.
42  *
43  * SimulationParameters is the Python type for data sources providing the
44  * simulation parameters aspect of input to simulation operations.
45  *
46  * \author M. Eric Irrgang <ericirrgang@gmail.com>
47  * \ingroup module_python
48  */
49
50 #include "gmxapi/exceptions.h"
51
52 #include "gmxapi/gmxapicompat.h"
53 #include "gmxapi/compat/mdparams.h"
54 #include "gmxapi/compat/tpr.h"
55
56 #include "module.h"
57
58 using gmxapi::GmxapiType;
59
60 namespace gmxpy
61 {
62
63
64 void detail::export_tprfile(pybind11::module& module)
65 {
66     namespace py = pybind11;
67     using gmxapicompat::GmxMdParams;
68     using gmxapicompat::readTprFile;
69     using gmxapicompat::TprReadHandle;
70
71     py::class_<GmxMdParams> mdparams(module, "SimulationParameters");
72     // We don't want Python users to create invalid params objects, so don't
73     // export a constructor until we can default initialize a valid one.
74     //    mdparams.def(py::init());
75     mdparams.def("extract",
76                  [](const GmxMdParams& self) {
77                      py::dict dictionary;
78                      for (const auto& key : gmxapicompat::keys(self))
79                      {
80                          try
81                          {
82                              // TODO: More complete typing and dispatching.
83                              // This only handles the two types described in the initial implementation.
84                              // Less trivial types (strings, maps, arrays) warrant additional
85                              // design discussion before being exposed through an interface
86                              // like this one.
87                              // Also reference https://gitlab.com/gromacs/gromacs/-/issues/2993
88
89                              // We can use templates and/or tag dispatch in a more complete
90                              // future implementation.
91                              const auto& paramType = gmxapicompat::mdParamToType(key);
92                              if (paramType == GmxapiType::FLOAT64)
93                              {
94                                  dictionary[key.c_str()] = extractParam(self, key, double());
95                              }
96                              else if (paramType == GmxapiType::INT64)
97                              {
98                                  dictionary[key.c_str()] = extractParam(self, key, int64_t());
99                              }
100                          }
101                          catch (const gmxapicompat::ValueError& e)
102                          {
103                              throw gmxapi::ProtocolError(std::string("Unknown parameter: ") + key);
104                          }
105                      }
106                      return dictionary;
107                  },
108                  "Get a dictionary of the parameters.");
109
110     // Overload a setter for each known type and None
111     mdparams.def("set",
112                  [](GmxMdParams* self, const std::string& key, int64_t value) {
113                      gmxapicompat::setParam(self, key, value);
114                  },
115                  py::arg("key").none(false),
116                  py::arg("value").none(false),
117                  "Use a dictionary to update simulation parameters.");
118     mdparams.def("set",
119                  [](GmxMdParams* self, const std::string& key, double value) {
120                      gmxapicompat::setParam(self, key, value);
121                  },
122                  py::arg("key").none(false),
123                  py::arg("value").none(false),
124                  "Use a dictionary to update simulation parameters.");
125     mdparams.def("set",
126                  [](GmxMdParams* self, const std::string& key, py::none) {
127                      // unsetParam(self, key);
128                  },
129                  py::arg("key").none(false),
130                  py::arg("value"),
131                  "Use a dictionary to update simulation parameters.");
132
133
134     py::class_<TprReadHandle> tprfile(module, "TprFile");
135     tprfile.def("params", [](const TprReadHandle& self) {
136         auto params = gmxapicompat::getMdParams(self);
137         return params;
138     });
139
140     module.def("read_tprfile",
141                &readTprFile,
142                py::arg("filename"),
143                "Get a handle to a TPR file resource for a given file name.");
144
145     module.def("write_tprfile",
146                [](std::string filename, const GmxMdParams& parameterObject) {
147                    auto tprReadHandle = gmxapicompat::getSourceFileHandle(parameterObject);
148                    auto params        = gmxapicompat::getMdParams(tprReadHandle);
149                    auto structure     = gmxapicompat::getStructureSource(tprReadHandle);
150                    auto state         = gmxapicompat::getSimulationState(tprReadHandle);
151                    auto topology      = gmxapicompat::getTopologySource(tprReadHandle);
152                    gmxapicompat::writeTprFile(filename, *params, *structure, *state, *topology);
153                },
154                py::arg("filename").none(false),
155                py::arg("parameters"),
156                "Write a new TPR file with the provided data.");
157
158     module.def("copy_tprfile",
159                [](const gmxapicompat::TprReadHandle& input, std::string outFile) {
160                    return gmxapicompat::copy_tprfile(input, outFile);
161                },
162                py::arg("source"),
163                py::arg("destination"),
164                "Copy a TPR file from ``source`` to ``destination``.");
165
166     module.def("rewrite_tprfile",
167                [](std::string input, std::string output, double end_time) {
168                    return gmxapicompat::rewrite_tprfile(input, output, end_time);
169                },
170                py::arg("source"),
171                py::arg("destination"),
172                py::arg("end_time"),
173                "Copy a TPR file from ``source`` to ``destination``, replacing `nsteps` with "
174                "``end_time``.");
175 }
176
177 } // end namespace gmxpy