Generalize constraints on MPI rank counts for tests
authorMark Abraham <mark.j.abraham@gmail.com>
Fri, 29 Oct 2021 11:44:45 +0000 (11:44 +0000)
committerAndrey Alekseenko <al42and@gmail.com>
Fri, 29 Oct 2021 11:44:45 +0000 (11:44 +0000)
When the number of ranks is unsuitable the test is skipped.

In some cases minor changes are made so that useful tests run on wider ranges of rank counts.

Closes #3785

docs/dev-manual/testutils.rst
src/gromacs/domdec/tests/haloexchange_mpi.cpp
src/gromacs/fft/tests/fft_mpi.cpp
src/gromacs/mdrunutility/tests/threadaffinity_mpi.cpp
src/gromacs/utility/tests/CMakeLists.txt
src/gromacs/utility/tests/physicalnodecommunicator_mpi.cpp
src/programs/mdrun/tests/multisimtest.cpp
src/programs/mdrun/tests/pmetest.cpp
src/testutils/include/testutils/mpitest.h
src/testutils/testinit.cpp
src/testutils/tests/mpitest.cpp

index 75341a296350daf4ffc5d81851ad58fcf5a7d80b..adfb09c55f67d1e2cd3278e1e5df9ad035c0608d 100644 (file)
@@ -22,9 +22,10 @@ Shared code used to implement the tests is in ``src/external/googletest/`` and
 ``src/testutils/`` (see below).
 
 The tests are built if ``BUILD_TESTING=ON`` (the default) and
-``GMX_BUILD_UNITTESTS=ON`` (the default) in CMake.
-Each module produces a separate unit test binary (:file:`{module}-test`) under
-``bin/``, which can execute all the tests for that module.
+``GMX_BUILD_UNITTESTS=ON`` (the default) in CMake. Each module
+produces at least one separate unit test binary
+(:file:`{module}-test`) under ``bin/``, which can execute tests for
+that module.
 
 The tests can be executed in a few different ways:
 
@@ -47,13 +48,20 @@ The tests can be executed in a few different ways:
   3. If unit tests and/or regression tests are not available, a message is
      printed.
 
-- Directly executing a test binary.  This provides the most useful output for
-  diagnosing failures, and allows debugging test failures.  The output
-  identifies the individual test(s) that fail, and shows the results of all
-  failing assertions.  Some tests also add extra information to failing
-  assertions to make it easier to identify the reason.  It is possible to
-  control which tests are run using command line options.  Execute the binary
-  with ``-h`` to get additional information.
+- The implementation of ``make check`` calls CTest via the ``ctest`` binary
+  to run all the individual test binaries. More fine-grained control is available
+  there, e.g. filtering by test name or label, or increasing verbosity.
+- Directly executing a test binary.  This provides the most useful
+  output for diagnosing failures, and allows debugging test failures.
+  The output identifies the individual test(s) that fail, and shows
+  the results of all failing assertions.  Some tests also add extra
+  information to failing assertions to make it easier to identify the
+  reason. Some tests are skipped because they cannot run with the
+  number of MPI ranks or GPU devices detected.  Explicit information
+  about such cases can be obtained by using the ``-echo-reasons`` flag
+  to the test binary.  It is possible to control which tests are run
+  using command line options.  Execute the binary with ``--help`` to
+  get additional information.
 
 When executed using CTest, the tests produce XML output in
 ``Testing/Temporary/``, containing the result of each test as well as failure
@@ -142,7 +150,7 @@ some pointers to find tests that use certain functionality:
   following it, plus headers required to make them compile.
 - The same file contains also simple tests using the reference framework to
   check line wrapping (the tests for ``gmx::TextLineWrapper``).  The test fixture
-  for these tests is in ``src/testutils/stringtest.h``/``.cpp``.  The string test
+  for these tests is in ``src/testutils/include/testutils/stringtest.h``/``.cpp``.  The string test
   fixture also demonstrates how to add a custom command line option to the
   test binary to influence the test execution.
 - ``src/gromacs/selection/tests/`` contains more complex use of the
