Allow MPI-enabled clients to provide the library communicator.
authorM. Eric Irrgang <ericirrgang@gmail.com>
Mon, 28 Sep 2020 12:55:50 +0000 (15:55 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 29 Sep 2020 06:09:40 +0000 (06:09 +0000)
Provide library helpers and template headers to allow clients to
createContext() with explitly assigned resources.

Relates to #3644

api/gmxapi/cpp/CMakeLists.txt
api/gmxapi/cpp/context.cpp
api/gmxapi/cpp/context_impl.h
api/gmxapi/cpp/multiprocessingresources.cpp [new file with mode: 0644]
api/gmxapi/include/gmxapi/context.h
python_packaging/src/gmxapi/pycontext.cpp
src/api/cpp/tests/CMakeLists.txt
src/api/cpp/tests/context.cpp [new file with mode: 0644]
src/api/cpp/tests/restraint.cpp
src/api/cpp/tests/runner.cpp
src/api/cpp/tests/stopsignaler.cpp

index 129ffb69fcb91f30b96b3267b7bcb3ed64a256e2..f9d5a983738b291cb638f7b3715a68afd95fca38 100644 (file)
@@ -35,6 +35,7 @@
 # This list file provides the Gromacs::gmxapi cmake target.
 
 target_sources(gmxapi PRIVATE
+               multiprocessingresources.cpp
                context.cpp
                exceptions.cpp
                gmxapi.cpp
index befb90ad59b95027e6d6d2c2505f4225f8350f16..e2efd4f812a88088dc83cdad4a69e923409f41eb 100644 (file)
 namespace gmxapi
 {
 
-MpiContextManager::MpiContextManager() :
-    communicator_(std::make_unique<MPI_Comm>(GMX_LIB_MPI ? MPI_COMM_WORLD : MPI_COMM_NULL))
+MpiContextManager::MpiContextManager(MPI_Comm communicator) :
+    communicator_(std::make_unique<MPI_Comm>(communicator))
 {
+    GMX_ASSERT(*communicator_ == MPI_COMM_NULL ? !GMX_LIB_MPI : GMX_LIB_MPI,
+               "Invalid communicator value for this library configuration.");
     // Safely increments the GROMACS MPI initialization counter after checking
     // whether the MPI library is already initialized. After this call, MPI_Init
     // or MPI_Init_thread has been called exactly once.
@@ -108,6 +110,11 @@ MpiContextManager::~MpiContextManager()
     }
 }
 
+MpiContextManager::MpiContextManager() :
+    MpiContextManager(GMX_LIB_MPI ? MPI_COMM_WORLD : MPI_COMM_NULL)
+{
+}
+
 MPI_Comm MpiContextManager::communicator() const
 {
     if (!communicator_)
@@ -119,6 +126,21 @@ MPI_Comm MpiContextManager::communicator() const
 
 ContextImpl::~ContextImpl() = default;
 
+Context createContext(const ResourceAssignment& resources)
+{
+    CommHandle handle;
+    resources.applyCommunicator(&handle);
+    if (handle.communicator == MPI_COMM_NULL)
+    {
+        throw UsageError("MPI-enabled Simulator contexts require a valid communicator.");
+    }
+    auto contextmanager = MpiContextManager(handle.communicator);
+    auto impl           = ContextImpl::create(std::move(contextmanager));
+    GMX_ASSERT(impl, "ContextImpl creation method should not be able to return null.");
+    auto context = Context(impl);
+    return context;
+}
+
 Context createContext()
 {
     MpiContextManager contextmanager;
@@ -239,9 +261,9 @@ std::shared_ptr<Session> ContextImpl::launch(const Workflow& work)
         LogFilePtr       logFileGuard     = nullptr;
         gmx_multisim_t*  ms               = simulationContext.multiSimulation_.get();
         std::tie(startingBehavior, logFileGuard) =
-                handleRestart(findIsSimulationMasterRank(ms, communicator), communicator, ms,
-                              options_.mdrunOptions.appendingBehavior, ssize(options_.filenames),
-                              options_.filenames.data());
+                handleRestart(findIsSimulationMasterRank(ms, simulationContext.simulationCommunicator_),
+                              communicator, ms, options_.mdrunOptions.appendingBehavior,
+                              ssize(options_.filenames), options_.filenames.data());
 
         auto builder = MdrunnerBuilder(std::move(mdModules),
                                        compat::not_null<SimulationContext*>(&simulationContext));
@@ -299,12 +321,6 @@ std::shared_ptr<Session> ContextImpl::launch(const Workflow& work)
     return launchedSession;
 }
 
-// As of gmxapi 0.0.3 there is only one Context type
-Context::Context() : Context{ ContextImpl::create(MpiContextManager()) }
-{
-    GMX_ASSERT(impl_, "Context requires a non-null implementation member.");
-}
-
 std::shared_ptr<Session> Context::launch(const Workflow& work)
 {
     return impl_->launch(work);
@@ -312,7 +328,10 @@ std::shared_ptr<Session> Context::launch(const Workflow& work)
 
 Context::Context(std::shared_ptr<ContextImpl> impl) : impl_{ std::move(impl) }
 {
-    GMX_ASSERT(impl_, "Context requires a non-null implementation member.");
+    if (!impl_)
+    {
+        throw UsageError("Context requires a non-null implementation member.");
+    }
 }
 
 void Context::setMDArgs(const MDArgs& mdArgs)
index 7d912f4b81b54444a455e963187be27b15fc180a..2f57c1404ef3014bdf1c94f0f0b0d067e20a2bcf 100644 (file)
@@ -90,6 +90,20 @@ public:
      */
     MpiContextManager();
 
+    /*!
+     * \brief Borrow a communicator and ensure that MPI is initialized, if applicable.
+     *
+     * \param communicator Optional communicator representing a client-managed MPI environment.
+     *
+     *
+     *
+     * Note that the communicator must be MPI_COMM_NULL if and only if GROMACS was built without an
+     * external MPI library.
+     *
+     * \todo (#3650?) Decide whether to exclude this from tMPI environments or find a sensible invariant.
+     */
+    explicit MpiContextManager(MPI_Comm communicator);
+
     ~MpiContextManager();
 
     /*!
diff --git a/api/gmxapi/cpp/multiprocessingresources.cpp b/api/gmxapi/cpp/multiprocessingresources.cpp
new file mode 100644 (file)
index 0000000..2c2f107
--- /dev/null
@@ -0,0 +1,72 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2020, 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.
+ */
+
+/*! \file
+ * \brief Provide a bridge to communication resources appropriate for the library.
+ *
+ * Define the helper functions that the library provides for to help the client
+ * implement the interfaces. (This is dependent on the library configuration.)
+ *
+ * \author "M. Eric Irrgang <ericirrgang@gmail.com"
+ */
+
+#include "gmxapi/mpi/multiprocessingresources.h"
+
+#include "config.h"
+
+#include "gromacs/utility/gmxmpi.h"
+
+#include "context_impl.h"
+
+namespace gmxapi
+{
+
+ResourceAssignment::~ResourceAssignment() = default;
+
+// Base implementation is only overridden by client-provided code in certain
+// combinations of library and client build configurations.
+void ResourceAssignment::applyCommunicator(CommHandle* dst) const
+{
+    dst->communicator = MPI_COMM_NULL;
+}
+
+#if GMX_LIB_MPI
+void offerComm(MPI_Comm src, CommHandle* dst)
+{
+    dst->communicator = src;
+}
+#endif
+
+} // end namespace gmxapi
index ec434059af419517dcccae5dc632a11ce2bd8a8d..296b24d7e3ec5007cb78373d3a93590828f02cdb 100644 (file)
@@ -123,7 +123,7 @@ public:
     /*!
      * \brief Get a handle to a new default context object.
      */
-    Context();
+    Context() = delete;
     ~Context();
 
     /*!
index 797c7a2c95256967d610fbc0691aba6e8e8be690..5b9878efbb3d6d7ea8064cfbf3d6750a43d6cf60 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019, by the GROMACS development team, led by
+ * Copyright (c) 2019,2020, 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.
@@ -74,7 +74,7 @@ std::shared_ptr<gmxapi::Context> PyContext::get() const
 }
 
 PyContext::PyContext() :
-    context_{ std::make_shared<gmxapi::Context>() },
+    context_{ std::make_shared<gmxapi::Context>(gmxapi::createContext()) },
     workNodes_{ std::make_shared<gmxapi::MDWorkSpec>() }
 {
     assert(context_);
index 893e736182b9c194b6ffd86d2b43ea9af63f7855..66380fd929fb2cb3d19111f32fcc0d7b19e87823 100644 (file)
@@ -51,6 +51,7 @@
 add_library(gmxapi-testsupport INTERFACE)
 target_include_directories(gmxapi-testsupport INTERFACE ${CMAKE_CURRENT_SOURCE_DIR})
 
+# TODO: Test tMPI build in a MPI-enabled client context.
 gmx_add_gtest_executable(gmxapi-test
     CPP_SOURCE_FILES
         restraint.cpp
@@ -84,6 +85,7 @@ set_tests_properties(GmxapiExternalInterfaceTests PROPERTIES
 if (GMX_MPI)
     gmx_add_gtest_executable(gmxapi-mpi-test MPI
                              CPP_SOURCE_FILES
+                             context.cpp
                              restraint.cpp
                              runner.cpp
                              status.cpp
diff --git a/src/api/cpp/tests/context.cpp b/src/api/cpp/tests/context.cpp
new file mode 100644 (file)
index 0000000..e873871
--- /dev/null
@@ -0,0 +1,108 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018,2019,2020, 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.
+ */
+
+/*! \file
+ * \brief Test gmxapi::Context
+ *
+ * Provides additional test coverage of template headers only used by client code.
+ */
+
+#include "config.h"
+
+#include <gtest/gtest.h>
+
+#include "gromacs/utility/gmxmpi.h"
+
+#include "gmxapi/context.h"
+#include "gmxapi/mpi/gmxapi_mpi.h"
+
+
+namespace gmxapi
+{
+
+namespace testing
+{
+
+namespace
+{
+
+TEST(GmxApiMpiTest, AllContext)
+{
+    // Default Implicit COMM_WORLD for MPI builds.
+    auto context = createContext();
+}
+
+#if GMX_LIB_MPI
+TEST(GmxApiMpiTest, NullContext)
+{
+    // Explicit COMM_NULL is not supported.
+    EXPECT_ANY_THROW(assignResource(MPI_COMM_NULL));
+}
+
+TEST(GmxApiMpiTest, MpiWorldContext)
+{
+    // Note that this test is only compiled when GMX_MPI is enabled for the
+    // build tree, so we cannot unit test the behavior of non-MPI GROMACS
+    // provided with MPI-enabled Context. For that, we defer to the Python
+    // package testing.
+    // Note also that the code should look the same for tMPI or regular MPI.
+
+    // Explicit COMM_WORLD.
+    auto resources = assignResource(MPI_COMM_WORLD);
+    EXPECT_TRUE(resources->size() != 0);
+
+    // Store the rank for debugging convenience.
+    [[maybe_unused]] auto rank = resources->rank();
+
+    auto context = createContext(*resources);
+}
+
+TEST(GmxApiMpiTest, MpiSplitContext)
+{
+    // Explicit sub-communicator.
+    MPI_Comm communicator = MPI_COMM_NULL;
+    int      rank{ 0 };
+    MPI_Comm_rank(MPI_COMM_WORLD, &rank);
+    // Run each rank as a separate ensemble member.
+    MPI_Comm_split(MPI_COMM_WORLD, rank, rank, &communicator);
+    auto context = createContext(*assignResource(communicator));
+}
+#endif
+
+} // end anonymous namespace
+
+} // end namespace testing
+
+} // end namespace gmxapi
index a992527aef73f8c723dcec42279f0f8d360fae1c..11d001db648d6bcf36769e8b4bf9adecbaaacbeb 100644 (file)
@@ -123,7 +123,7 @@ TEST_F(GmxApiTest, ApiRunnerRestrainedMD)
     auto system = gmxapi::fromTprFile(runner_.tprFileName_);
 
     {
-        auto           context = std::make_shared<gmxapi::Context>();
+        auto           context = std::make_shared<gmxapi::Context>(gmxapi::createContext());
         gmxapi::MDArgs args    = makeMdArgs();
 
         context->setMDArgs(args);
index b03903b65c4221d33fe9d2d185f664e4e60bdbd4..d96760836ee2d1906f67d89bc365d8201213f4e7 100644 (file)
@@ -61,7 +61,7 @@ TEST_F(GmxApiTest, RunnerBasicMD)
     auto system = gmxapi::fromTprFile(runner_.tprFileName_);
 
     {
-        auto           context = std::make_shared<gmxapi::Context>();
+        auto           context = std::make_shared<gmxapi::Context>(gmxapi::createContext());
         gmxapi::MDArgs args    = makeMdArgs();
         // TODO the command line arguments should be set through the
         // usual command line options settings for the tests
@@ -82,7 +82,7 @@ TEST_F(GmxApiTest, RunnerBasicMD)
  */
 TEST_F(GmxApiTest, RunnerReinitialize)
 {
-    auto           context = std::make_shared<gmxapi::Context>();
+    auto           context = std::make_shared<gmxapi::Context>(gmxapi::createContext());
     gmxapi::MDArgs args    = makeMdArgs();
 
     makeTprFile(20);
@@ -140,7 +140,7 @@ TEST_F(GmxApiTest, RunnerContinuedMD)
     auto system = gmxapi::fromTprFile(runner_.tprFileName_);
 
     {
-        auto context = std::make_shared<gmxapi::Context>();
+        auto context = std::make_shared<gmxapi::Context>(gmxapi::createContext());
 
         {
             EXPECT_TRUE(context != nullptr);
index 29d13e6b838bb35fade5f4c1494e1c300159fbe7..5cfcf21d224ded061e01dbbf4982bddc6c9b4df4 100644 (file)
@@ -180,7 +180,7 @@ TEST_F(GmxApiTest, ApiRunnerStopSignalClient)
     const int nsteps = 1;
     makeTprFile(nsteps);
     auto system  = gmxapi::fromTprFile(runner_.tprFileName_);
-    auto context = std::make_shared<gmxapi::Context>();
+    auto context = std::make_shared<gmxapi::Context>(gmxapi::createContext());
 
     // Check assumptions about basic simulation behavior.
     {