Provide callbacks/notifications for MDModules
authorChristian Blau <cblau@gwdg.de>
Tue, 23 Jul 2019 15:25:06 +0000 (17:25 +0200)
committerChristian Blau <cblau@gwdg.de>
Thu, 8 Aug 2019 07:47:18 +0000 (09:47 +0200)
Adds functionality for MdModules to subscribe to be called back during
the simulation. Within the run, subscribed modules are notified of
events, that are distinguished by the function argument by the call back
function.

Implements the callbacks based on storing function pointers, following
the discussion in https://gerrit.gromacs.org/c/gromacs/+/10942

refs #2945

Change-Id: I61589215fa9beb79825f0b5261ed50b4116046ff

docs/doxygen/lib/mdmodules.md
src/gromacs/mdrun/mdmodulenotification.h [new file with mode: 0644]
src/gromacs/mdrun/mdmodules.cpp
src/gromacs/mdrun/mdmodules.h
src/gromacs/mdrun/runner.cpp
src/programs/mdrun/tests/CMakeLists.txt
src/programs/mdrun/tests/mdmodulenotification.cpp [new file with mode: 0644]

index b0fa0b17478fc6908b622bc092ab41f5c8659e5a..3da16c19b3a8163f622fa18c0430e6943baa7e93 100644 (file)
@@ -44,7 +44,7 @@ class, which is expressed with the final keyword in the class
 definition. Generally, modules will implement different flavours of
 functionality, perhaps based on user choices, or available computing
 resources. This should generally be implemented by providing variable
-behaviour for the methods that are called through the above
+behavior for the methods that are called through the above
 interfaces. Either code should branch at run time upon some data
 contained by the module (e.g. read from the mdp options), or that the
 module class should contain a pointer to an internal interface class
@@ -152,3 +152,42 @@ for initMdpTransform() to support specifying version-to-version conversions for
 any changed options, and applying the necessary conversions in sequence.  The
 main challenge is keeping track of the versions to know which conversions to
 apply.
+
+Callbacks to modules during setup and simulation
+------------------------------------------------
+
+During setup and simulation, modules receive required information like topologies
+and local atom sets by subscribing to callback functions.
+
+To include a notification for your module,  
+
+* Add the function signature for the callback function to the `notifier_type`
+  in the MdModules class
+
+  ```C++
+    notifier_type = registerMdModuleNotification<...,
+                      YourCallbackSignature,
+                      ...,
+  ```
+
+  (keep alphabetical order for ease of git merge)
+
+* Hand the notifier_ member of the MdModules Implementation class to your
+  builder createYourModule(&notifier_)
+
+* Add the function you want to subscribe with in the builder,
+  `notifier->subscribe(yourFunction)`
+  
+  * To subscribe class member functions of your module, you can use lambda expressions
+
+  ```C++
+    notifier->subscribe([modulePointer = yourModule.get()]
+      (YourCallbackSignature argument){modulePointer(argument);});
+  ```
+
+* During setup in , e.g., within `Mdrunner` use
+
+  ```C++
+    YourCallbackSignature argument();
+    mdModules_.notifier().notify(argument);
+  ```
diff --git a/src/gromacs/mdrun/mdmodulenotification.h b/src/gromacs/mdrun/mdmodulenotification.h
new file mode 100644 (file)
index 0000000..dcb273a
--- /dev/null
@@ -0,0 +1,199 @@
+/*
+ * 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.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares gmx::MdModuleNotification.
+ *
+ * \author Christian Blau <blau@kth.se>
+ * \inlibraryapi
+ * \ingroup module_mdrun
+ */
+
+#ifndef GMX_MDRUN_MDMODULENOTIFICATION_H
+#define GMX_MDRUN_MDMODULENOTIFICATION_H
+
+#include <functional>
+#include <vector>
+
+#include "gromacs/utility/basedefinitions.h"
+
+struct t_commrec;
+
+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>;
+};
+
+} // namespace gmx
+
+#endif
index ad2efe339afc9511f7cc3bac9e39925b0430c1eb..afc285a141c779106d468269977d5a83da056163 100644 (file)
@@ -102,6 +102,9 @@ class MDModules::Impl : public IMDOutputProvider
          * \todo include field_ in modules_
          */
         std::vector< std::shared_ptr<IMDModule> > modules_;
+
+        //! Manages resources and notifies the MD modules when available
+        MDModules::notifier_type notifier_;
 };
 
 MDModules::MDModules() : impl_(new Impl)
