SimulatorBuilder.add() for StopHandlerBuilder.
authorM. Eric Irrgang <ericirrgang@gmail.com>
Mon, 22 Jun 2020 13:27:52 +0000 (16:27 +0300)
committerJoe Jordan <ejjordan12@gmail.com>
Fri, 26 Jun 2020 10:54:08 +0000 (10:54 +0000)
Refs #3567

src/gromacs/mdrun/runner.cpp
src/gromacs/mdrun/simulatorbuilder.cpp
src/gromacs/mdrun/simulatorbuilder.h

index 5a2e7a97fa915ae91851dd4ce791a222d47617eb..1378888adfe2c18b97149725caea18e3b38c47dc 100644 (file)
@@ -1641,6 +1641,8 @@ int Mdrunner::mdrunner()
         GMX_ASSERT(stopHandlerBuilder_, "Runner must provide StopHandlerBuilder to simulator.");
         SimulatorBuilder simulatorBuilder;
 
+        simulatorBuilder.add(std::move(stopHandlerBuilder_));
+
         // build and run simulator object based on user-input
         auto simulator = simulatorBuilder.build(
                 useModularSimulator, fplog, cr, ms, mdlog, static_cast<int>(filenames.size()),
@@ -1649,7 +1651,7 @@ int Mdrunner::mdrunner()
                 mdModules_->outputProvider(), mdModules_->notifier(), inputrec, imdSession.get(),
                 pull_work, swap, &mtop, globalState.get(), &observablesHistory, mdAtoms.get(),
                 &nrnb, wcycle, fr, &enerd, &ekind, &runScheduleWork, replExParams, membed,
-                walltime_accounting, std::move(stopHandlerBuilder_), doRerun);
+                walltime_accounting, doRerun);
         simulator->run();
 
         if (fr->pmePpCommGpu)
index 0013091df029ab2a544d26480a612532173d8ab4..7f52cc13295e6eca99ab90299861e43713036a31 100644 (file)
@@ -89,9 +89,13 @@ std::unique_ptr<ISimulator> SimulatorBuilder::build(bool                     use
                                                     const ReplicaExchangeParameters& replExParams,
                                                     gmx_membed_t*                    membed,
                                                     gmx_walltime_accounting* walltime_accounting,
-                                                    std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder,
-                                                    bool                                doRerun)
+                                                    bool                     doRerun)
 {
+    if (!stopHandlerBuilder_)
+    {
+        throw APIError("You must add a StopHandlerBuilder before calling build().");
+    }
+
     if (useModularSimulator)
     {
         // NOLINTNEXTLINE(modernize-make-unique): make_unique does not work with private constructor
@@ -100,7 +104,7 @@ std::unique_ptr<ISimulator> SimulatorBuilder::build(bool                     use
                 constr, enforcedRotation, deform, outputProvider, mdModulesNotifier, inputrec,
                 imdSession, pull_work, swap, top_global, state_global, observablesHistory, mdAtoms,
                 nrnb, wcycle, fr, enerd, ekind, runScheduleWork, replExParams, membed,
-                walltime_accounting, std::move(stopHandlerBuilder), doRerun));
+                walltime_accounting, std::move(stopHandlerBuilder_), doRerun));
     }
     // NOLINTNEXTLINE(modernize-make-unique): make_unique does not work with private constructor
     return std::unique_ptr<LegacySimulator>(new LegacySimulator(
@@ -108,7 +112,7 @@ std::unique_ptr<ISimulator> SimulatorBuilder::build(bool                     use
             enforcedRotation, deform, outputProvider, mdModulesNotifier, inputrec, imdSession,
             pull_work, swap, top_global, state_global, observablesHistory, mdAtoms, nrnb, wcycle,
             fr, enerd, ekind, runScheduleWork, replExParams, membed, walltime_accounting,
-            std::move(stopHandlerBuilder), doRerun));
+            std::move(stopHandlerBuilder_), doRerun));
 }
 
 } // namespace gmx
index 1e51550ddb7b39aa6eac94282d1dfae515877a71..6a12a3695719cdff68f6e9de84be38bbd7775124 100644 (file)
@@ -86,58 +86,65 @@ struct MdrunOptions;
  * \brief Class preparing the creation of Simulator objects
  *
  * Objects of this class build Simulator objects, which in turn are used to
