Split off implementation details of MdModuleNotification
authorChristian Blau <cblau@gwdg.de>
Wed, 4 Mar 2020 11:22:34 +0000 (12:22 +0100)
committerChristian Blau <cblau@gerrit.gromacs.org>
Thu, 5 Mar 2020 09:59:35 +0000 (10:59 +0100)
For better readability and testing, moved the MdModuleNotification and
MdModulesNotifier into separate files, hiding away the implementation
details of the complex general callback-infrastructure and exposing more
clearly the callbacks that MdModules can sign up to.

Moved tests into correct directory.

Moved life-line explaining the notification machinery in MDModules from
MdModulesNotification to MdModuleNotifier, where it is better suited.

Change-Id: I0a73f99e450873903a888bcc04228737d0b814f3

src/gromacs/utility/mdmodulenotification-impl.h [new file with mode: 0644]
src/gromacs/utility/mdmodulenotification.h
src/gromacs/utility/tests/CMakeLists.txt
src/gromacs/utility/tests/mdmodulenotification-impl.cpp [moved from src/programs/mdrun/tests/mdmodulenotification.cpp with 97% similarity]
src/programs/mdrun/tests/CMakeLists.txt

diff --git a/src/gromacs/utility/mdmodulenotification-impl.h b/src/gromacs/utility/mdmodulenotification-impl.h
new file mode 100644 (file)
index 0000000..a402504
--- /dev/null
@@ -0,0 +1,167 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2020, 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.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares gmx::MdModuleNotification.
+ *
+ * \author Christian Blau <blau@kth.se>
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+
+#ifndef GMX_UTILITY_MDMODULENOTIFICATION_IMPL_H
+#define GMX_UTILITY_MDMODULENOTIFICATION_IMPL_H
+
+#include <functional>
+#include <vector>
+
+namespace gmx
+{
+
+/*! \libinternal \brief
+ * Subscribe and trigger notification functions.
+ *
+ * Extends MdModuleNotificationBase with new notification function and routine
+ * to subscribe new listeners.
+ *
+ * To create a class of this type that provides callbacks, e.g., for events
+ * EventA, and EventB use registerMdModuleNotification<EventA, EventB>::type.
+ *
+ * \tparam CallParameter of the function to be notified
+ * \tparam MdModuleNotificationBase class to be extended with a notification
+ *                                  with CallParameter
+ *
+ * \note All added subscribers are required to out-live the MdModuleNotification
+ *
+ */
+template<class CallParameter, class MdModuleNotificationBase>
+class MdModuleNotification : public MdModuleNotificationBase
+{
+public:
+    //! Make base class notification trigger available to this class
+    using MdModuleNotificationBase::notify;
+    //! Make base class subscription available to this class
+    using MdModuleNotificationBase::subscribe;
+
+    /*! \brief Trigger the subscribed notifications.
+     * \param[in] callParameter of the function to be called back
+     */
+    void notify(CallParameter callParameter) const
+    {
+        for (auto& callBack : callBackFunctions_)
+        {
+            callBack(callParameter);
+        }
+    }
+
+    /*! \brief
+     * Add callback function to be called when notification is triggered.
+     *
+     * Notifications are distinguished by their call signature.
+     *
+     * \param[in] callBackFunction to be called from this class
+     */
+    void subscribe(std::function<void(CallParameter)> callBackFunction)
+    {
+        callBackFunctions_.emplace_back(callBackFunction);
+    }
+
+private:
+    std::vector<std::function<void(CallParameter)>> callBackFunctions_;
+};
+
+/*! \internal
+ * \brief Aide to avoid nested MdModuleNotification definition.
+ *
+ * Instead of
+ * MdModuleNotification<CallParameterA, MdModuleNotification<CallParameterB, etc ... >>
+ * this allows to write
+ * registerMdModuleNotification<CallParameterA, CallParameterB, ...>::type
+ *
+ * \tparam CallParameter all the event types to be registered
+ */
+template<class... CallParameter>
+struct registerMdModuleNotification;
+
+/*! \internal \brief Template specialization to end parameter unpacking recursion.
+ */
+template<>
+struct registerMdModuleNotification<>
+{
+    /*! \internal
+     * \brief Do nothing but be base class of MdModuleNotification.
+     *
+     * Required so that using MdModuleNotificationBase::notify and
+     * MdModuleNotificationBase::subscribe are valid in derived class.
+     */
+    class NoCallParameter
+    {
+    public:
+        //! Do nothing but provide MdModuleNotification::notify to derived class
+        void notify() {}
+        //! Do nothing but provide MdModuleNotification::subscribe to derived class
+        void subscribe() {}
+    };
+    /*! \brief Defines a type if no notifications are managed.
+     *
+     * This ensures that code works with MdModuleCallParameterManagement that
+     * does not manage any notifications.
+     */
+    using type = NoCallParameter;
+};
+
+/*! \libinternal
+ * \brief Template specialization to assemble MdModuleNotification.
+ *
+ * Assembly of MdModuleNotification is performed by recursively taking off the
+ * front of the CallParameter parameter pack and constructing the nested type
+ * definition of MdModuleNotification base classes.
+ *
+ * \tparam CurrentCallParameter front of the template parameter pack
+ * \tparam CallParameter rest of the event types
+ */
+template<class CurrentCallParameter, class... CallParameter>
+struct registerMdModuleNotification<CurrentCallParameter, CallParameter...>
+{
+    // private:
+    //! The next type with rest of the arguments with the front parameter removed.
+    using next_type = typename registerMdModuleNotification<CallParameter...>::type;
+    //! The type of the MdModuleNotification
+    using type = MdModuleNotification<CurrentCallParameter, next_type>;
+};
+
+} // namespace gmx
+
+#endif
index 1f4b593651892ee35a478c76eee19b0ece79e5f4..b5df1e334096f7213d7e146870077140e390b936 100644 (file)
  */
 /*! \libinternal \file
  * \brief
- * Declares gmx::MdModuleNotification.
+ * Declares gmx::MdModulesNotifier.
  *
  * \author Christian Blau <blau@kth.se>
  * \inlibraryapi
  * \ingroup module_utility
  */
 