@@ -170,4 +173,9 @@ void MDModules::add(std::shared_ptr<gmx::IMDModule> module)
     impl_->modules_.emplace_back(std::move(module));
 }
 
+const MDModules::notifier_type &MDModules::notifier()
+{
+    return impl_->notifier_;
+}
+
 } // namespace gmx
index c069bbb4fa56cda37a46854a654fd32afd960099..e3c816f8abf9693b95d4b2b2d9039596391d78a3 100644 (file)
@@ -43,6 +43,7 @@
 #ifndef GMX_MDRUN_MDMODULES_H
 #define GMX_MDRUN_MDMODULES_H
 
+#include "gromacs/mdrun/mdmodulenotification.h"
 #include "gromacs/utility/classhelpers.h"
 
 struct ForceProviders;
@@ -59,6 +60,17 @@ class IKeyValueTreeTransformRules;
 class IMDOutputProvider;
 class KeyValueTreeObject;
 class IMDModule;
+class LocalAtomSetManager;
+
+/*! \libinternal \brief
+ * \brief Signals that the communication record is set up and provides this record.
+ */
+struct CommunicationIsSetup
+{
+    //! The communication record that is set up.
+    const t_commrec &communicationRecord_;
+};
+
 
 /*! \libinternal \brief
  * Manages the collection of all modules used for mdrun.
@@ -95,6 +107,12 @@ class MDModules
         MDModules();
         ~MDModules();
 
+        //! Register callback function types for MDModules
+        using notifier_type = registerMdModuleNotification<
+                    CommunicationIsSetup,
+                    LocalAtomSetManager *
+                    >::type;
+
         /*! \brief
          * Initializes a transform from mdp values to sectioned options.
          *
@@ -169,6 +187,10 @@ class MDModules
          */
         void add(std::shared_ptr<gmx::IMDModule> module);
 
+        /*! \brief Return a handle to the callbacks.
+         */
+        const notifier_type &notifier();
+
     private:
         class Impl;
 
index deef8f5b9e72b132f9ebcdc49315cfffb039b177..6e01ac3c470eadfe02d5cbeae0b8c9b6fb3bd4f0 100644 (file)
@@ -98,6 +98,7 @@
 #include "gromacs/mdlib/qmmm.h"
 #include "gromacs/mdlib/sighandler.h"
 #include "gromacs/mdlib/stophandler.h"
+#include "gromacs/mdrun/mdmodulenotification.h"
 #include "gromacs/mdrun/mdmodules.h"
 #include "gromacs/mdrun/simulationcontext.h"
 #include "gromacs/mdrunutility/handlerestart.h"
@@ -986,6 +987,7 @@ int Mdrunner::mdrunner()
                                            &mtop, inputrec,
                                            box, positionsFromStatePointer(globalState.get()),
                                            &atomSets);
+        mdModules_->notifier().notify(&atomSets);
         // Note that local state still does not exist yet.
     }
     else
