Merge branch 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 "module.h"
53 #include "mdparams.h"
54 #include "tprfile.h"
55
56 namespace gmxpy
57 {
58
59
60 void detail::export_tprfile(pybind11::module &module)
61 {
62     namespace py = pybind11;
63     using gmxapicompat::GmxMdParams;
64     using gmxapicompat::TprReadHandle;
65     using gmxapicompat::readTprFile;
66
67     py::class_<GmxMdParams> mdparams(module, "SimulationParameters");
68     // We don't want Python users to create invalid params objects, so don't
69     // export a constructor until we can default initialize a valid one.
70     //    mdparams.def(py::init());
71     mdparams.def("extract",
72                  [](const GmxMdParams &self)
73                  {
74                      py::dict dictionary;
75                      for (const auto &key: gmxapicompat::keys(self))
76                      {
77                          try
78                          {
79                              // TODO: More complete typing and dispatching.
80                              // This only handles the two types described in the initial implementation.
81                              // Less trivial types (strings, maps, arrays) warrant additional
82                              // design discussion before being exposed through an interface
83                              // like this one.
84                              // Also reference https://redmine.gromacs.org/issues/2993
85
86                              // We can use templates and/or tag dispatch in a more complete
87                              // future implementation.
88                              const auto &paramType = gmxapicompat::mdParamToType(key);
89                              if (gmxapicompat::isFloat(paramType))
90                              {
91                                  dictionary[key.c_str()] = extractParam(self, key, double());
92                              }
93                              else if (gmxapicompat::isInt(paramType))
94                              {
95                                  dictionary[key.c_str()] = extractParam(self, key, int64_t());
96                              }
97                          }
98                          catch (const gmxapicompat::ValueError &e)
99                          {
100                              throw gmxapi::ProtocolError(std::string("Unknown parameter: ") + key);
101                          }
102                      }
103                      return dictionary;
104                  },
105                  "Get a dictionary of the parameters.");
106
107     // Overload a setter for each known type and None
108     mdparams.def("set",
109                  [](GmxMdParams* self, const std::string &key, py::int_ value)
110                  {
111                      gmxapicompat::setParam(self, key, py::cast<int64_t>(value));
112                  },
113                  py::arg("key").none(false),
114                  py::arg("value").none(false),
115                  "Use a dictionary to update simulation parameters.");
116     mdparams.def("set",
117                  [](GmxMdParams* self, const std::string &key, py::float_ value)
118                  {
119                      gmxapicompat::setParam(self, key, py::cast<double>(value));
120                  },
121                  py::arg("key").none(false),
122                  py::arg("value").none(false),
123                  "Use a dictionary to update simulation parameters.");
124     mdparams.def("set",
125                  [](GmxMdParams* self, const std::string &key, py::none)
126                  {
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",
136                 [](const TprReadHandle &self)
137                 {
138                     auto params = gmxapicompat::getMdParams(self);
139                     return params;
140                 });
141
142     module.def("read_tprfile",
143                &readTprFile,
144                py::arg("filename"),
145                "Get a handle to a TPR file resource for a given file name.");
146
147     module.def("write_tprfile",
148                [](std::string filename, const GmxMdParams &parameterObject)
149                {
150                    auto tprReadHandle = gmxapicompat::getSourceFileHandle(parameterObject);
151                    auto params = gmxapicompat::getMdParams(tprReadHandle);
152                    auto structure = gmxapicompat::getStructureSource(tprReadHandle);
153                    auto state = gmxapicompat::getSimulationState(tprReadHandle);
154                    auto topology = gmxapicompat::getTopologySource(tprReadHandle);
155                    gmxapicompat::writeTprFile(filename, params, structure, state, topology);
156                },
157                py::arg("filename").none(false),
158                py::arg("parameters"),
159                "Write a new TPR file with the provided data.");
160
161     module.def("copy_tprfile",
162                [](const gmxapicompat::TprReadHandle &input, std::string outFile)
163                {
164                    return gmxpy::copy_tprfile(input, outFile);
165                },
166                py::arg("source"),
167                py::arg("destination"),
168                "Copy a TPR file from `source` to `destination`."
169                );
170
171     module.def("rewrite_tprfile",
172                [](std::string input, std::string output, double end_time)
173                {
174                    return gmxpy::rewrite_tprfile(input, output, end_time);
175                },
176                py::arg("source"),
177                py::arg("destination"),
178                py::arg("end_time"),
179                "Copy a TPR file from `source` to `destination`, replacing `nsteps` with `end_time`.");
180 }
181
182 } // end namespace gmxpy