Split MdModuleNotifications into topics
authorChristian Blau <cblau@gwdg.de>
Fri, 17 Jan 2020 13:15:14 +0000 (14:15 +0100)
committerChristian Blau <cblau@gerrit.gromacs.org>
Tue, 18 Feb 2020 10:44:16 +0000 (11:44 +0100)
Use different notifers to distinguish MdModule callbacks to

 - preprocessing
 - simulation setup
 - checkpointing

refs #3076

Change-Id: I08b2ff7e90e42fce44043d7c0f2a9304d404445c

src/gromacs/applied_forces/densityfitting.cpp
src/gromacs/fileio/checkpoint.cpp
src/gromacs/gmxpreprocess/grompp.cpp
src/gromacs/gmxpreprocess/readir.cpp
src/gromacs/mdlib/energyoutput.cpp
src/gromacs/mdrun/runner.cpp
src/gromacs/utility/mdmodulenotification.h

index d5d8ea90afbaabeb1b494866f8af9e21e14607ac..5a134e4327ab0a1d9c929d5511a895cb07dae96f 100644 (file)
@@ -176,12 +176,11 @@ public:
      * allow the MdModuleNotifier to statically distinguish the callback
      * function type from other 'int' function callbacks.
      *
-     * \param[in] pbc MdModuleNotification class that contains a variable
-     *                that enumerates the periodic boundary condition.
+     * \param[in] pbcType enumerates the periodic boundary condition.
      */
-    void setPeriodicBoundaryConditionType(PeriodicBoundaryConditionType pbc)
+    void setPeriodicBoundaryConditionType(const PbcType& pbcType)
     {
-        pbcType_ = std::make_unique<PbcType>(pbc.pbcType);
+        pbcType_ = std::make_unique<PbcType>(pbcType);
     }
 
     //! Get the periodic boundary conditions
@@ -230,8 +229,8 @@ public:
      * \param[in] notifier allows the module to subscribe to notifications from MdModules.
      *
      * The density fitting code subscribes to these notifications:
-     *   - setting atom group indices in the densityFittingOptions_ by
-     *     taking a parmeter const IndexGroupsAndNames &
+     *   - setting atom group indices in the densityFittingOptions_ from an
+     *     index group string by taking a parmeter const IndexGroupsAndNames &
      *   - storing its internal parameters in a tpr file by writing to a
      *     key-value-tree during pre-processing by a function taking a
      *     KeyValueTreeObjectBuilder as parameter
@@ -254,73 +253,73 @@ public:
         // and subscribed, and will be dispatched correctly at run time
         // based on the type of the parameter required by the lambda.
 
-        // Setting atom group indices
+        // Setting the atom group indices from index group string
         const auto setFitGroupIndicesFunction = [this](const IndexGroupsAndNames& indexGroupsAndNames) {
             densityFittingOptions_.setFitGroupIndices(indexGroupsAndNames);
         };
-        notifier->notifier_.subscribe(setFitGroupIndicesFunction);
+        notifier->preProcessingNotifications_.subscribe(setFitGroupIndicesFunction);
 
         // Writing internal parameters during pre-processing
         const auto writeInternalParametersFunction = [this](KeyValueTreeObjectBuilder treeBuilder) {
             densityFittingOptions_.writeInternalParametersToKvt(treeBuilder);
         };
-        notifier->notifier_.subscribe(writeInternalParametersFunction);
+        notifier->preProcessingNotifications_.subscribe(writeInternalParametersFunction);
 
         // Reading internal parameters during simulation setup
         const auto readInternalParametersFunction = [this](const KeyValueTreeObject& tree) {
             densityFittingOptions_.readInternalParametersFromKvt(tree);
         };
-        notifier->notifier_.subscribe(readInternalParametersFunction);
+        notifier->simulationSetupNotifications_.subscribe(readInternalParametersFunction);
 
         // Checking for consistency with all .mdp options
         const auto checkEnergyCaluclationFrequencyFunction =
                 [this](EnergyCalculationFrequencyErrors* energyCalculationFrequencyErrors) {
                     densityFittingOptions_.checkEnergyCaluclationFrequency(energyCalculationFrequencyErrors);
                 };
-        notifier->notifier_.subscribe(checkEnergyCaluclationFrequencyFunction);
+        notifier->preProcessingNotifications_.subscribe(checkEnergyCaluclationFrequencyFunction);
 
         // constructing local atom sets during simulation setup
         const auto setLocalAtomSetFunction = [this](LocalAtomSetManager* localAtomSetManager) {
             this->constructLocalAtomSet(localAtomSetManager);
         };
