Improve naming and docs.
authorM. Eric Irrgang <ericirrgang@gmail.com>
Tue, 15 Sep 2020 19:18:22 +0000 (22:18 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Fri, 18 Sep 2020 12:42:44 +0000 (12:42 +0000)
Follow up from https://gitlab.com/gromacs/gromacs/-/merge_requests/541

Fixes #3673

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

index 509b009dcb0cf2ee6b447673f3973d8da0c3f344..301a02ffe2e1be16d253cda408c178809087bda3 100644 (file)
@@ -352,12 +352,12 @@ 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;
-    newRunner.stopHandlerBuilder_ = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
-    newRunner.inputHolder_        = inputHolder_;
+    newRunner.libraryWorldCommunicator = MPI_COMM_WORLD;
+    newRunner.simulationCommunicator   = MPI_COMM_WORLD;
+    newRunner.ms                       = ms;
+    newRunner.startingBehavior         = startingBehavior;
+    newRunner.stopHandlerBuilder_      = std::make_unique<StopHandlerBuilder>(*stopHandlerBuilder_);
+    newRunner.inputHolder_             = inputHolder_;
 
     threadMpiMdrunnerAccessBarrier();
 
@@ -399,8 +399,8 @@ void Mdrunner::spawnThreads(int numThreadsToLaunch)
 
     // Give the master thread the newly created valid communicator for
     // the simulation.
-    worldCommunicator = MPI_COMM_WORLD;
-    communicator      = MPI_COMM_WORLD;
+    libraryWorldCommunicator = MPI_COMM_WORLD;
+    simulationCommunicator   = MPI_COMM_WORLD;
     threadMpiMdrunnerAccessBarrier();
 #else
     GMX_UNUSED_VALUE(numThreadsToLaunch);
@@ -763,7 +763,7 @@ int Mdrunner::mdrunner()
     {
         fplog = gmx_fio_getfp(logFileHandle);
     }
-    const bool       isSimulationMasterRank = findIsSimulationMasterRank(ms, communicator);
+    const bool isSimulationMasterRank = findIsSimulationMasterRank(ms, simulationCommunicator);
     gmx::LoggerOwner logOwner(buildLogger(fplog, isSimulationMasterRank));
     gmx::MDLogger    mdlog(logOwner.logger());
 
@@ -775,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(worldCommunicator, gmx_physicalnode_id_hash());
+    PhysicalNodeCommunicator physicalNodeComm(libraryWorldCommunicator, gmx_physicalnode_id_hash());
     hwinfo = gmx_detect_hardware(mdlog, physicalNodeComm);
 
     gmx_print_detected_hardware(fplog, isSimulationMasterRank && isMasterSim(ms), mdlog, hwinfo);
@@ -846,12 +846,13 @@ 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(worldCommunicator, gmx_physicalnode_id_hash());
+        physicalNodeComm =
+                PhysicalNodeCommunicator(libraryWorldCommunicator, gmx_physicalnode_id_hash());
     }
 
-    GMX_RELEASE_ASSERT(ms || communicator == MPI_COMM_WORLD,
-                       "Must have valid world communicator unless running a multi-simulation");
-    CommrecHandle crHandle = init_commrec(communicator);
+    GMX_RELEASE_ASSERT(ms || simulationCommunicator != MPI_COMM_NULL,
+                       "Must have valid communicator unless running a multi-simulation");
+    CommrecHandle crHandle = init_commrec(simulationCommunicator);
     t_commrec*    cr       = crHandle.get();
     GMX_RELEASE_ASSERT(cr != nullptr, "Must have valid commrec");
 
