Rearrange exposure of gmxapicompat simulation input tools.
authorM. Eric Irrgang <ericirrgang@gmail.com>
Fri, 23 Aug 2019 18:22:47 +0000 (21:22 +0300)
committerEric Irrgang <ericirrgang@gmail.com>
Mon, 26 Aug 2019 12:03:42 +0000 (14:03 +0200)
* Add gmxapi/gmxapicompat.h installed header. Move the most benign API
  aspects to this header, including exceptions. Retain some extra
  details in gmxapi/compat/tpr.h and gmxapi/compat/mdparams.h.
* Move the provisional typing enum to gmxapi namespace.
* Change usages of several types to unique_ptr handles to allow opaque
  pointers in gmxapi/gmxapicompat.h

Change-Id: I9c0ae50250b6bc1a1ad3e1885356a5ae718c610e

python_packaging/src/gmxapi/export_tprfile.cpp
src/api/CMakeLists.txt
src/api/cpp/include/gmxapi/compat/exceptions.h [deleted file]
src/api/cpp/include/gmxapi/compat/mdparams.h
src/api/cpp/include/gmxapi/compat/tpr.h
src/api/cpp/include/gmxapi/gmxapi.h
src/api/cpp/include/gmxapi/gmxapicompat.h [new file with mode: 0644]
src/api/cpp/include/gmxapi/gromacsfwd.h
src/api/cpp/tpr.cpp

index 802c43c68ea2647451bb473f31fa1f87c27f263b..8134ba9f68c92df2799fc7ce407e9c0c1d6c7e3a 100644 (file)
 
 #include "gmxapi/exceptions.h"
 
-#include "gmxapi/compat/exceptions.h"
+#include "gmxapi/gmxapicompat.h"
 #include "gmxapi/compat/mdparams.h"
 #include "gmxapi/compat/tpr.h"
 
 #include "module.h"
 
+using gmxapi::GmxapiType;
+
 namespace gmxpy
 {
 
@@ -88,11 +90,11 @@ void detail::export_tprfile(pybind11::module &module)
                              // We can use templates and/or tag dispatch in a more complete
                              // future implementation.
                              const auto &paramType = gmxapicompat::mdParamToType(key);
-                             if (paramType == gmxapicompat::GmxapiType::FLOAT64)
+                             if (paramType == GmxapiType::FLOAT64)
                              {
                                  dictionary[key.c_str()] = extractParam(self, key, double());
                              }
-                             else if (paramType == gmxapicompat::GmxapiType::INT64)
+                             else if (paramType == GmxapiType::INT64)
                              {
                                  dictionary[key.c_str()] = extractParam(self, key, int64_t());
                              }
@@ -154,7 +156,7 @@ void detail::export_tprfile(pybind11::module &module)
                    auto structure = gmxapicompat::getStructureSource(tprReadHandle);
                    auto state = gmxapicompat::getSimulationState(tprReadHandle);
                    auto topology = gmxapicompat::getTopologySource(tprReadHandle);
-                   gmxapicompat::writeTprFile(filename, params, structure, state, topology);
+                   gmxapicompat::writeTprFile(filename, *params, *structure, *state, *topology);
                },
                py::arg("filename").none(false),
                py::arg("parameters"),
index 0d8aea499c1a66acde831e0bb272cbde60e92b84..c3b3f3ba51922125fa1cc16a2344505d80409baf 100644 (file)
@@ -57,18 +57,20 @@ set(GMXAPI_RELEASE ${GMXAPI_MAJOR}.${GMXAPI_MINOR}.${GMXAPI_PATCH})
 # absolute paths in the build tree. Path component before `gmxapi` can be
 # stripped as appropriate by consumers of this list.
 set(GMXAPI_PUBLIC_HEADERS
-    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/compat/exceptions.h
-    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/compat/mdparams.h
-    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/compat/tpr.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/context.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/exceptions.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/gmxapi.h
+    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/gmxapicompat.h
+    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/gromacsfwd.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/md.h
-    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/md/mdmodule.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/session.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/status.h
     ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/system.h
-    ${CMAKE_CURRENT_BINARY_DIR}/cpp/include/gmxapi/version.h)
+    ${CMAKE_CURRENT_BINARY_DIR}/cpp/include/gmxapi/version.h
+    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/compat/mdparams.h
+    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/compat/tpr.h
+    ${CMAKE_CURRENT_SOURCE_DIR}/cpp/include/gmxapi/md/mdmodule.h
+    )
 
 add_subdirectory(cpp)
 