@@ -1280,6 +1282,7 @@ int Mdrunner::mdrunner()
     t_nrnb nrnb;
     if (thisRankHasDuty(cr, DUTY_PP))
     {
+        mdModules_->notifier().notify(CommunicationIsSetup {*cr});
         /* Initiate forcerecord */
         fr                 = new t_forcerec;
         fr->forceProviders = mdModules_->initForceProviders();
index c7179463e2f342ccfd7720fe0293cdf181dce050..676fe6f25d03490cc03246a001635366f8f5e14e 100644 (file)
@@ -38,6 +38,7 @@ gmx_add_unit_test_library(mdrun_test_infrastructure
     energyreader.cpp
     energycomparison.cpp
     moduletest.cpp
+    mdmodulenotification.cpp
     terminationhelper.cpp
     trajectorycomparison.cpp
     trajectoryreader.cpp
diff --git a/src/programs/mdrun/tests/mdmodulenotification.cpp b/src/programs/mdrun/tests/mdmodulenotification.cpp
new file mode 100644 (file)
index 0000000..da70de7
--- /dev/null
@@ -0,0 +1,186 @@
+/*
+ * 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.
+ */
+/*! \internal \file
+ * \brief
+ * Tests MdModuleNotification
+ *
+ * \author Christian Blau <blau@kth.se>
+ * \ingroup module_mdrun_integration_tests
+ */
+#include "gmxpre.h"
+
+#include <gmock/gmock.h>
+
+#include "gromacs/mdrun/mdmodulenotification.h"
+
+namespace gmx
+{
+
+namespace
+{
+
+struct EventA{};
+struct EventB{};
+
+class EventACallee final
+{
+    public:
+        void callback(EventA /*a*/)
+        {
+            notifiedEventA_ = true;
+        };
+
+        bool notifiedEventA() { return notifiedEventA_; }
+
+    private:
+        bool notifiedEventA_ = false;
+};
+
+class EventBCallee final
+{
+    public:
+        void callback(EventB * /* bPointer */)
+        {
+            notifiedEventB_ = true;
+        };
+
+        bool notifiedEventB() { return notifiedEventB_; }
+
+    private:
+        bool notifiedEventB_ = false;
+};
+
+class EventAandBCallee final
+{
+    public:
+        void notify(EventB * /* bPointer */)
+        {
+            notifiedEventB_ = true;
+        };
+
+        void callback(EventA /* a */)
+        {
+            notifiedEventA_ = true;
+        };
+
+        bool notifiedEventB() { return notifiedEventB_; }
+        bool notifiedEventA() { return notifiedEventA_; }
+
+    private:
+        bool notifiedEventB_ = false;
+        bool notifiedEventA_ = false;
+};
+
+TEST(MDModuleNotificationTest, addConsumer)
+{
+    registerMdModuleNotification<EventA>::type notifications;
+    EventACallee eventACallee;
+
+    EXPECT_FALSE(eventACallee.notifiedEventA());
+
+    notifications.subscribe([&eventACallee](EventA eventA) { eventACallee.callback(eventA); });
+    notifications.notify(EventA {});
+
+    EXPECT_TRUE(eventACallee.notifiedEventA());
+}
+
+TEST(MDModuleNotificationTest, addConsumerWithPointerParameter)
+{
+    registerMdModuleNotification<EventB *>::type notifications;
+    EventBCallee eventBCallee;
+
+    EXPECT_FALSE(eventBCallee.notifiedEventB());
+
+    notifications.subscribe([&eventBCallee](EventB * eventB) { eventBCallee.callback(eventB); });
+    EventB * eventBPointer = nullptr;
+    notifications.notify(eventBPointer);
+
+    EXPECT_TRUE(eventBCallee.notifiedEventB());
+}
+
+TEST(MDModuleNotificationTest, addTwoDifferentConsumers)
+{
+    registerMdModuleNotification<EventA, EventB *>::type notifications;
+    EventBCallee eventBCallee;
+    EventACallee eventACallee;
+
+    EXPECT_FALSE(eventACallee.notifiedEventA());
+    EXPECT_FALSE(eventBCallee.notifiedEventB());
+
+    notifications.subscribe([&eventBCallee](EventB *eventB) { eventBCallee.callback(eventB); });
+    notifications.subscribe([&eventACallee](EventA eventA) { eventACallee.callback(eventA); });
+
+    EventB * eventBPointer = nullptr;
+    notifications.notify(eventBPointer);
+
+    EXPECT_FALSE(eventACallee.notifiedEventA());
+    EXPECT_TRUE(eventBCallee.notifiedEventB());
+
+    notifications.notify(EventA {});
+
+    EXPECT_TRUE(eventACallee.notifiedEventA());
+    EXPECT_TRUE(eventBCallee.notifiedEventB());
+}
+
+TEST(MDModuleNotificationTest, consumerOfTwoResources)
+{
+    registerMdModuleNotification<EventA, EventB *>::type notifications;
+
+    EventAandBCallee     callee;
+
+    EXPECT_FALSE(callee.notifiedEventB());
+    EXPECT_FALSE(callee.notifiedEventA());
+
+    // requires a template parameter here, because call is ambiguous otherwise
+    notifications.subscribe([&callee](EventA msg) { callee.callback(msg); });
+    notifications.subscribe([&callee](EventB *msg) { callee.notify(msg); });
+
+    EventB * eventBp = nullptr;
+
+    notifications.notify(eventBp);
+
+    EXPECT_FALSE(callee.notifiedEventA());
+    EXPECT_TRUE(callee.notifiedEventB());
+
+    notifications.notify(EventA {});
+
+    EXPECT_TRUE(callee.notifiedEventA());
+    EXPECT_TRUE(callee.notifiedEventB());
+}
+
+
+} // namespace
+
+} // namespace gmx