Merge remote-tracking branch 'origin/release-2019'
[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, 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::TprReadHandle;
69     using gmxapicompat::readTprFile;
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                  {
78                      py::dict dictionary;
79                      for (const auto &key: gmxapicompat::keys(self))
80                      {
81                          try
82                          {
83                              // TODO: More complete typing and dispatching.
84                              // This only handles the two types described in the initial implementation.
85                              // Less trivial types (strings, maps, arrays) warrant additional
86                              // design discussion before being exposed through an interface
87                              // like this one.
88                              // Also reference https://redmine.gromacs.org/issues/2993
89
90                              // We can use templates and/or tag dispatch in a more complete
91                              // future implementation.
92                              const auto &paramType = gmxapicompat::mdParamToType(key);
93                              if (paramType == GmxapiType::FLOAT64)
94                              {
95                                  dictionary[key.c_str()] = extractParam(self, key, double());
96                              }
97                              else if (paramType == GmxapiType::INT64)
98                              {
99                                  dictionary[key.c_str()] = extractParam(self, key, int64_t());
100                              }
101                          }
102                          catch (const gmxapicompat::ValueError &e)
103                          {
104                              throw gmxapi::ProtocolError(std::string("Unknown parameter: ") + key);
105                          }
106                      }
107                      return dictionary;
108                  },
109                  "Get a dictionary of the parameters.");
110
111     // Overload a setter for each known type and None
112     mdparams.def("set",
113                  [](GmxMdParams* self, const std::string &key, py::int_ value)
114                  {
115                      gmxapicompat::setParam(self, key, py::cast<int64_t>(value));
116                  },
117                  py::arg("key").none(false),
118                  py::arg("value").none(false),
119                  "Use a dictionary to update simulation parameters.");
120     mdparams.def("set",
121                  [](GmxMdParams* self, const std::string &key, py::float_ value)
122                  {
123                      gmxapicompat::setParam(self, key, py::cast<double>(value));
124                  },
125                  py::arg("key").none(false),
126                  py::arg("value").none(false),
127                  "Use a dictionary to update simulation parameters.");
128     mdparams.def("set",
129                  [](GmxMdParams* self, const std::string &key, py::none)
130                  {
131                      // unsetParam(self, key);
132                  },
133                  py::arg("key").none(false),
134                  py::arg("value"),
135                  "Use a dictionary to update simulation parameters.");
136
137
138     py::class_<TprReadHandle> tprfile(module, "TprFile");
139     tprfile.def("params",
140                 [](const TprReadHandle &self)
141                 {
142                     auto params = gmxapicompat::getMdParams(self);
143                     return params;
144                 });
145
146     module.def("read_tprfile",
147                &readTprFile,
148                py::arg("filename"),
149                "Get a handle to a TPR file resource for a given file name.");
150
151     module.def("write_tprfile",
152                [](std::string filename, const GmxMdParams &parameterObject)
153                {
154                    auto tprReadHandle = gmxapicompat::getSourceFileHandle(parameterObject);
155                    auto params = gmxapicompat::getMdParams(tprReadHandle);
156                    auto structure = gmxapicompat::getStructureSource(tprReadHandle);
157                    auto state = gmxapicompat::getSimulationState(tprReadHandle);
158                    auto topology = gmxapicompat::getTopologySource(tprReadHandle);
159                    gmxapicompat::writeTprFile(filename, *params, *structure, *state, *topology);
160                },
161                py::arg("filename").none(false),
162                py::arg("parameters"),
163                "Write a new TPR file with the provided data.");
164
165     module.def("copy_tprfile",
166                [](const gmxapicompat::TprReadHandle &input, std::string outFile)
167                {
168                    return gmxapicompat::copy_tprfile(input, outFile);
169                },
170                py::arg("source"),
171                py::arg("destination"),
172                "Copy a TPR file from `source` to `destination`."
173                );
174
175     module.def("rewrite_tprfile",
176                [](std::string input, std::string output, double end_time)
177                {
178                    return gmxapicompat::rewrite_tprfile(input, output, end_time);
179                },
180                py::arg("source"),
181                py::arg("destination"),
182                py::arg("end_time"),
183                "Copy a TPR file from `source` to `destination`, replacing `nsteps` with `end_time`.");
184 }
185
186 } // end namespace gmxpy