@@ -179,3 +187,14 @@ Here are some things to keep in mind when working with the unit tests:
 .. _Google Mock: http://code.google.com/p/googlemock/
 
 .. include:: /fragments/doxygen-links.rst
+
+MPI tests
+---------
+
+If your test makes specific requirements on the number of MPI ranks,
+or needs a communicator as part of its implementation, then there are
+GROMACS-specific extensions that make normal-looking GoogleTests work
+well in these cases. Use ``GMX_TEST_MPI(RankRequirement)`` and declare
+the test with ``gmx_add_mpi_unit_test`` to teach ``CTest`` how to run
+the test regardless of whether the build is with thread-MPI or real
+MPI. See ``src/testutils/include/mpitest.h`` for details.
index 45b8364f8087a5a81373f9748ceb95d53baa446e..fac6cc09abc64ab40da09c11e394436caf7ac87f 100644 (file)
@@ -192,8 +192,9 @@ void define1dRankTopology(gmx_domdec_t* dd)
     int rank;
     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
 
-    dd->neighbor[0][0] = (rank + 1) % 4;
-    dd->neighbor[0][1] = (rank == 0) ? 3 : rank - 1;
+    const int numRanks = getNumberOfTestMpiRanks();
+    dd->neighbor[0][0] = (rank + 1) % numRanks;
+    dd->neighbor[0][1] = (rank == 0) ? (numRanks - 1) : rank - 1;
 }
 
 /*! \brief Define 2D rank topology with 4 MPI tasks
@@ -510,7 +511,7 @@ void checkResults2dHaloWith2PulsesInDim1(const RVec* x, const gmx_domdec_t* dd,
 
 TEST(HaloExchangeTest, Coordinates1dHaloWith1Pulse)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireRankCount<4>);
 
     // Set up atom data
     const int        numHomeAtoms  = 10;
@@ -566,7 +567,7 @@ TEST(HaloExchangeTest, Coordinates1dHaloWith1Pulse)
 
 TEST(HaloExchangeTest, Coordinates1dHaloWith2Pulses)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireRankCount<4>);
 
     // Set up atom data
     const int        numHomeAtoms  = 10;
@@ -623,7 +624,7 @@ TEST(HaloExchangeTest, Coordinates1dHaloWith2Pulses)
 
 TEST(HaloExchangeTest, Coordinates2dHaloWith1PulseInEachDim)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireRankCount<4>);
 
     // Set up atom data
     const int        numHomeAtoms  = 10;
@@ -679,7 +680,7 @@ TEST(HaloExchangeTest, Coordinates2dHaloWith1PulseInEachDim)
 
 TEST(HaloExchangeTest, Coordinates2dHaloWith2PulsesInDim1)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireRankCount<4>);
 
     // Set up atom data
     const int        numHomeAtoms  = 10;
index 17721021961a306aee94a5133664ded2bc443646..25b87cee0512796615a5e135456911f1c31ff4fa 100644 (file)
@@ -216,14 +216,14 @@ public:
 
 TEST_P(GpuFftTest3D, GpuFftDecomposition)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireRankCount<4>);
     GpuFftTestParams params = GetParam();
     runTest(params);
 }
 
 std::vector<GpuFftTestParams> const inputs{
-    { IVec{ 5, 6, 9 }, 4, 1, FftBackend::HeFFTe_CUDA}, // slab decomposition
-    { IVec{ 5, 6, 9 }, 2, 2, FftBackend::HeFFTe_CUDA} // pencil decomposition
+    { IVec{ 5, 6, 9 }, 4, 1, FftBackend::HeFFTe_CUDA }, // slab decomposition
+    { IVec{ 5, 6, 9 }, 2, 2, FftBackend::HeFFTe_CUDA }  // pencil decomposition
 };
 
 INSTANTIATE_TEST_SUITE_P(GpuFft, GpuFftTest3D, ::testing::ValuesIn(inputs));
index 9cd90ccd1e549ed281d59eac3bf2c76b54f1687a..1800a34d219a042602ba3a909f47e19606ef1bd7 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2016,2017,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2016,2017,2019,2020,2021, 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.
 
 #include "threadaffinitytest.h"
 
+namespace gmx
+{
+namespace test
+{
 namespace
 {
 
 using gmx::test::ThreadAffinityTestHelper;
 
+//! Helper to establish test requirements on MPI ranks
+class RequireEvenRankCountWithAtLeastFourRanks
+{
+public:
+    //! Function to require even ranks with at least four ranks
+    static bool conditionSatisfied(const int numRanks)
+    {
+        return (numRanks > 2) && (numRanks % 2 == 0);
+    }
+    //! Text to echo when skipping a test that does not satisfy the requirement
+    inline static const char* s_skipReason = "an even rank count of at least four is required";
+};
+
 TEST(ThreadAffinityMultiRankTest, PinsWholeNode)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     ThreadAffinityTestHelper helper;
-    helper.setLogicalProcessorCount(4);
+    helper.setLogicalProcessorCount(getNumberOfTestMpiRanks());
     helper.expectPinningMessage(false, 1);
     helper.expectAffinitySet(gmx_node_rank());
     helper.setAffinity(1);
@@ -61,11 +78,11 @@ TEST(ThreadAffinityMultiRankTest, PinsWholeNode)
 
 TEST(ThreadAffinityMultiRankTest, PinsWithOffsetAndStride)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::On);
     helper.setOffsetAndStride(1, 2);
-    helper.setLogicalProcessorCount(8);
+    helper.setLogicalProcessorCount(2 * getNumberOfTestMpiRanks());
     helper.expectWarningMatchingRegex("Applying core pinning offset 1");
     helper.expectPinningMessage(true, 2);
     helper.expectAffinitySet(1 + 2 * gmx_node_rank());
@@ -74,7 +91,7 @@ TEST(ThreadAffinityMultiRankTest, PinsWithOffsetAndStride)
 
 TEST(ThreadAffinityMultiRankTest, PinsTwoNodes)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
     helper.setPhysicalNodeId(gmx_node_rank() / 2);
     helper.setLogicalProcessorCount(2);
@@ -85,30 +102,32 @@ TEST(ThreadAffinityMultiRankTest, PinsTwoNodes)
 
 TEST(ThreadAffinityMultiRankTest, DoesNothingWhenDisabled)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::Off);
-    helper.setLogicalProcessorCount(4);
+    helper.setLogicalProcessorCount(getNumberOfTestMpiRanks());
     helper.setAffinity(1);
 }
 
 TEST(ThreadAffinityMultiRankTest, HandlesTooManyThreadsWithAuto)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     ThreadAffinityTestHelper helper;
-    helper.setLogicalProcessorCount(6);
+    const int                threadsPerRank = 2;
+    helper.setLogicalProcessorCount(threadsPerRank * getNumberOfTestMpiRanks() - 1);
     helper.expectWarningMatchingRegex("Oversubscribing the CPU");
-    helper.setAffinity(2);
+    helper.setAffinity(threadsPerRank);
 }
 
 TEST(ThreadAffinityMultiRankTest, HandlesTooManyThreadsWithForce)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     ThreadAffinityTestHelper helper;
+    const int                threadsPerRank = 2;
     helper.setAffinityOption(ThreadAffinity::On);
-    helper.setLogicalProcessorCount(6);
+    helper.setLogicalProcessorCount(threadsPerRank * getNumberOfTestMpiRanks() - 1);
     helper.expectWarningMatchingRegex("Oversubscribing the CPU");
-    helper.setAffinity(2);
+    helper.setAffinity(threadsPerRank);
 }
 
 class ThreadAffinityHeterogeneousNodesTest : public ::testing::Test
@@ -118,11 +137,11 @@ public:
     static int  indexInNode() { return gmx_node_rank() % 2; }
     static bool isMaster() { return gmx_node_rank() == 0; }
 
-    static void setupNodes(ThreadAffinityTestHelper* helper, std::array<int, 2> cores)
+    static void setupNodes(ThreadAffinityTestHelper* helper, int coresOnNodeZero, int coresOnOtherNodes)
     {
         const int node = currentNode();
         helper->setPhysicalNodeId(node);
-        helper->setLogicalProcessorCount(cores[node]);
+        helper->setLogicalProcessorCount(node == 0 ? coresOnNodeZero : coresOnOtherNodes);
     }
     static void expectNodeAffinitySet(ThreadAffinityTestHelper* helper, int node, int core)
     {
@@ -135,10 +154,10 @@ public:
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, PinsOnMasterOnly)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::On);
-    setupNodes(&helper, { { 2, 1 } });
+    setupNodes(&helper, 2, 1);
     helper.expectWarningMatchingRegexIf("Oversubscribing the CPU", isMaster() || currentNode() == 1);
     if (currentNode() == 0)
     {
@@ -150,25 +169,25 @@ TEST_F(ThreadAffinityHeterogeneousNodesTest, PinsOnMasterOnly)
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, PinsOnNonMasterOnly)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::On);
-    setupNodes(&helper, { { 1, 2 } });
+    setupNodes(&helper, 1, 2);
     helper.expectWarningMatchingRegexIf("Oversubscribing the CPU", currentNode() == 0);
-    if (currentNode() == 1)
+    if (currentNode() >= 1)
     {
         helper.expectPinningMessage(false, 1);
+        expectNodeAffinitySet(&helper, currentNode(), indexInNode());
     }
-    expectNodeAffinitySet(&helper, 1, indexInNode());
     helper.setAffinity(1);
 }
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, HandlesUnknownHardwareOnNonMaster)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::On);
-    setupNodes(&helper, { { 2, 0 } });
+    setupNodes(&helper, 2, 0);
     helper.expectWarningMatchingRegexIf("No information on available cores",
                                         isMaster() || currentNode() == 1);
     if (currentNode() == 0)
@@ -181,9 +200,9 @@ TEST_F(ThreadAffinityHeterogeneousNodesTest, HandlesUnknownHardwareOnNonMaster)
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, PinsAutomaticallyOnMasterOnly)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
-    setupNodes(&helper, { { 2, 1 } });
+    setupNodes(&helper, 2, 1);
     helper.expectWarningMatchingRegexIf("Oversubscribing the CPU", isMaster() || currentNode() == 1);
     if (currentNode() == 0)
     {
@@ -195,27 +214,27 @@ TEST_F(ThreadAffinityHeterogeneousNodesTest, PinsAutomaticallyOnMasterOnly)
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, PinsAutomaticallyOnNonMasterOnly)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
-    setupNodes(&helper, { { 1, 2 } });
+    setupNodes(&helper, 1, 2);
     helper.expectWarningMatchingRegexIf("Oversubscribing the CPU", currentNode() == 0);
-    if (currentNode() == 1)
+    if (currentNode() >= 1)
     {
         helper.expectPinningMessage(false, 1);
+        expectNodeAffinitySet(&helper, currentNode(), indexInNode());
     }
-    expectNodeAffinitySet(&helper, 1, indexInNode());
     helper.setAffinity(1);
 }
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, HandlesInvalidOffsetOnNonMasterOnly)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::On);
     helper.setOffsetAndStride(2, 0);
-    setupNodes(&helper, { { 4, 2 } });
+    setupNodes(&helper, 4, 2);
     helper.expectWarningMatchingRegex("Applying core pinning offset 2");
-    helper.expectWarningMatchingRegexIf("Requested offset too large", isMaster() || currentNode() == 1);
+    helper.expectWarningMatchingRegexIf("Requested offset too large", isMaster() || currentNode() >= 1);
     if (currentNode() == 0)
     {
         helper.expectPinningMessage(false, 1);
@@ -226,11 +245,11 @@ TEST_F(ThreadAffinityHeterogeneousNodesTest, HandlesInvalidOffsetOnNonMasterOnly
 
 TEST_F(ThreadAffinityHeterogeneousNodesTest, HandlesInvalidStrideOnNonMasterOnly)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(RequireEvenRankCountWithAtLeastFourRanks);
     ThreadAffinityTestHelper helper;
     helper.setAffinityOption(ThreadAffinity::On);
     helper.setOffsetAndStride(0, 2);
-    setupNodes(&helper, { { 4, 2 } });
+    setupNodes(&helper, 4, 2);
     helper.expectWarningMatchingRegexIf("Requested stride too large", isMaster() || currentNode() == 1);
     if (currentNode() == 0)
     {
@@ -241,3 +260,5 @@ TEST_F(ThreadAffinityHeterogeneousNodesTest, HandlesInvalidStrideOnNonMasterOnly
 }
 
 } // namespace
+} // namespace test
+} // namespace gmx
index 142d4c19848f44455940b7bee976c16a93c01e33..d2a35adda283bab249f93a31be15eee09d1446ee 100644 (file)
@@ -63,7 +63,7 @@ gmx_add_unit_test(UtilityUnitTests utility-test
 # TODO: Explicitly link specific modules.
 target_link_libraries(utility-test PRIVATE legacy_modules)
 
-gmx_add_mpi_unit_test(UtilityMpiUnitTests utility-mpi-test 4
+gmx_add_mpi_unit_test(UtilityMpiUnitTests utility-mpi-test 2
     CPP_SOURCE_FILES
         physicalnodecommunicator_mpi.cpp
         )
index a97e344a88e618f22fc1de730ecc99844c69bfc5..215c243efe2bbaf4d711d425af9a543e3e0ac6f7 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018,2019, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,2021, 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.
 
 namespace gmx
 {
+namespace test
+{
 namespace
 {
 
 TEST(PhysicalNodeCommunicatorTest, CanConstruct)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     PhysicalNodeCommunicator comm(MPI_COMM_WORLD, 0);
 }
 
 TEST(PhysicalNodeCommunicatorTest, CanCallBarrier)
 {
-    GMX_MPI_TEST(4);
+    GMX_MPI_TEST(AllowAnyRankCount);
     PhysicalNodeCommunicator comm(MPI_COMM_WORLD, 0);
     comm.barrier();
 }
 
 } // namespace
+} // namespace test
 } // namespace gmx
index 26e5de90b52d421e1d730c1b906a3c69433b7d4d..bd84d04db35ae1ea4208b31d1b9c5be78bdd7c74 100644 (file)
@@ -60,6 +60,7 @@
 #include "gromacs/utility/stringutil.h"
 
 #include "testutils/cmdlinetest.h"
+#include "testutils/mpitest.h"
 
 #include "moduletest.h"
 #include "terminationhelper.h"
index a4e4aa37b532b3829d53afa1f89edc5929c0eda7..7dadd5d0ea9ca03032a09491bfe0d9801f661966 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2016,2017,2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2016,2017,2018,2019,2020,2021, 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.
@@ -59,6 +59,7 @@
 #include "gromacs/hardware/device_management.h"
 #include "gromacs/hardware/hw_info.h"
 #include "gromacs/trajectory/energyframe.h"
+#include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/gmxmpi.h"
 #include "gromacs/utility/physicalnodecommunicator.h"
index 37d543937cbfb64c74fcf46f14728e3dc59edda8..04b041b2b8220dc0b712bd6e1b39e95c86e83dd9 100644 (file)
 #include "config.h"
 
 #include <functional>
+#include <string>
 #include <type_traits>
 
-#include "gromacs/utility/basenetwork.h"
-
 namespace gmx
 {
 namespace test
@@ -64,6 +63,7 @@ namespace test
  * \ingroup module_testutils
  */
 int getNumberOfTestMpiRanks();
