Fix task assignment with multiruns
authorBerk Hess <hess@kth.se>
Tue, 15 Sep 2020 10:17:46 +0000 (10:17 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 15 Sep 2020 10:17:46 +0000 (10:17 +0000)
Change e0f481ae broke the physical node hardware detection and task
assignment by passing the simulation instead of the world communicator
to the creation of the physical node communicator. Now the world
communicator is stored and passed.

Fixes #3666

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

index 7f8dbc5a8617b98861cedfd080cc21188a7b03c5..781ccea41c4566bfa13f489ecee1b07dee5c1bf2 100644 (file)
@@ -352,6 +352,7 @@ Mdrunner Mdrunner::cloneOnSpawnedThread() const
     newRunner.pforce          = pforce;
     // Give the spawned thread the newly created valid communicator
     // for the simulation.
+    newRunner.worldCommunicator   = MPI_COMM_WORLD;
     newRunner.communicator        = MPI_COMM_WORLD;
     newRunner.ms                  = ms;
     newRunner.startingBehavior    = startingBehavior;
@@ -398,7 +399,8 @@ void Mdrunner::spawnThreads(int numThreadsToLaunch)
 
     // Give the master thread the newly created valid communicator for
     // the simulation.
-    communicator = MPI_COMM_WORLD;
+    worldCommunicator = MPI_COMM_WORLD;
+    communicator      = MPI_COMM_WORLD;
     threadMpiMdrunnerAccessBarrier();
 #else
     GMX_UNUSED_VALUE(numThreadsToLaunch);
@@ -773,7 +775,7 @@ int Mdrunner::mdrunner()
     // this is expressed, e.g. by expressly running detection only the
     // master rank for thread-MPI, rather than relying on the mutex
     // and reference count.
-    PhysicalNodeCommunicator physicalNodeComm(communicator, gmx_physicalnode_id_hash());
+    PhysicalNodeCommunicator physicalNodeComm(worldCommunicator, gmx_physicalnode_id_hash());
     hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
 
     gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
@@ -844,7 +846,7 @@ int Mdrunner::mdrunner()
         spawnThreads(hw_opt.nthreads_tmpi);
         // The spawned threads enter mdrunner() and execution of
         // master and spawned threads joins at the end of this block.
-        physicalNodeComm = PhysicalNodeCommunicator(communicator, gmx_physicalnode_id_hash());
+        physicalNodeComm = PhysicalNodeCommunicator(worldCommunicator, gmx_physicalnode_id_hash());
     }
 
     GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
@@ -1908,6 +1910,9 @@ private:
     //! Command-line override for the duration of a neighbor list with the Verlet scheme.
     int nstlist_ = 0;
 
+    //! World communicator, used for hardware detection and task assignment
+    MPI_Comm worldCommunicator_ = MPI_COMM_NULL;
+
     //! Multisim communicator handle.
     gmx_multisim_t* multiSimulation_;
 
@@ -1957,8 +1962,9 @@ Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules
                                                        compat::not_null<SimulationContext*> context) :
     mdModules_(std::move(mdModules))
 {
-    communicator_    = context->communicator_;
-    multiSimulation_ = context->multiSimulation_.get();
+    worldCommunicator_ = context->worldCommunicator_;
+    communicator_      = context->simulationCommunicator_;
+    multiSimulation_   = context->multiSimulation_.get();
 }
 
 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
@@ -2008,6 +2014,8 @@ Mdrunner Mdrunner::BuilderImplementation::build()
 
     newRunner.filenames = filenames_;
 
+    newRunner.worldCommunicator = worldCommunicator_;
+
     newRunner.communicator = communicator_;
 
     // nullptr is a valid value for the multisim handle
index 21f8317c790741ad68fd6b97c76310b6e0fcb762..9f5790b9004bafa147b90039c4ac5d8c1ef3868d 100644 (file)
@@ -261,6 +261,13 @@ private:
     //! \brief Non-owning handle to file used for logging.
     t_fileio* logFileHandle = nullptr;
 