-        notifier->notifier_.subscribe(setLocalAtomSetFunction);
+        notifier->simulationSetupNotifications_.subscribe(setLocalAtomSetFunction);
 
         // constructing local atom sets during simulation setup
-        const auto setPeriodicBoundaryContionsFunction = [this](PeriodicBoundaryConditionType pbc) {
+        const auto setPeriodicBoundaryContionsFunction = [this](const PbcType& pbc) {
             this->densityFittingSimulationParameters_.setPeriodicBoundaryConditionType(pbc);
         };
-        notifier->notifier_.subscribe(setPeriodicBoundaryContionsFunction);
+        notifier->simulationSetupNotifications_.subscribe(setPeriodicBoundaryContionsFunction);
 
         // setting the simulation time step
         const auto setSimulationTimeStepFunction = [this](const SimulationTimeStep& simulationTimeStep) {
             this->densityFittingSimulationParameters_.setSimulationTimeStep(simulationTimeStep.delta_t);
         };
-        notifier->notifier_.subscribe(setSimulationTimeStepFunction);
+        notifier->simulationSetupNotifications_.subscribe(setSimulationTimeStepFunction);
 
         // adding output to energy file
         const auto requestEnergyOutput =
                 [this](MdModulesEnergyOutputToDensityFittingRequestChecker* energyOutputRequest) {
                     this->setEnergyOutputRequest(energyOutputRequest);
                 };
-        notifier->notifier_.subscribe(requestEnergyOutput);
+        notifier->simulationSetupNotifications_.subscribe(requestEnergyOutput);
 
         // writing checkpoint data
         const auto checkpointDataWriting = [this](MdModulesWriteCheckpointData checkpointData) {
             this->writeCheckpointData(checkpointData);
         };
-        notifier->notifier_.subscribe(checkpointDataWriting);
+        notifier->checkpointingNotifications_.subscribe(checkpointDataWriting);
 
         // reading checkpoint data
         const auto checkpointDataReading = [this](MdModulesCheckpointReadingDataOnMaster checkpointData) {
             this->readCheckpointDataOnMaster(checkpointData);
         };
-        notifier->notifier_.subscribe(checkpointDataReading);
+        notifier->checkpointingNotifications_.subscribe(checkpointDataReading);
 
         // broadcasting checkpoint data
         const auto checkpointDataBroadcast = [this](MdModulesCheckpointReadingBroadcast checkpointData) {
             this->broadcastCheckpointData(checkpointData);
         };
-        notifier->notifier_.subscribe(checkpointDataBroadcast);
+        notifier->checkpointingNotifications_.subscribe(checkpointDataBroadcast);
     }
 
     //! From IMDModule
index 6388c181a07413b9ebfcd43e4407ff7566158aac..86fa1695f0f5709838395c15d7060529a91ee69e 100644 (file)
@@ -2089,7 +2089,7 @@ static void do_cpt_mdmodules(int                           fileVersion,
         gmx::MdModulesCheckpointReadingDataOnMaster mdModuleCheckpointReadingDataOnMaster = {
             mdModuleCheckpointParameterTree, fileVersion
         };
-        mdModulesNotifier.notifier_.notify(mdModuleCheckpointReadingDataOnMaster);
+        mdModulesNotifier.checkpointingNotifications_.notify(mdModuleCheckpointReadingDataOnMaster);
     }
 }
 