+
 //! \cond internal
 /*! \brief
  * Helper function for GMX_MPI_TEST().
@@ -73,58 +73,114 @@ int getNumberOfTestMpiRanks();
 bool threadMpiTestRunner(std::function<void()> testBody);
 //! \endcond
 
-/*! \brief
- * Declares that this test is an MPI-enabled unit test.
+/*! \brief Implementation of MPI test runner for thread-MPI
  *
- * \param[in] expectedRankCount Expected number of ranks for this test.
- *     The test will fail if run with unsupported number of ranks.
+ * See documentation GMX_MPI_TEST */
+#if GMX_THREAD_MPI
+#    define GMX_MPI_TEST_INNER                                                                  \
+        do                                                                                      \
+        {                                                                                       \
+            using MyTestClass = std::remove_reference_t<decltype(*this)>;                       \
+            if (!::gmx::test::threadMpiTestRunner([this]() { this->MyTestClass::TestBody(); })) \
+            {                                                                                   \
+                return;                                                                         \
+            }                                                                                   \
+        } while (0)
+#else
+#    define GMX_MPI_TEST_INNER
+#endif
+
+/*! \brief Declares that this test is an MPI-enabled unit test and
+ * expresses the conditions under which it can run.
  *
  * To write unit tests that run under MPI, you need to do a few things:
- *  - Put GMX_MPI_TEST() as the first statement in your test body and
- *    specify the number of ranks this test expects.
+ *  - Put GMX_MPI_TEST(RankRequirement) as the first statement in your
+ *    test body and either declare or use a suitable class as its
+ *    argument to express what requirements exist on the number of MPI
+ *    ranks for this test.
  *  - Declare your unit test in CMake with gmx_add_mpi_unit_test().
- *    Note that all tests in the binary should fulfill the conditions above,
- *    and work with the same number of ranks.
- * TODO: Figure out a mechanism for mixing tests with different rank counts in
- * the same binary (possibly, also MPI and non-MPI tests).
+ *    Note that all tests in the binary should fulfill the conditions above.
  *
  * When you do the above, the following will happen:
  *  - The test will get compiled only if thread-MPI or real MPI is enabled.
- *  - The test will get executed on the number of ranks specified.
- *    If you are using real MPI, the whole test binary is run under MPI and
- *    test execution across the processes is synchronized (GMX_MPI_TEST()
- *    actually has no effect in this case, the synchronization is handled at a
- *    higher level).
- *    If you are using thread-MPI, GMX_MPI_TEST() is required and it
- *    initializes thread-MPI with the specified number of threads and runs the
- *    rest of the test on each of the threads.
+ *  - The test will get executed only when the specified condition on
+ *    the the number of ranks is satisfied.
+ *  - If you are using real MPI, the whole test binary is run under
+ *    MPI and test execution across the processes is synchronized
+ *    (GMX_MPI_TEST() actually has no effect in this case, the
+ *    synchronization is handled at a higher level).
+ *  - If you are using thread-MPI, GMX_MPI_TEST() is required and it
+ *    initializes thread-MPI with the specified number of threads and
+ *    runs the rest of the test on each of the threads.
+ *
+ * \param[in] RankRequirement Class that expresses the necessary
+ *     conditions on the number of MPI ranks for the test to continue.
+ *     If run with unsupported number of ranks, the remainder of the
+ *     test body is skipped, and the GTEST_SKIP() mechanism used to
+ *     report the reason why the number of MPI ranks is unsuitable.
+ *
+ * The RankRequirement class must have two static members; a static
+ * method \c bool conditionSatisfied(const int) that can be passed the
+ * number of ranks present at run time and return whether the test can
+ * run with that number of ranks, and a static const string \c
+ * s_skipReason describing the reason why the test cannot be run, when
+ * that is the case.
  *
  * You need to be extra careful for variables in the test fixture, if you use
  * one: when run under thread-MPI, these will be shared across all the ranks,
  * while under real MPI, these are naturally different for each process.
  * Local variables in the test body are private to each rank in both cases.
  *
- * Currently, it is not possible to specify the number of ranks as one, because
- * that will lead to problems with (at least) thread-MPI, but such tests can be
- * written as serial tests anyways.
+ * Currently, it is not possible to require the use of a single MPI
+ * rank, because that will lead to problems with (at least)
+ * thread-MPI, but such tests can be written as serial tests anyway.
  *
  * \ingroup module_testutils
  */