diff --git a/src/api/cpp/include/gmxapi/compat/exceptions.h b/src/api/cpp/include/gmxapi/compat/exceptions.h
deleted file mode 100644 (file)
index 8b74899..0000000
+++ /dev/null
@@ -1,108 +0,0 @@
-/*
- * This file is part of the GROMACS molecular simulation package.
- *
- * Copyright (c) 2019, by the GROMACS development team, led by
- * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
- * and including many others, as listed in the AUTHORS file in the
- * top-level source directory and at http://www.gromacs.org.
- *
- * GROMACS is free software; you can redistribute it and/or
- * modify it under the terms of the GNU Lesser General Public License
- * as published by the Free Software Foundation; either version 2.1
- * of the License, or (at your option) any later version.
- *
- * GROMACS is distributed in the hope that it will be useful,
- * but WITHOUT ANY WARRANTY; without even the implied warranty of
- * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
- * Lesser General Public License for more details.
- *
- * You should have received a copy of the GNU Lesser General Public
- * License along with GROMACS; if not, see
- * http://www.gnu.org/licenses, or write to the Free Software Foundation,
- * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
- *
- * If you want to redistribute modifications to GROMACS, please
- * consider that scientific software is very special. Version
- * control is crucial - bugs must be traceable. We will be happy to
- * consider code for inclusion in the official distribution, but
- * derived work must not be called official GROMACS. Details are found
- * in the README & COPYING files - if they are missing, get the
- * official version at http://www.gromacs.org.
- *
- * To help us fund GROMACS development, we humbly ask that you cite
- * the research papers on the package. Check out http://www.gromacs.org.
- */
-/*! \file
- * \brief C++ exceptions for the gmxapi compatibility tools.
- *
- * Used internally for the gmxapi compatibility helpers that manage type
- * mappings of older GROMACS structures. The long-term disposition of this
- * code is uncertain, but the headers are not likely to be public. If they
- * do persist in some form, we can integrate the exception hierarchy into
- * whatever module takes ownership of this code.
- *
- * Exceptions defined here should only be caught by code that understands the
- * implementation details of these compatibility tools. Exposure of these
- * exceptions outside of the installed object files should be treated as a bug.
- *
- * \author M. Eric Irrgang <ericirrgang@gmail.com>
- * \ingroup gmxapi_compat
- */
-
-#ifndef GMXAPICOMPAT_EXCEPTIONS_H
-#define GMXAPICOMPAT_EXCEPTIONS_H
-
-#include <exception>
-#include <string>
-
-namespace gmxapicompat
-{
-
-/*!
- * \brief Generic exception class for gmxapicompat.
- */
-class Exception : public std::exception
-{
-    public:
-        using std::exception::exception;
-
-        explicit Exception(const std::string &message) :
-            message_ {message}
-        {}
-        explicit Exception(const char* message) : Exception(std::string(message)) {}
-
-        const char *what() const noexcept override
-        {
-            return message_.c_str();
-        }
-
-    private:
-        std::string message_;
-};
-
-/*!
- * \brief The key name provided for a key-value mapping is invalid.
- */
-class KeyError : public Exception
-{
-    using Exception::Exception;
-};
-
-/*!
- * \brief The value provided for a key-value mapping is invalid.
- */
-class ValueError : public Exception
-{
-    using Exception::Exception;
-};
-
-/*!
- * \brief Value provided for a key-value mapping is of an incompatible type.
- */
-class TypeError : public Exception
-{
-    using Exception::Exception;
-};
-
-}      // end namespace gmxapicompat
-#endif //GMXAPICOMPAT_EXCEPTIONS_H
index 4b5f532d57ce3b1dba0c668eec0e811a3fe4850a..cbf07f2a6f699460bacad3837bc3bd3b4ca098ec 100644 (file)
  */
 
 /*! \file
- * \brief Compatibility header for functionality differences in gmxapi releases.
+ * \brief Compatibility header for simulation parameters.
  *
- * Also handle the transitioning installed headers from GROMACS 2019 moving forward.
- *
- * \todo Configure for gmxapi 0.0.7, 0.0.8, GROMACS 2019, GROMACS master...
- *
- * \defgroup gmxapi_compat
  * \author M. Eric Irrgang <ericirrgang@gmail.com>
  * \ingroup gmxapi_compat
  */
 #include <map>
 #include <memory>
 #include <string>