- * run molecular simulations. Currently, this only has a single public
- * `build` function which takes all arguments needed to build the
- * `LegacySimulator`.
+ * run molecular simulations.
  */
 class SimulatorBuilder
 {
 public:
+    void add(std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder)
+    {
+        stopHandlerBuilder_ = std::move(stopHandlerBuilder);
+    }
+
     /*! \brief Build a Simulator object based on input data
      *
      * Return a pointer to a simulation object. The use of a parameter
      * pack insulates the builder from changes to the arguments of the
      * Simulator objects.
      *
-     * @return  Unique pointer to a Simulator object
+     * \throws gmx::APIError if expected set-up methods have not been called before build()
+     *
+     * \return  Unique pointer to a Simulator object
      */
-    std::unique_ptr<ISimulator> build(bool                                useModularSimulator,
-                                      FILE*                               fplog,
-                                      t_commrec*                          cr,
-                                      const gmx_multisim_t*               ms,
-                                      const MDLogger&                     mdlog,
-                                      int                                 nfile,
-                                      const t_filenm*                     fnm,
-                                      const gmx_output_env_t*             oenv,
-                                      const MdrunOptions&                 mdrunOptions,
-                                      StartingBehavior                    startingBehavior,
-                                      VirtualSitesHandler*                vsite,
-                                      Constraints*                        constr,
-                                      gmx_enfrot*                         enforcedRotation,
-                                      BoxDeformation*                     deform,
-                                      IMDOutputProvider*                  outputProvider,
-                                      const MdModulesNotifier&            mdModulesNotifier,
-                                      t_inputrec*                         inputrec,
-                                      ImdSession*                         imdSession,
-                                      pull_t*                             pull_work,
-                                      t_swap*                             swap,
-                                      gmx_mtop_t*                         top_global,
-                                      t_state*                            state_global,
-                                      ObservablesHistory*                 observablesHistory,
-                                      MDAtoms*                            mdAtoms,
-                                      t_nrnb*                             nrnb,
-                                      gmx_wallcycle*                      wcycle,
-                                      t_forcerec*                         fr,
-                                      gmx_enerdata_t*                     enerd,
-                                      gmx_ekindata_t*                     ekind,
-                                      MdrunScheduleWorkload*              runScheduleWork,
-                                      const ReplicaExchangeParameters&    replExParams,
-                                      gmx_membed_t*                       membed,
-                                      gmx_walltime_accounting*            walltime_accounting,
-                                      std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder,
-                                      bool                                doRerun);
-};
+    std::unique_ptr<ISimulator> build(bool                     useModularSimulator,
+                                      FILE*                    fplog,
+                                      t_commrec*               cr,
+                                      const gmx_multisim_t*    ms,
+                                      const MDLogger&          mdlog,
+                                      int                      nfile,
+                                      const t_filenm*          fnm,
+                                      const gmx_output_env_t*  oenv,
+                                      const MdrunOptions&      mdrunOptions,
+                                      StartingBehavior         startingBehavior,
+                                      VirtualSitesHandler*     vsite,
+                                      Constraints*             constr,
+                                      gmx_enfrot*              enforcedRotation,
+                                      BoxDeformation*          deform,
+                                      IMDOutputProvider*       outputProvider,
+                                      const MdModulesNotifier& mdModulesNotifier,
+                                      t_inputrec*              inputrec,
+                                      ImdSession*              imdSession,
+                                      pull_t*                  pull_work,
+                                      t_swap*                  swap,
+                                      gmx_mtop_t*              top_global,
 
+                                      t_state*                         state_global,
+                                      ObservablesHistory*              observablesHistory,
+                                      MDAtoms*                         mdAtoms,
+                                      t_nrnb*                          nrnb,
+                                      gmx_wallcycle*                   wcycle,
+                                      t_forcerec*                      fr,
+                                      gmx_enerdata_t*                  enerd,
+                                      gmx_ekindata_t*                  ekind,
+                                      MdrunScheduleWorkload*           runScheduleWork,
+                                      const ReplicaExchangeParameters& replExParams,
+                                      gmx_membed_t*                    membed,
+                                      gmx_walltime_accounting*         walltime_accounting,
+                                      bool                             doRerun);
+
+private:
+    std::unique_ptr<StopHandlerBuilder> stopHandlerBuilder_;
+};
 } // namespace gmx
 
 #endif // GMX_MDRUN_SIMULATORBUILDER_SIMULATORBUILDER_H