@@ -2377,7 +2377,7 @@ void write_checkpoint(const char*                   fn,
         gmx::KeyValueTreeBuilder          builder;
         gmx::MdModulesWriteCheckpointData mdModulesWriteCheckpoint = { builder.rootObject(),
                                                                        headerContents.file_version };
-        mdModulesNotifier.notifier_.notify(mdModulesWriteCheckpoint);
+        mdModulesNotifier.checkpointingNotifications_.notify(mdModulesWriteCheckpoint);
         auto                     tree = builder.build();
         gmx::FileIOXdrSerializer serializer(fp);
         gmx::serializeKeyValueTree(tree, &serializer);
@@ -2819,7 +2819,7 @@ void load_checkpoint(const char*                   fn,
     {
         gmx_bcast(sizeof(headerContents.step), &headerContents.step, cr);
         gmx::MdModulesCheckpointReadingBroadcast broadcastCheckPointData = { *cr, headerContents.file_version };
-        mdModulesNotifier.notifier_.notify(broadcastCheckPointData);
+        mdModulesNotifier.checkpointingNotifications_.notify(broadcastCheckPointData);
     }
     ir->bContinuation = TRUE;
     // TODO Should the following condition be <=? Currently if you
index eb09d636e817548011cdd7cc216d05d5600a3f24..29ac09631610d8b4b04f3f547ead77fa79853b51 100644 (file)
@@ -2394,7 +2394,7 @@ int gmx_grompp(int argc, char* argv[])
 
     {
         gmx::KeyValueTreeBuilder internalParameterBuilder;
-        mdModules.notifier().notifier_.notify(internalParameterBuilder.rootObject());
+        mdModules.notifier().preProcessingNotifications_.notify(internalParameterBuilder.rootObject());
         ir->internalParameters =
                 std::make_unique<gmx::KeyValueTreeObject>(internalParameterBuilder.build());
     }
index 79a79bbf75c76de6782c59e8cbf9c94576aeb2fc..8f9cc8112f0ae6b66500126b0d81be294f61b4fd 100644 (file)
@@ -562,7 +562,7 @@ void check_ir(const char*                   mdparin,
         // Inquire all MdModules, if their parameters match with the energy
         // calculation frequency
         gmx::EnergyCalculationFrequencyErrors energyCalculationFrequencyErrors(ir->nstcalcenergy);
-        mdModulesNotifier.notifier_.notify(&energyCalculationFrequencyErrors);
+        mdModulesNotifier.preProcessingNotifications_.notify(&energyCalculationFrequencyErrors);
 
         // Emit all errors from the energy calculation frequency checks
         for (const std::string& energyFrequencyErrorMessage :
@@ -3615,7 +3615,7 @@ void do_index(const char*                   mdparin,
 
     gmx::IndexGroupsAndNames defaultIndexGroupsAndNames(
             *defaultIndexGroups, gmx::arrayRefFromArray(gnames, defaultIndexGroups->nr));
-    notifier.notifier_.notify(defaultIndexGroupsAndNames);
+    notifier.preProcessingNotifications_.notify(defaultIndexGroupsAndNames);
 
     auto accelerations          = gmx::splitString(is->acc);
     auto accelerationGroupNames = gmx::splitString(is->accgrps);
index 01ef7b2cc9860e27ea03a6fae670a8d4eeae7274..345b8438fd9ac15b72e48125aaf50cd4ba1e94ef 100644 (file)
@@ -242,7 +242,7 @@ EnergyOutput::EnergyOutput(ener_file*               fp_ene,
     bEner_[F_COM_PULL]   = ((ir->bPull && pull_have_potential(pull_work)) || ir->bRot);
 
     MdModulesEnergyOutputToDensityFittingRequestChecker mdModulesAddOutputToDensityFittingFieldRequest;
-    mdModulesNotifier.notifier_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
+    mdModulesNotifier.simulationSetupNotifications_.notify(&mdModulesAddOutputToDensityFittingFieldRequest);
 
     bEner_[F_DENSITYFITTING] = mdModulesAddOutputToDensityFittingFieldRequest.energyOutputToDensityFitting_;
 
index f68714ce70b21e40f0aa19d13d2b6cc91d3d76af..5c4ed8a49c0662b730ac608d5af6b9a04926c077 100644 (file)
@@ -917,7 +917,7 @@ int Mdrunner::mdrunner()
 
     // TODO: Error handling
     mdModules_->assignOptionsToModules(*inputrec->params, nullptr);
-    const auto& mdModulesNotifier = mdModules_->notifier().notifier_;
+    const auto& mdModulesNotifier = mdModules_->notifier().simulationSetupNotifications_;
 
     if (inputrec->internalParameters != nullptr)
     {
@@ -1335,7 +1335,7 @@ int Mdrunner::mdrunner()
     {
         mdModulesNotifier.notify(*cr);
         mdModulesNotifier.notify(&atomSets);
-        mdModulesNotifier.notify(PeriodicBoundaryConditionType{ inputrec->pbcType });
+        mdModulesNotifier.notify(inputrec->pbcType);
         mdModulesNotifier.notify(SimulationTimeStep{ inputrec->delta_t });
         /* Initiate forcerecord */
         fr                 = new t_forcerec;
index 891ec575c24164bcfdf90574838f89cc4663e9db..1f4b593651892ee35a478c76eee19b0ece79e5f4 100644 (file)
@@ -202,10 +202,6 @@ class IndexGroupsAndNames;
 struct MdModulesCheckpointReadingDataOnMaster;
 struct MdModulesCheckpointReadingBroadcast;
 struct MdModulesWriteCheckpointData;
-struct PeriodicBoundaryConditionType
-{
-    PbcType pbcType;
-};
 
 struct MdModulesEnergyOutputToDensityFittingRequestChecker
 {
@@ -253,21 +249,75 @@ struct SimulationTimeStep
     const double delta_t;
 };
 
+/*! \libinternal
+ * \brief Collection of callbacks to MDModules at differnt run-times.
+ *
+ * MDModules use members of this struct to sign up for callback functionality.
+ *
+ * The members of the struct represent callbacks at these run-times:
+ *
+ *  When pre-processing the simulation data
+ *  When reading and writing check-pointing data
+ *  When setting up simulation after reading in the tpr file
+ *
+ * The template arguments to the members of this struct directly reflect
+ * the callback function signature. Arguments passed as pointers are always
+ * meant to be modified, but never meant to be stored (in line with the policy
+ * everywhere else).
+ */
 struct MdModulesNotifier
 {
-    //! Register callback function types for MdModule
-    registerMdModuleNotification<const t_commrec&,
-                                 EnergyCalculationFrequencyErrors*,
-                                 IndexGroupsAndNames,
-                                 KeyValueTreeObjectBuilder,
-                                 const KeyValueTreeObject&,
+    /*! \brief Pre-processing callback functions.
+     *
+     * EnergyCalculationFrequencyErrors* allows modules to check if they match
+     *                                   their required calculation frequency
+     *                                   and add their error message if needed
+     *                                   to the collected error messages
+     * IndexGroupsAndNames provides modules with atom indices and their names
+     * KeyValueTreeObjectBuilder enables writing of module internal data to
+     *                           .tpr files.
+     */
+    registerMdModuleNotification<EnergyCalculationFrequencyErrors*, IndexGroupsAndNames, KeyValueTreeObjectBuilder>::type preProcessingNotifications_;
+
+    /*! \brief Checkpointing callback functions.
+     *
+     * MdModulesCheckpointReadingDataOnMaster provides modules with their
+     *                                        checkpointed data on the master
+     *                                        node and checkpoint file version
+     * MdModulesCheckpointReadingBroadcast provides modules with a communicator
+     *                                     and the checkpoint file version to
+     *                                     distribute their data
+     * MdModulesWriteCheckpointData provides the modules with a key-value-tree
+     *                              builder to store their checkpoint data and
+     *                              the checkpoint file version
+     */
+    registerMdModuleNotification<MdModulesCheckpointReadingDataOnMaster,
+                                 MdModulesCheckpointReadingBroadcast,
+                                 MdModulesWriteCheckpointData>::type checkpointingNotifications_;
+
+    /*! \brief Callbacks during simulation setup.
+     *
+     * const KeyValueTreeObject& provides modules with the internal data they
+     *                           wrote to .tpr files
+     * LocalAtomSetManager* enables modules to add atom indices to local atom sets
+     *                      to be managed
+     * MdModulesEnergyOutputToDensityFittingRequestChecker* enables modules to
+     *                      report if they want to write their energy output
+     *                      to the density fitting field in the energy files
+     * const PbcType& provides modules with the periodic boundary condition type
+     *                that is used during the simulation
+     * const SimulationTimeStep& provides modules with the simulation time-step
+     *                           that allows them to interconvert between step
+     *                           time information
+     * const t_commrec& provides a communicator to the modules during simulation
+     *                  setup
+     */
+    registerMdModuleNotification<const KeyValueTreeObject&,
                                  LocalAtomSetManager*,
                                  MdModulesEnergyOutputToDensityFittingRequestChecker*,
-                                 MdModulesCheckpointReadingDataOnMaster,
-                                 MdModulesCheckpointReadingBroadcast,
-                                 MdModulesWriteCheckpointData,
-                                 PeriodicBoundaryConditionType,
-                                 const SimulationTimeStep&>::type notifier_;
+                                 const PbcType&,
+                                 const SimulationTimeStep&,
+                                 const t_commrec&>::type simulationSetupNotifications_;
 };
 
 } // namespace gmx