-#if GMX_THREAD_MPI
-#    define GMX_MPI_TEST(expectedRankCount)                                                     \
-        do                                                                                      \
-        {                                                                                       \
-            ASSERT_EQ(expectedRankCount, ::gmx::test::getNumberOfTestMpiRanks());               \
-            using MyTestClass = std::remove_reference_t<decltype(*this)>;                       \
-            if (!::gmx::test::threadMpiTestRunner([this]() { this->MyTestClass::TestBody(); })) \
-            {                                                                                   \
-                return;                                                                         \
-            }                                                                                   \
-        } while (0)
-#else
-#    define GMX_MPI_TEST(expectedRankCount) \
-        ASSERT_EQ(expectedRankCount, ::gmx::test::getNumberOfTestMpiRanks())
-#endif
+#define GMX_MPI_TEST(RankRequirement)                                                         \
+    const int numRanks = ::gmx::test::getNumberOfTestMpiRanks();                              \
+    if (!RankRequirement::conditionSatisfied(numRanks))                                       \
+    {                                                                                         \
+        GTEST_SKIP() << std::string("Test skipped because ") + RankRequirement::s_skipReason; \
+        return;                                                                               \
+    }                                                                                         \
+    GMX_MPI_TEST_INNER;
+
+//! Helper for GMX_MPI_TEST to permit any rank count
+class AllowAnyRankCount
+{
+public:
+    /*! \brief Function called by GMX_MPI_CONDITIONAL_TEST to see
+     * whether the test conditions are satisifed */
+    static bool conditionSatisfied(const int /* numRanks */) { return true; }
+    //! Reason to echo when skipping the test
+    inline static const char* s_skipReason = "UNUSED - any rank count satisfies";
+};
+
+//! Helper for GMX_MPI_TEST to permit only a specific rank count
+template<int requiredNumRanks>
+class RequireRankCount
+{
+public:
+    //! Function to require a specific number of ranks
+    static bool conditionSatisfied(const int numRanks) { return numRanks == requiredNumRanks; }
+    //! Text to echo when skipping a test that does not satisfy the requirement
+    inline static const std::string s_skipReason =
+            std::to_string(requiredNumRanks) + " ranks are required";
+};
+
+//! Helper for GMX_MPI_TEST to permit only a specific rank count
+template<int minimumNumRanks>
+class RequireMinimumRankCount
+{
+public:
+    //! Function to require at least the minimum number of ranks
+    static bool conditionSatisfied(const int numRanks) { return numRanks >= minimumNumRanks; }
+    //! Text to echo when skipping a test that does not satisfy the requirement
+    inline static const std::string s_skipReason =
+            std::to_string(minimumNumRanks) + " or more ranks are required";
+};
+
 
 } // namespace test
 } // namespace gmx