+#include <vector>
 
-struct t_inputrec;
+#include "gmxapi/gmxapicompat.h"
+#include "gmxapi/gmxapi.h"
 
-/*!
- * \brief Compatibility code for features that may not be in gmxapi yet.
- */
 namespace gmxapicompat
 {
 
-
-/*!
- * \brief Label the types recognized by gmxapi.
- *
- * Provide an enumeration to aid in translating data between languages, APIs,
- * and storage formats.
- *
- * \todo The spec should explicitly map these to types in APIs already used.
- * e.g. MPI, Python, numpy, GROMACS, JSON, etc.
- * \todo Actually check the size of the types.
- *
- * \see https://redmine.gromacs.org/issues/2993 for discussion.
- */
-enum class GmxapiType
-{
-    NULLTYPE,     //! Reserved
-    MAP,          //! Mapping of key name (string) to a value of some MdParamType
-    BOOL,         //! Boolean logical type
-    INT64,        //! 64-bit integer type
-    FLOAT64,      //! 64-bit float type
-    STRING,       //! string with metadata
-    NDARRAY,      //! multi-dimensional array with metadata
-};
-
-
-/*!
- * \brief Static map of GROMACS MDP user input to normalized "type".
- *
- * Note that only fields present in the TPR file are named. Additional names
- * may be accepted as mdp file entries, but we cannot discern which parameter
- * name was used from inspection of the TPR file and this is an interim solution
- * that does not need to support a complete MDP file converter.
- */
-std::map<std::string, GmxapiType> simulationParameterTypeMap();
-
-std::map<std::string, bool t_inputrec::*> boolParams();
-std::map<std::string, int t_inputrec::*> int32Params();
-std::map<std::string, float t_inputrec::*> float32Params();
-std::map<std::string, double t_inputrec::*> float64Params();
-std::map<std::string, int64_t t_inputrec::*> int64Params();
-
-/*!
- * \brief Static mapping of parameter names to gmxapi types for GROMACS.
- *
- * \param name MDP entry name.
- * \return enumeration value for known parameters.
- *
- * \throws gmxapi_compat::ValueError for parameters with no mapping.
- */
-GmxapiType mdParamToType(const std::string &name);
-
 // Forward declaration for private implementation class for GmxMdParams
 class GmxMdParamsImpl;
 
-/*!
- * \brief Handle / manager for GROMACS molecular computation input parameters.
- *
- * Interface should be consistent with MDP file entries, but data maps to TPR
- * file interface. For type safety and simplicity, we don't have generic operator
- * accessors. Instead, we have templated accessors that throw exceptions when
- * there is trouble.
- *
- * When MDP input is entirely stored in a key-value tree, this class can be a
- * simple adapter or wrapper. Until then, we need a manually maintained mapping
- * of MDP entries to TPR data.
- *
- * Alternatively, we could update the infrastructure used by list_tpx to provide
- * more generic output, but our efforts may be better spent in updating the
- * infrastructure for the key-value tree input system.
- */
 class GmxMdParams
 {
     public:
@@ -140,27 +67,21 @@ class GmxMdParams
         GmxMdParams(GmxMdParams &&) noexcept;
         GmxMdParams &operator=(GmxMdParams &&) noexcept;
 
+        explicit GmxMdParams(std::unique_ptr<GmxMdParamsImpl> &&impl);
+
         std::unique_ptr<GmxMdParamsImpl> params_;
 };
 
 /*!
- * \brief A set of overloaded functions to fetch parameters of the indicated type, if possible.
+ * \brief Get the list of parameter key words.
  *
- * \param params Handle to a parameters structure from which to extract.
- * \param name Parameter name
- * \param (tag) type for dispatch
+ * \param params molecular simulation parameters object reference.
+ * \return A new vector of parameter names.
  *
- * Could be used for dispatch and/or some sort of templating in the future, but
- * invoked directly for now.
+ * \note The returned data is a copy. Modifying the return value has no affect on
+ * the original object inspected.
  */
-int extractParam(const gmxapicompat::GmxMdParams &params, const std::string &name, int);
-int64_t extractParam(const gmxapicompat::GmxMdParams& params, const std::string& name, int64_t);
-float extractParam(const gmxapicompat::GmxMdParams &params, const std::string &name, float);
-double extractParam(const gmxapicompat::GmxMdParams &params, const std::string &name, double);
-
-void setParam(gmxapicompat::GmxMdParams* params, const std::string &name, double value);
-void setParam(gmxapicompat::GmxMdParams* params, const std::string &name, int64_t value);
-// TODO: unsetParam
+std::vector<std::string> keys(const GmxMdParams &params);
 
 }      // end namespace gmxapicompat
 