@@ -1185,7 +1186,7 @@ int Mdrunner::mdrunner()
 
     // Produce the task assignment for this rank - done after DD is constructed
     GpuTaskAssignments gpuTaskAssignments = GpuTaskAssignmentsBuilder::build(
-            gpuIdsToUse, userGpuTaskAssignment, *hwinfo, communicator, physicalNodeComm,
+            gpuIdsToUse, userGpuTaskAssignment, *hwinfo, simulationCommunicator, physicalNodeComm,
             nonbondedTarget, pmeTarget, bondedTarget, updateTarget, useGpuForNonbonded,
             useGpuForPme, thisRankHasDuty(cr, DUTY_PP),
             // TODO cr->duty & DUTY_PME should imply that a PME
@@ -1923,13 +1924,13 @@ private:
     int nstlist_ = 0;
 
     //! World communicator, used for hardware detection and task assignment
-    MPI_Comm worldCommunicator_ = MPI_COMM_NULL;
+    MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
 
     //! Multisim communicator handle.
     gmx_multisim_t* multiSimulation_;
 
     //! mdrun communicator
-    MPI_Comm communicator_ = MPI_COMM_NULL;
+    MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
 
     //! Print a warning if any force is larger than this (in kJ/mol nm).
     real forceWarningThreshold_ = -1;
@@ -1981,9 +1982,9 @@ Mdrunner::BuilderImplementation::BuilderImplementation(std::unique_ptr<MDModules
                                                        compat::not_null<SimulationContext*> context) :
     mdModules_(std::move(mdModules))
 {
-    worldCommunicator_ = context->worldCommunicator_;
-    communicator_      = context->simulationCommunicator_;
-    multiSimulation_   = context->multiSimulation_.get();
+    libraryWorldCommunicator_ = context->libraryWorldCommunicator_;
+    simulationCommunicator_   = context->simulationCommunicator_;
+    multiSimulation_          = context->multiSimulation_.get();
 }
 
 Mdrunner::BuilderImplementation::~BuilderImplementation() = default;
@@ -2033,9 +2034,9 @@ Mdrunner Mdrunner::BuilderImplementation::build()
 
     newRunner.filenames = filenames_;
 
-    newRunner.worldCommunicator = worldCommunicator_;
+    newRunner.libraryWorldCommunicator = libraryWorldCommunicator_;
 
-    newRunner.communicator = communicator_;
+    newRunner.simulationCommunicator = simulationCommunicator_;
 
     // nullptr is a valid value for the multisim handle
     newRunner.ms = multiSimulation_;
index 9f5790b9004bafa147b90039c4ac5d8c1ef3868d..2dd38531bb86ce11037c759501ff5be5dab30a5d 100644 (file)
@@ -266,14 +266,14 @@ private:
      * 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;
+    MPI_Comm libraryWorldCommunicator = MPI_COMM_NULL;
 
-    /*! \brief Non-owning handle to communication data structure.
+    /*! \brief Non-owning handle to communication data structure for the current simulation.
      *
      * 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 communicator = MPI_COMM_NULL;
+    MPI_Comm simulationCommunicator = MPI_COMM_NULL;
 
     //! \brief Non-owning handle to multi-simulation handler.
     gmx_multisim_t* ms = nullptr;
index 6beb0e9b0d42c1aa0f42460d7f80763063c2b620..229da88281b78f7b9ffebf4f4522ed0f2695f3cb 100644 (file)
 namespace gmx
 {
 //! \cond libapi
-SimulationContext::SimulationContext(MPI_Comm                          worldCommunicator,
-                                     const ArrayRef<const std::string> multiSimDirectoryNames) :
-    worldCommunicator_(worldCommunicator)
+SimulationContext::SimulationContext(MPI_Comm                    communicator,
+                                     ArrayRef<const std::string> multiSimDirectoryNames) :
+    libraryWorldCommunicator_(communicator)
 {
-    GMX_RELEASE_ASSERT((GMX_LIB_MPI && (worldCommunicator != MPI_COMM_NULL))
-                               || (!GMX_LIB_MPI && (worldCommunicator == MPI_COMM_NULL)),
+    GMX_RELEASE_ASSERT((GMX_LIB_MPI && (communicator != MPI_COMM_NULL))
+                               || (!GMX_LIB_MPI && (communicator == MPI_COMM_NULL)),
                        "With real MPI, a non-null communicator is required. "
                        "Without it, the communicator must be null.");
     if (multiSimDirectoryNames.empty())
     {
-        simulationCommunicator_ = worldCommunicator;
+        simulationCommunicator_ = communicator;
     }
     else
     {
-        multiSimulation_ = buildMultiSimulation(worldCommunicator, multiSimDirectoryNames);
+        multiSimulation_ = buildMultiSimulation(communicator, multiSimDirectoryNames);
         // Use the communicator resulting from the split for the multi-simulation.
         simulationCommunicator_ = multiSimulation_->simulationComm_;
     }
index 717b4f30a9f943488c153184fe6bdcdcee7549d2..0fb81b5865f9f99b613d4986f3d901ecd78a8f66 100644 (file)
@@ -108,45 +108,67 @@ public:
     /*!
      * \brief Construct
      *
-     * \param worldCommunicator       MPI communicator for this simulation context
+     * \param communicator Communicator for all collaborating simulation processes.
      * \param multiSimDirectoryNames  Names of any directories used with -multidir
+     *
+     * Caller is responsible for keeping *communicator* valid for the life of SimulationContext,
+     * and then freeing the communicator, if appropriate.
+     *
+     * With an external MPI library (GMX_LIB_MPI==1), the client promises that
+     * gmx::init() has called MPI_Init and the provided communicator is valid to use.
+     * This communicator is "borrowed" (not duplicated) from the caller.
+     * Additional communicators may be split from the provided communicator
+     * during the life of the SimulationContext or its consumers.
+     *
+     * With thread-MPI, *communicator* must be MPI_COMM_NULL.
+     * 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.)
+     *
+     * \todo Refine distribution of environment management.
+     *       There should be a context level at which only the current simulator directory matters,
+     *       and a level above which encapsulates multisim details in a specialized type.
      */
-    explicit SimulationContext(MPI_Comm                    worldCommunicator,
-                               ArrayRef<const std::string> multiSimDirectoryNames);
+    explicit SimulationContext(MPI_Comm communicator, ArrayRef<const std::string> multiSimDirectoryNames);
 
     /*!
-     * \brief MPI communicator object for this GROMACS instance.
+     * \brief MPI resources for the entire simulation context.
      *
-     * 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 an external MPI library (GMX_LIB_MPI==1),
+     * gmx::init() has called MPI_Init and the provided communicator is valid to use.
+     * The communicator is "borrowed" (not duplicated) from the caller.
      *
-     * With thread-MPI in both cases, the communicator is set up later
+     * With thread-MPI, 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;
+    MPI_Comm libraryWorldCommunicator_ = MPI_COMM_NULL;
 
     /*!
-     * \brief MPI communicator object for this simulation object.
+     * \brief MPI communicator object for this simulation.
      *
-     * 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 an external MPI library (GMX_LIB_MPI==1),
+     * gmx::init() has called MPI_Init and the provided communicator is valid to use.
+     * The communicator is "borrowed" (not duplicated) from the caller.
      *
-     * With thread-MPI in both cases, the communicator is set up later
+     * With thread-MPI, 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.)
+     * ranks.
      */
     MPI_Comm simulationCommunicator_ = MPI_COMM_NULL;
 
     /*!
      * \brief Multi-sim handler (if required by e.g. gmx mdrun
      * -multidir; only supported with real MPI)
+     *
+     * 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) communicator splitting
+     *       can be undertaken during multi-sim setup (acquisition of the multisim resource).
+     *
+     * Multi-simulation is not supported with thread-MPI.
      */
     std::unique_ptr<gmx_multisim_t> multiSimulation_;
 };