index 834aa5cffa6df9b34ba280e3ed37bb5d048200e0..e04b49bec5da4fbc7e58cb724933ef3c1c998f02 100644 (file)
@@ -68,6 +68,7 @@
 #include "gromacs/utility/programcontext.h"
 #include "gromacs/utility/textwriter.h"
 
+#include "testutils/mpitest.h"
 #include "testutils/refdata.h"
 #include "testutils/test_hardware_environment.h"
 #include "testutils/testfilemanager.h"
@@ -88,7 +89,7 @@ namespace
  *
  * This context overrides the installationPrefix() implementation to always
  * load data files from the source directory, as the test binaries are never
- * installed.  It also support overriding the directory through a command-line
+ * installed.  It also supports overriding the directory through a command-line
  * option to the test binary.
  *
  * \ingroup module_testutils
@@ -142,6 +143,29 @@ void printHelp(const Options& options)
 // Never releases ownership.
 std::unique_ptr<TestProgramContext> g_testContext;
 
+/*! \brief Makes GoogleTest non-failures more verbose
+ *
+ * By default, GoogleTest does not echo messages appended to explicit
+ * assertions of success with SUCCEEDED() e.g.
+ *
+ *    GTEST_SKIP() << "reason why";
+ *
+ * produces no output. This test listener changes that behavior, so
+ * that the message is echoed.
+ *
+ * When run with multiple ranks, only the master rank should use this
+ * listener, else the output can be very noisy. */
+class SuccessListener : public testing::EmptyTestEventListener
+{
+    void OnTestPartResult(const testing::TestPartResult& result) override
+    {
+        if (result.type() == testing::TestPartResult::kSuccess)
+        {
+            printf("%s\n", result.message());
+        }
+    }
+};
+
 } // namespace
 
 //! \cond internal