-#ifndef GMX_MDRUNUTILITY_MDMODULENOTIFICATION_H
-#define GMX_MDRUNUTILITY_MDMODULENOTIFICATION_H
+#ifndef GMX_UTILITY_MDMODULENOTIFICATION_H
+#define GMX_UTILITY_MDMODULENOTIFICATION_H
 
-#include <functional>
 #include <string>
 #include <vector>
 
+#include "gromacs/utility/mdmodulenotification-impl.h"
+
 struct t_commrec;
 enum class PbcType : int;
 
 namespace gmx
 {
 
-/*! \libinternal \brief
- * Subscribe and trigger notification functions.
- *
- * Extends MdModuleNotificationBase with new notification function and routine
- * to subscribe new listeners.
- *
- * To create a class of this type that provides callbacks, e.g., for events
- * EventA, and EventB use registerMdModuleNotification<EventA, EventB>::type.
- *
- * \tparam CallParameter of the function to be notified
- * \tparam MdModuleNotificationBase class to be extended with a notification
- *                                  with CallParameter
- *
-   \msc
-   wordwraparcs=true,
-   hscale="2";
-
-   runner [label="runner:\nMdrunner"],
-   CallParameter [label = "eventA:\nCallParameter"],
-   MOD [label = "mdModules_:\nMdModules"],
-   ModuleA [label="moduleA"],
-   ModuleB [label="moduleB"],
-   MdModuleNotification [label="notifier_:\nMdModuleNotification"];
-
-   MOD box MdModuleNotification [label = "mdModules_ owns notifier_ and moduleA/B"];
-   MOD =>> ModuleA [label="instantiates(notifier_)"];
-   ModuleA =>> MdModuleNotification [label="subscribe(otherfunc)"];
-   ModuleA =>> MOD;
-   MOD =>> ModuleB [label="instantiates(notifier_)"];
-   ModuleB =>> MdModuleNotification [label="subscribe(func)"];
-   ModuleB =>> MOD;
-   runner =>> CallParameter [label="instantiate"];
-   CallParameter =>> runner ;
-   runner =>> MOD [label="notify(eventA)"];
-   MOD =>> MdModuleNotification [label="notify(eventA)"];
-   MdModuleNotification =>> ModuleA [label="notify(eventA)"];
-   ModuleA -> ModuleA [label="func(eventA)"];
-   MdModuleNotification =>> ModuleB [label="notify(eventA)"];
-   ModuleB -> ModuleB [label="otherfunc(eventA)"];
-
-   \endmsc
- *
- * \note All added subscribers are required to out-live the MdModuleNotification
- *
- */
-template<class CallParameter, class MdModuleNotificationBase>
-class MdModuleNotification : public MdModuleNotificationBase
-{
-public:
-    //! Make base class notification trigger available to this class
-    using MdModuleNotificationBase::notify;
-    //! Make base class subscription available to this class
-    using MdModuleNotificationBase::subscribe;
-
-    /*! \brief Trigger the subscribed notifications.
-     * \param[in] callParameter of the function to be called back
-     */
-    void notify(CallParameter callParameter) const
-    {
-        for (auto& callBack : callBackFunctions_)
-        {
-            callBack(callParameter);
-        }
-    }
-
-    /*! \brief
-     * Add callback function to be called when notification is triggered.
-     *
-     * Notifications are distinguished by their call signature.
-     *
-     * \param[in] callBackFunction to be called from this class
-     */
-    void subscribe(std::function<void(CallParameter)> callBackFunction)
-    {
-        callBackFunctions_.emplace_back(callBackFunction);
-    }
-
-private:
-    std::vector<std::function<void(CallParameter)>> callBackFunctions_;
-};
-
-/*! \internal
- * \brief Aide to avoid nested MdModuleNotification definition.
- *
- * Instead of
- * MdModuleNotification<CallParameterA, MdModuleNotification<CallParameterB, etc ... >>
- * this allows to write
- * registerMdModuleNotification<CallParameterA, CallParameterB, ...>::type
- *
- * \tparam CallParameter all the event types to be registered
- */
-template<class... CallParameter>
-struct registerMdModuleNotification;
-
-/*! \internal \brief Template specialization to end parameter unpacking recursion.
- */
-template<>
-struct registerMdModuleNotification<>
-{
-    /*! \internal
-     * \brief Do nothing but be base class of MdModuleNotification.
-     *
-     * Required so that using MdModuleNotificationBase::notify and
-     * MdModuleNotificationBase::subscribe are valid in derived class.
-     */
-    class NoCallParameter
-    {
-    public:
-        //! Do nothing but provide MdModuleNotification::notify to derived class
-        void notify() {}
-        //! Do nothing but provide MdModuleNotification::subscribe to derived class
-        void subscribe() {}
-    };
-    /*! \brief Defines a type if no notifications are managed.
-     *
-     * This ensures that code works with MdModuleCallParameterManagement that
-     * does not manage any notifications.
-     */
-    using type = NoCallParameter;
-};
-
-/*! \libinternal
- * \brief Template specialization to assemble MdModuleNotification.
- *
- * Assembly of MdModuleNotification is performed by recursively taking off the
- * front of the CallParameter parameter pack and constructing the nested type
- * definition of MdModuleNotification base classes.
- *
- * \tparam CurrentCallParameter front of the template parameter pack
- * \tparam CallParameter rest of the event types
- */
-template<class CurrentCallParameter, class... CallParameter>
-struct registerMdModuleNotification<CurrentCallParameter, CallParameter...>
-{
-    // private:
-    //! The next type with rest of the arguments with the front parameter removed.
-    using next_type = typename registerMdModuleNotification<CallParameter...>::type;
-    //! The type of the MdModuleNotification
-    using type = MdModuleNotification<CurrentCallParameter, next_type>;
-};
-
 class KeyValueTreeObject;
 class KeyValueTreeObjectBuilder;
 class LocalAtomSetManager;
@@ -203,8 +63,13 @@ struct MdModulesCheckpointReadingDataOnMaster;
 struct MdModulesCheckpointReadingBroadcast;
 struct MdModulesWriteCheckpointData;
 
+/*! \libinternal \brief Check if module outputs energy to a specific field.
+ *
+ * Ensures that energy is output for this module.
+ */
 struct MdModulesEnergyOutputToDensityFittingRequestChecker
 {
+    //! Trigger output to density fitting energy field
     bool energyOutputToDensityFitting_ = false;
 };
 
@@ -243,6 +108,8 @@ private:
     std::vector<std::string> errorMessages_;
 };
 