index 7a576a7fadd03451810ce12a4169bf5df5bd82c6..0479307ab876b961b2f03356b22fc5aea18d335f 100644 (file)
 #include <memory>
 #include <vector>
 
+#include "gmxapi/gmxapicompat.h"
 #include "gmxapi/compat/mdparams.h"
 
 namespace gmxapicompat
 {
 
-/*!
- * \brief Facade for objects that can provide atomic data for a configuration.
- */
-class StructureSource;
-
-/*!
- * \brief Facade for objects that can provide molecular topology information for a structure.
- */
-class TopologySource;
-
-/*!
- * \brief Proxy to simulation state data.
- */
-class SimulationState;
-
 /*!
  * \brief Manager for TPR file resources.
  *
@@ -86,18 +72,6 @@ class SimulationState;
  */
 class TprContents;
 
-/*!
- * \brief Handle for a TPR data resource.
- *
- * Can provide StructureSource, TopologySource, GmxMdParams, and SimulationState.
- *
- * This is the type of object we allow Python clients to hold references to, though
- * we don't expose any methods to Python. Python clients should acquire access
- * to TPR file contents with read_tpr().
- *
- * \todo gmxapi C++ API should provide mechanisms for subscribing to simulation
- *       input data from various sources.
- */
 class TprReadHandle
 {
     public:
@@ -123,29 +97,6 @@ class TprReadHandle
         std::shared_ptr<TprContents> tprContents_;
 };
 
-/*!
- * \brief Open a TPR file and retrieve a handle.
- *
- * \param filename Path of file to read.
- * \return handle that may share ownership of TPR file resource.
- */
-TprReadHandle readTprFile(const std::string &filename);
-
-/*!
- * \brief Write a new TPR file to the filesystem with the provided contents.
- *
- * \param filename output file path
- * \param params simulation parameters
- * \param structure system structure (atomic configuration)
- * \param state simulation state
- * \param topology molecular topology
- */
-void writeTprFile(const std::string     &filename,
-                  const GmxMdParams     &params,
-                  const StructureSource &structure,
-                  const SimulationState &state,
-                  const TopologySource  &topology);
-
 /*!
  * \brief Helper function for early implementation.
  *
@@ -155,41 +106,6 @@ void writeTprFile(const std::string     &filename,
  */
 TprReadHandle getSourceFileHandle(const GmxMdParams &params);
 
-/*!
- * \brief Get a topology source from the TPR contents collection.
- * \param handle
- * \return
- *
- * \todo replace with a helper template on T::topologySource() member function existence.
- */
-
-TopologySource getTopologySource(const TprReadHandle &handle);
-
-/*!
- * \brief Get a source of simulation state from the TPR contents collection.
- * \param handle
- * \return
- *
- * \todo template on T::simulationState() member function existence.
- */
-SimulationState getSimulationState(const TprReadHandle &handle);
-
-/*!
- * \brief Get a source of atomic structure from the TPR contents collection.
- * \param handle
- * \return
- */
-StructureSource getStructureSource(const TprReadHandle &handle);
-
-/*!
- * \brief Get an initialized parameters structure.
- * \param handle
- * \return
- */
-GmxMdParams getMdParams(const TprReadHandle &handle);
-
-std::vector<std::string> keys(const GmxMdParams &params);
-
 class StructureSource
 {
     public:
@@ -215,7 +131,7 @@ class SimulationState
  * \param outFile output TPR file name
  * \return true if successful. else false.
  */
-bool copy_tprfile(const gmxapicompat::TprReadHandle &input, std::string outFile);
+bool copy_tprfile(const gmxapicompat::TprReadHandle &input, const std::string &outFile);
 
 /*!
  * \brief Copy and possibly update TPR file by name.
@@ -225,7 +141,7 @@ bool copy_tprfile(const gmxapicompat::TprReadHandle &input, std::string outFile)
  * \param endTime Replace `nsteps` in infile with `endTime/dt`
  * \return true if successful, else false
  */
-bool rewrite_tprfile(std::string inFile, std::string outFile, double endTime);
+bool rewrite_tprfile(const std::string &inFile, const std::string &outFile, double endTime);
 
 }      // end namespace gmxapicompat
 
index bc21459b1b1021f1ef3f78035127a449f4df4b06..e38ecd73b1272b38a16f712f11e80d0c36f29562 100644 (file)
@@ -499,6 +499,28 @@ class MDHolder
         /*! \endcond */
 };
 