@@ -198,6 +222,7 @@ void initTestUtils(const char* dataPath,
 
         bool        bHelp = false;
         std::string sourceRoot;
+        bool        echoReasons = false;
         Options     options;
         // TODO: A single option that accepts multiple names would be nicer.
         // Also, we recognize -help, but GTest doesn't, which leads to a bit
@@ -209,6 +234,8 @@ void initTestUtils(const char* dataPath,
         // TODO: Make this into a FileNameOption (or a DirectoryNameOption).
         options.addOption(
                 StringOption("src-root").store(&sourceRoot).description("Override source tree location (for data files)"));
+        options.addOption(
+                BooleanOption("echo-reasons").store(&echoReasons).description("When succeeding or skipping a test, echo the reason"));
         // The potential MPI test event listener must be initialized first,
         // because it should appear in the start of the event listener list,
         // before other event listeners that may generate test failures
@@ -240,6 +267,11 @@ void initTestUtils(const char* dataPath,
             g_testContext->overrideSourceRoot(sourceRoot);
             TestFileManager::setInputDataDirectory(Path::join(sourceRoot, dataPath));
         }
+        // Echo success messages only from the master MPI rank
+        if (echoReasons && (gmx_node_rank() == 0))
+        {
+            testing::UnitTest::GetInstance()->listeners().Append(new SuccessListener);
+        }
     }
     catch (const std::exception& ex)
     {
index d41ee8d9281162648e952c26a2a47cea4f06d757..f3d8c0f4eca09bef4bea9c4eaa7090aec062614b 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2016,2019, by the GROMACS development team, led by
+ * Copyright (c) 2016,2019,2021, 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.
 
 #include "config.h"
 
+#include <vector>
+
 #include <gtest/gtest.h>
+#include <gmock/gmock.h>
 
 #include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/gmxmpi.h"
 
+namespace gmx
+{
+namespace test
+{
 namespace
 {
 
 class MpiSelfTest : public ::testing::Test
 {
 public:
-    MpiSelfTest() : reached{ 0, 0 } {}
-
-    int reached[2];
+    //! Whether each rank participated, relevant only on rank 0
+    std::vector<int> reached_;
 };
 
 TEST_F(MpiSelfTest, Runs)
 {
-    GMX_MPI_TEST(2);
-#if GMX_THREAD_MPI
-    reached[gmx_node_rank()] = 1;
-    MPI_Barrier(MPI_COMM_WORLD);
-#else
+    GMX_MPI_TEST(RequireMinimumRankCount<2>);
+    if (gmx_node_rank() == 0)
+    {
+        reached_.resize(getNumberOfTestMpiRanks(), 0);
+    }
     int value = 1;
-    MPI_Gather(&value, 1, MPI_INT, reached, 1, MPI_INT, 0, MPI_COMM_WORLD);
-#endif
+    MPI_Gather(&value, 1, MPI_INT, reached_.data(), 1, MPI_INT, 0, MPI_COMM_WORLD);
     if (gmx_node_rank() == 0)
     {
-        EXPECT_EQ(1, reached[0]);
-        EXPECT_EQ(1, reached[1]);
+        EXPECT_THAT(reached_, testing::Each(value));
     }
 }
 
 } // namespace
+} // namespace test
+} // namespace gmx