+    /*! \brief Non-owning handle to world communication data structure for task assigment.
+     *
+     * With real MPI, gets a value from the SimulationContext
+     * supplied to the MdrunnerBuilder. With thread-MPI gets a
+     * value after threads have been spawned. */
+    MPI_Comm worldCommunicator = MPI_COMM_NULL;
+
     /*! \brief Non-owning handle to communication data structure.
      *
      * With real MPI, gets a value from the SimulationContext
index 0c421213bc0f4b10d3bbbb1f79c15faa5b4938ac..6beb0e9b0d42c1aa0f42460d7f80763063c2b620 100644 (file)
 namespace gmx
 {
 //! \cond libapi
-SimulationContext::SimulationContext(MPI_Comm                          communicator,
+SimulationContext::SimulationContext(MPI_Comm                          worldCommunicator,
                                      const ArrayRef<const std::string> multiSimDirectoryNames) :
-    communicator_(communicator)
+    worldCommunicator_(worldCommunicator)
 {
-    GMX_RELEASE_ASSERT((GMX_LIB_MPI && (communicator != MPI_COMM_NULL))
-                               || (!GMX_LIB_MPI && (communicator == MPI_COMM_NULL)),
+    GMX_RELEASE_ASSERT((GMX_LIB_MPI && (worldCommunicator != MPI_COMM_NULL))
+                               || (!GMX_LIB_MPI && (worldCommunicator == MPI_COMM_NULL)),
                        "With real MPI, a non-null communicator is required. "
                        "Without it, the communicator must be null.");
-    if (!multiSimDirectoryNames.empty())
+    if (multiSimDirectoryNames.empty())
     {
-        multiSimulation_ = buildMultiSimulation(communicator, multiSimDirectoryNames);
+        simulationCommunicator_ = worldCommunicator;
+    }
+    else
+    {
+        multiSimulation_ = buildMultiSimulation(worldCommunicator, multiSimDirectoryNames);
         // Use the communicator resulting from the split for the multi-simulation.
-        communicator_ = multiSimulation_->simulationComm_;
+        simulationCommunicator_ = multiSimulation_->simulationComm_;
     }
 }
 
index 6833d8d231c1c84a11bc4f9d3cbf57bbd2bf06f6..717b4f30a9f943488c153184fe6bdcdcee7549d2 100644 (file)
@@ -108,10 +108,26 @@ public:
     /*!
      * \brief Construct
      *
-     * \param communicator            MPI communicator for this (set of) simulations
+     * \param worldCommunicator       MPI communicator for this simulation context
      * \param multiSimDirectoryNames  Names of any directories used with -multidir
      */
-    explicit SimulationContext(MPI_Comm communicator, ArrayRef<const std::string> multiSimDirectoryNames);
+    explicit SimulationContext(MPI_Comm                    worldCommunicator,
+                               ArrayRef<const std::string> multiSimDirectoryNames);
+
+    /*!
+     * \brief MPI communicator object for this GROMACS instance.
+     *
+     * With real MPI,
+     *   the gmx wrapper binary has called MPI_Init, thus
+     *     MPI_COMM_WORLD is now valid to use, and
+     *   (in future) the gmxapi runner will handle such details
+     *     (e.g. via mpi4py) before creating its SimulationContext.
+     *
+     * With thread-MPI in both cases, the communicator is set up later
+     * during the process of spawning the threads that will be the MPI
+     * ranks. (Multi-simulation is not supported with thread-MPI.)
+     */
+    MPI_Comm worldCommunicator_ = MPI_COMM_NULL;
 
     /*!
      * \brief MPI communicator object for this simulation object.
@@ -121,18 +137,12 @@ public:
      *     MPI_COMM_WORLD is now valid to use, and
      *   (in future) the gmxapi runner will handle such details
      *     (e.g. via mpi4py) before creating its SimulationContext.
-     * In both cases, if a multi-simulation is in use, then its
-     * communicator(s) are found in multiSimulation_. This
-     * communicator is that of all ranks from all simulations, and
-     * will later be split into one for each simulation.
-     * TODO Perhaps (for simplicity) that communicator splitting
-     * task can be undertaken during multi-sim setup.
      *
      * With thread-MPI in both cases, the communicator is set up later
      * during the process of spawning the threads that will be the MPI
      * ranks. (Multi-simulation is not supported with thread-MPI.)
      */
-    MPI_Comm communicator_ = MPI_COMM_NULL;
+    MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
 
     /*!
      * \brief Multi-sim handler (if required by e.g. gmx mdrun