-}      // end namespace gmxapi
+/*!
+ * \brief Label the types recognized by gmxapi.
+ *
+ * Provide an enumeration to aid in translating data between languages, APIs,
+ * and storage formats.
+ *
+ * \todo The spec should explicitly map these to types in APIs already used.
+ * e.g. MPI, Python, numpy, GROMACS, JSON, etc.
+ * \todo Actually check the size of the types.
+ *
+ * \see https://redmine.gromacs.org/issues/2993 for discussion.
+ */
+enum class GmxapiType
+{
+    NULLTYPE, //! Reserved
+    MAP,      //! Mapping of key name (string) to a value of some MdParamType
+    BOOL,     //! Boolean logical type
+    INT64,    //! 64-bit integer type
+    FLOAT64,  //! 64-bit float type
+    STRING,   //! string with metadata
+    NDARRAY,  //! multi-dimensional array with metadata
+};
+}             // end namespace gmxapi
 
-#endif // header guard
+#endif        // header guard
diff --git a/src/api/cpp/include/gmxapi/gmxapicompat.h b/src/api/cpp/include/gmxapi/gmxapicompat.h
new file mode 100644 (file)
index 0000000..f9dad25
--- /dev/null
@@ -0,0 +1,242 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2019, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+
+/*! \file
+ * \brief Compatibility header for functionality differences in gmxapi releases.
+ *
+ * Also handle the transitioning installed headers from GROMACS 2019 moving forward.
+ *
+ * \todo Configure for gmxapi 0.0.7, 0.0.8, GROMACS 2019, GROMACS master...
+ *
+ * \defgroup gmxapi_compat
+ * \author M. Eric Irrgang <ericirrgang@gmail.com>
+ * \ingroup gmxapi_compat
+ */
+
+#ifndef GMXAPICOMPAT_H
+#define GMXAPICOMPAT_H
+
+#include <map>
+#include <string>
+
+#include "gmxapi/exceptions.h"
+#include "gmxapi/gmxapi.h"
+#include "gmxapi/gromacsfwd.h"
+
+/*!
+ * \brief Compatibility code for features that may not be in the gmxapi specification yet.
+ *
+ * \ingroup gmxapi_compat
+ */
+namespace gmxapicompat
+{
+
+/*!
+ * \brief The key name provided for a key-value mapping is invalid.
+ */
+class KeyError : public gmxapi::BasicException<KeyError>
+{
+    using gmxapi::BasicException<KeyError>::BasicException;
+};
+
+/*!
+ * \brief The value provided for a key-value mapping is invalid.
+ */
+class ValueError : public gmxapi::BasicException<ValueError>
+{
+    using gmxapi::BasicException<ValueError>::BasicException;
+};
+
+/*!
+ * \brief Value provided for a key-value mapping is of an incompatible type.
+ */
+class TypeError : public gmxapi::BasicException<TypeError>
+{
+    using gmxapi::BasicException<TypeError>::BasicException;
+};
+
+/*!
+ * \brief Static map of GROMACS MDP user input to normalized "type".
+ *
+ * Note that only fields present in the TPR file are named. Additional names
+ * may be accepted as mdp file entries, but we cannot discern which parameter
+ * name was used from inspection of the TPR file and this is an interim solution
+ * that does not need to support a complete MDP file converter.
+ */
+std::map<std::string, gmxapi::GmxapiType> simulationParameterTypeMap();
+
+std::map<std::string, bool t_inputrec::*> boolParams();
+std::map<std::string, int t_inputrec::*> int32Params();
+std::map<std::string, float t_inputrec::*> float32Params();
+std::map<std::string, double t_inputrec::*> float64Params();
+std::map<std::string, int64_t t_inputrec::*> int64Params();
+
+/*!
+ * \brief Static mapping of parameter names to gmxapi types for GROMACS.
+ *
+ * \param name MDP entry name.
+ * \return enumeration value for known parameters.
+ *
+ * \throws gmxapi_compat::ValueError for parameters with no mapping.
+ */
+gmxapi::GmxapiType mdParamToType(const std::string &name);
+
+/*!
+ * \brief Facade for objects that can provide atomic data for a configuration.
+ */
+class StructureSource;
+
+/*!
+ * \brief Facade for objects that can provide molecular topology information for a structure.
+ */
+class TopologySource;
+
+/*!
+ * \brief Proxy to simulation state data.
+ */
+class SimulationState;
+
+/*!
+ * \brief Handle / manager for GROMACS molecular computation input parameters.
+ *
+ * Interface should be consistent with MDP file entries, but data maps to TPR
+ * file interface. For type safety and simplicity, we don't have generic operator
+ * accessors. Instead, we have templated accessors that throw exceptions when
+ * there is trouble.
+ *
+ * When MDP input is entirely stored in a key-value tree, this class can be a
+ * simple adapter or wrapper. Until then, we need a manually maintained mapping
+ * of MDP entries to TPR data.
+ *
+ * Alternatively, we could update the infrastructure used by list_tpx to provide
+ * more generic output, but our efforts may be better spent in updating the
+ * infrastructure for the key-value tree input system.
+ */
+class GmxMdParams;
+
+/*!
+ * \brief Handle for a TPR data resource.
+ *
+ * Can provide StructureSource, TopologySource, GmxMdParams, and SimulationState.
+ *
+ * This is the type of object we allow Python clients to hold references to, though
+ * we don't expose any methods to Python. Python clients should acquire access
+ * to TPR file contents with read_tpr().
+ *
+ * \todo gmxapi C++ API should provide mechanisms for subscribing to simulation
+ *       input data from various sources.
+ */
+class TprReadHandle;
+
+/*!
+ * \brief Open a TPR file and retrieve a handle.
+ *
+ * \param filename Path of file to read.
+ * \return handle that may share ownership of TPR file resource.
+ */
+std::unique_ptr<TprReadHandle> readTprFile(const std::string &filename);
+
+/*!
+ * \brief Write a new TPR file to the filesystem with the provided contents.
+ *
+ * \param filename output file path
+ * \param params simulation parameters
+ * \param structure system structure (atomic configuration)
+ * \param state simulation state
+ * \param topology molecular topology
+ *
+ * \throws ValueError for invalid or irreconcilable input.
+ */
+void writeTprFile(const std::string     &filename,
+                  const GmxMdParams     &params,
+                  const StructureSource &structure,
+                  const SimulationState &state,
+                  const TopologySource  &topology);
+
+/*!
+ * \brief Get a topology source from the TPR contents collection.
+ * \param handle
+ * \return
+ *
+ * \todo replace with a helper template on T::topologySource() member function existence.
+ */
+
+std::unique_ptr<TopologySource> getTopologySource(const TprReadHandle &handle);
+
+/*!
+ * \brief Get a source of simulation state from the TPR contents collection.
+ * \param handle
+ * \return
+ *
+ * \todo template on T::simulationState() member function existence.
+ */
+std::unique_ptr<SimulationState> getSimulationState(const TprReadHandle &handle);
+
+/*!
+ * \brief Get a source of atomic structure from the TPR contents collection.
+ * \param handle
+ * \return
+ */
+std::unique_ptr<StructureSource> getStructureSource(const TprReadHandle &handle);
+
+/*!
+ * \brief Get an initialized parameters structure.
+ * \param handle
+ * \return
+ */
+std::unique_ptr<GmxMdParams> getMdParams(const TprReadHandle &handle);
+
+/*!
+ * \brief A set of overloaded functions to fetch parameters of the indicated type, if possible.
+ *
+ * \param params Handle to a parameters structure from which to extract.
+ * \param name Parameter name
+ * \param (tag) type for dispatch
+ *
+ * Could be used for dispatch and/or some sort of templating in the future, but
+ * invoked directly for now.
+ */
+int extractParam(const GmxMdParams &params, const std::string &name, int);
+int64_t extractParam(const GmxMdParams&params, const std::string& name, int64_t);
+float extractParam(const GmxMdParams &params, const std::string &name, float);
+double extractParam(const GmxMdParams &params, const std::string &name, double);
+
+void setParam(GmxMdParams* params, const std::string &name, double value);
+void setParam(GmxMdParams* params, const std::string &name, int64_t value);
+// TODO: unsetParam
+
+}      // end namespace gmxapicompat
+
+#endif //GMXAPICOMPAT_H
index 5a8f66d19857ead30455e3df4bb55e733658ec56..e46f183f914d5145e0182d6da03d6fecf43d1bcd 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
  * Basic API clients only need to compile
  * and link against the gmxapi target, but some gmxapi classes use opaque pointers to
  * library classes that are forward-declared here.