+/*! \libinternal \brief Provides the simulation time step in ps.
+ */
 struct SimulationTimeStep
 {
     //! Time step (ps)
@@ -259,6 +126,35 @@ struct SimulationTimeStep
  *  When pre-processing the simulation data
  *  When reading and writing check-pointing data
  *  When setting up simulation after reading in the tpr file
+ *
+   \msc
+   wordwraparcs=true,
+   hscale="2";
+
+   runner [label="runner:\nMdrunner"],
+   CallParameter [label = "eventA:\nCallParameter"],
+   MOD [label = "mdModules_:\nMdModules"],
+   ModuleA [label="moduleA"],
+   ModuleB [label="moduleB"],
+   MdModuleNotification [label="notifier_:\nMdModuleNotification"];
+
+   MOD box MdModuleNotification [label = "mdModules_ owns notifier_ and moduleA/B"];
+   MOD =>> ModuleA [label="instantiates(notifier_)"];
+   ModuleA =>> MdModuleNotification [label="subscribe(otherfunc)"];
+   ModuleA =>> MOD;
+   MOD =>> ModuleB [label="instantiates(notifier_)"];
+   ModuleB =>> MdModuleNotification [label="subscribe(func)"];
+   ModuleB =>> MOD;
+   runner =>> CallParameter [label="instantiate"];
+   CallParameter =>> runner ;
+   runner =>> MOD [label="notify(eventA)"];
+   MOD =>> MdModuleNotification [label="notify(eventA)"];
+   MdModuleNotification =>> ModuleA [label="notify(eventA)"];
+   ModuleA -> ModuleA [label="func(eventA)"];
+   MdModuleNotification =>> ModuleB [label="notify(eventA)"];
+   ModuleB -> ModuleB [label="otherfunc(eventA)"];
+
+   \endmsc
  *
  * The template arguments to the members of this struct directly reflect
  * the callback function signature. Arguments passed as pointers are always
index 701b1f6ffbce20722c21c2d8a05e66bb747518c9..5c3fb9161434dadcb0c78e3605f73e1447fe8f50 100644 (file)
@@ -46,6 +46,7 @@ gmx_add_unit_test(UtilityUnitTests utility-test
                   keyvaluetreetransform.cpp
                   listoflists.cpp
                   logger.cpp
+                  mdmodulenotification-impl.cpp
                   mutex.cpp
                   path.cpp
                   physicalnodecommunicator.cpp
similarity index 97%
rename from src/programs/mdrun/tests/mdmodulenotification.cpp
rename to src/gromacs/utility/tests/mdmodulenotification-impl.cpp
index a2a32ecc8363d4e24f126fa25d9299bbf8f6721c..d77b3848e46fd88f5be49458e2a7f121b006f024 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020, 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.
  * Tests MdModuleNotification
  *
  * \author Christian Blau <blau@kth.se>
- * \ingroup module_mdrun_integration_tests
+ * \ingroup module_utility
  */
 #include "gmxpre.h"
 
-#include <gmock/gmock.h>
+#include "gromacs/utility/mdmodulenotification-impl.h"
 
-#include "gromacs/utility/mdmodulenotification.h"
+#include <gmock/gmock.h>
 
 namespace gmx
 {
index 24af506d980463cd454b414bf83cd17ea49fbb41..aa8d74e9514296d07b36dad2aa17c9b47ec155ce 100644 (file)
@@ -39,7 +39,6 @@ gmx_add_unit_test_library(mdrun_test_infrastructure
     energyreader.cpp
     energycomparison.cpp
     moduletest.cpp
-    mdmodulenotification.cpp
     terminationhelper.cpp
     trajectorycomparison.cpp
     trajectoryreader.cpp