- *
- * We don't want to include ::gmx headers if we don't have to, but we need to declare
- * some things in the ::gmx namespace somewhere. These are forward declarations for
- * opaque pointers in libgromacs for client code building against libgmxapi.
- * Client code that is
- * more entwined with libgromacs can include headers from there.
+ * Client code should not need to include this header directly.
  *
  * For maximal compatibility with other libgmxapi clients (such as third-party
  * Python modules), client code should use the wrappers and protocols in the
- * gmxapi.h header. Note that there is a separate CMake target to build the full
- * developer documentation for gmxapi.
+ * gmxapi.h header.
  *
+ * Note that there is a separate CMake target to build the full
+ * developer documentation for gmxapi.
  * Refer to GMXAPI developer docs for the protocols that map gmxapi interfaces to
  * GROMACS library interfaces.
  * Refer to the GROMACS developer
  * documentation for details on library interfaces forward-declared in this header.
  *
- * \todo It would be nice to include links to the documentation for these classes, too.
+ * \todo Improve documentation cross-linking.
  */
 
+// Forward declaration for src/gromacs/mdtypes/inputrec.h
+struct t_inputrec;
+
 namespace gmx
 {
 
+// Forward declaration for libgromacs header gromacs/restraint/restraintpotential.h
 class IRestraintPotential;
 
 }      // end namespace gmx
index 5caf67d97e254ad0d86d13aadfe32d22b3024442..965544804dc8ab727e088ac7c86edba07633045c 100644 (file)
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/programcontext.h"
 
-#include "gmxapi/compat/exceptions.h"
+#include "gmxapi/gmxapi.h"
+#include "gmxapi/gmxapicompat.h"
 #include "gmxapi/compat/mdparams.h"
 
+using gmxapi::GmxapiType;
+
 namespace gmxapicompat
 {
 
@@ -222,7 +225,7 @@ std::map<std::string, GmxapiType> simulationParameterTypeMap()
 //            TBD
 
     };
-};
+}
 
 /*
  * TODO: Visitor for predetermined known types.
@@ -519,7 +522,7 @@ class GmxMdParamsImpl final
                 if (source_)
                 {
                     auto memberPointer                    = float32Params().at(key);
-                    source_->inputRecord().*memberPointer = value;
+                    source_->inputRecord().*memberPointer = static_cast<float>(value);
                 }
 
             }
@@ -747,41 +750,41 @@ std::vector<std::string> keys(const GmxMdParams &params)
 }
 
 
-TprReadHandle readTprFile(const std::string &filename)
+std::unique_ptr<TprReadHandle> readTprFile(const std::string &filename)
 {
     auto tprfile = gmxapicompat::TprContents(filename);
-    auto handle  = gmxapicompat::TprReadHandle(std::move(tprfile));
+    auto handle  = std::make_unique<gmxapicompat::TprReadHandle>(std::move(tprfile));
     return handle;
 }
 
-GmxMdParams getMdParams(const TprReadHandle &handle)
+std::unique_ptr<GmxMdParams> getMdParams(const TprReadHandle &handle)
 {
     auto tprfile = handle.get();
     // TODO: convert to exception / decide whether null handles are allowed.
     assert(tprfile);
-    GmxMdParams params;
-    params.params_ = std::make_unique<GmxMdParamsImpl>(tprfile);
+    auto params_impl = std::make_unique<GmxMdParamsImpl>(tprfile);
+    auto params      = std::make_unique<GmxMdParams>(std::move(params_impl));
     return params;
 }
 
-TopologySource getTopologySource(const TprReadHandle &handle)
+std::unique_ptr<TopologySource> getTopologySource(const TprReadHandle &handle)
 {
-    TopologySource source;
-    source.tprFile_ = handle.get();
+    auto source = std::make_unique<TopologySource>();
+    source->tprFile_ = handle.get();
     return source;
 }
 
-SimulationState getSimulationState(const TprReadHandle &handle)
+std::unique_ptr<SimulationState> getSimulationState(const TprReadHandle &handle)
 {
-    SimulationState source;
-    source.tprFile_ = handle.get();
+    auto source = std::make_unique<SimulationState>();
+    source->tprFile_ = handle.get();
     return source;
 }
 
-StructureSource getStructureSource(const TprReadHandle &handle)
+std::unique_ptr<StructureSource> getStructureSource(const TprReadHandle &handle)
 {
-    StructureSource source;
-    source.tprFile_ = handle.get();
+    auto source = std::make_unique<StructureSource>();
+    source->tprFile_ = handle.get();
     return source;
 }
 
@@ -846,24 +849,32 @@ GmxMdParams::GmxMdParams(GmxMdParams &&) noexcept = default;
 
 GmxMdParams &GmxMdParams::operator=(GmxMdParams &&) noexcept = default;
 
+GmxMdParams::GmxMdParams(std::unique_ptr<GmxMdParamsImpl> &&impl)
+{
+    // We use swap instead of move construction so that we don't have
+    // to worry about the restrictions on Deleters.
+    // Ref: https://en.cppreference.com/w/cpp/memory/unique_ptr/unique_ptr
+    params_.swap(impl);
+};
+
 // maybe this should return a handle to the new file?
-bool copy_tprfile(const gmxapicompat::TprReadHandle &input, std::string outFile)
+bool copy_tprfile(const gmxapicompat::TprReadHandle &input, const std::string &outFile)
 {
     if (!input.get())
     {
         return false;
     }
     gmxapicompat::writeTprFile(outFile,
-                               gmxapicompat::getMdParams(input),
-                               gmxapicompat::getStructureSource(input),
-                               gmxapicompat::getSimulationState(input),
-                               gmxapicompat::getTopologySource(input));
+                               *gmxapicompat::getMdParams(input),
+                               *gmxapicompat::getStructureSource(input),
+                               *gmxapicompat::getSimulationState(input),
+                               *gmxapicompat::getTopologySource(input));
     return true;
 }
 
-bool rewrite_tprfile(std::string inFile, std::string outFile, double endTime)
+bool rewrite_tprfile(const std::string &inFile, const std::string &outFile, double endTime)
 {
-    bool              success = false;
+    auto              success = false;
 
     const char      * top_fn = inFile.c_str();
 
@@ -887,7 +898,7 @@ bool rewrite_tprfile(std::string inFile, std::string outFile, double endTime)
 
     double  run_t    = irInstance.init_step*irInstance.delta_t + irInstance.init_t;
 
-    irInstance.nsteps = static_cast<int64_t>((endTime - run_t) / irInstance.delta_t + 0.5);
+    irInstance.nsteps = lround((endTime - run_t) / irInstance.delta_t);
 
     write_tpx_state(outFile.c_str(), &irInstance, &state, &mtop);