Add integration test for empty domain
authorMark Abraham <mark.j.abraham@gmail.com>
Fri, 15 May 2015 14:40:17 +0000 (16:40 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 21 May 2015 09:55:50 +0000 (11:55 +0200)
This box is large enough that the default 2D DD will contain an
empty domain (at least initially).

Current Jenkins will run this with two domains when run with real MPI.
Extended the CMakery to make this also work with thread-MPI.

Refs #1734

Change-Id: I7b7269cdd1faaa562afeb7a5dd3f75fc19ceff85

src/programs/mdrun/tests/CMakeLists.txt
src/programs/mdrun/tests/domain_decomposition.cpp [new file with mode: 0644]
src/programs/mdrun/tests/moduletest.cpp
src/programs/mdrun/tests/moduletest.h
src/programs/mdrun/tests/spc2_and_vacuum.gro [new file with mode: 0644]
src/testutils/TestMacros.cmake

index 8467cf6f5fff3323c6bd3ebd889f02a3c280f6c0..c5eb3730d479c25160aa7f8837100ff43dc61ddb 100644 (file)
@@ -64,6 +64,7 @@ gmx_build_unit_test(
     multisim.cpp
     multisimtest.cpp
     replicaexchange.cpp
+    domain_decomposition.cpp
     # files with code for test fixtures
     moduletest.cpp
     # pseudo-library for code for mdrun
diff --git a/src/programs/mdrun/tests/domain_decomposition.cpp b/src/programs/mdrun/tests/domain_decomposition.cpp
new file mode 100644 (file)
index 0000000..5dd3fa3
--- /dev/null
@@ -0,0 +1,77 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2015, 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.
+ */
+
+/*! \internal \file
+ * \brief
+ * Tests special cases in domain decomposition
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ * \ingroup module_mdrun_integration_tests
+ */
+#include "gmxpre.h"
+
+#include <string>
+
+#include <gtest/gtest.h>
+
+#include "testutils/cmdlinetest.h"
+
+#include "moduletest.h"
+
+namespace
+{
+
+//! Test fixture for domain decomposition special cases
+class DomainDecompositionSpecialCasesTest : public gmx::test::MdrunTestFixture
+{
+};
+
+//! When run with 2+ domains, ensures an empty cell, to make sure that zero-sized things work
+TEST_F(DomainDecompositionSpecialCasesTest, AnEmptyDomainWorks)
+{
+    runner_.useStringAsMdpFile("cutoff-scheme = Verlet\n");
+    runner_.useTopGroAndNdxFromDatabase("spc2");
+    /* Over-ride the use of spc2.gro with one with both waters close
+     * together.
+     *
+     * TODO Maybe refactor to avoid over-writing the .gro filename set
+     * above. */
+    runner_.useGroFromDatabase("spc2_and_vacuum");
+    ASSERT_EQ(0, runner_.callGrompp());
+
+    ASSERT_EQ(0, runner_.callMdrun());
+}
+
+} // namespace
index 718e4c647bc30eca1b20d60aa1c4a1d1c2021759..845aa58ac0fadebb683eebd647bb366c4567b342 100644 (file)
@@ -51,6 +51,7 @@
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/file.h"
+#include "gromacs/utility/gmxmpi.h"
 #include "programs/mdrun/mdrun_main.h"
 
 #include "testutils/cmdlinetest.h"
@@ -119,6 +120,7 @@ SimulationRunner::SimulationRunner(IntegrationTestFixture *fixture) :
 void
 SimulationRunner::useEmptyMdpFile()
 {
+    // TODO When removing the group scheme, update actual and potential users of useEmptyMdpFile
     useStringAsMdpFile("cutoff-scheme = Group\n");
 }
 
@@ -148,6 +150,12 @@ SimulationRunner::useTopGroAndNdxFromDatabase(const char *name)
     ndxFileName_ = fixture_->fileManager_.getInputFilePath((std::string(name) + ".ndx").c_str());
 }
 
+void
+SimulationRunner::useGroFromDatabase(const char *name)
+{
+    groFileName_ = fixture_->fileManager_.getInputFilePath((std::string(name) + ".gro").c_str());
+}
+
 int
 SimulationRunner::callGromppOnThisRank()
 {
@@ -167,15 +175,23 @@ SimulationRunner::callGromppOnThisRank()
 int
 SimulationRunner::callGrompp()
 {
+    int returnValue = 0;
 #ifdef GMX_LIB_MPI
-    // When compiled with external MPI, only call one instance of the
-    // grompp function
-    if (0 != gmx_node_rank())
+    // When compiled with external MPI, we're trying to run mdrun with
+    // MPI, but we need to make sure that we only do grompp on one
+    // rank
+    if (0 == gmx_node_rank())
+#endif
     {
-        return 0;
+        returnValue = callGromppOnThisRank();
     }
+#ifdef GMX_LIB_MPI
+    // Make sure rank zero has written the .tpr file before other
+    // ranks try to read it. Thread-MPI and serial do this just fine
+    // on their own.
+    MPI_Barrier(MPI_COMM_WORLD);
 #endif
-    return callGromppOnThisRank();
+    return returnValue;
 }
 
 int
@@ -185,6 +201,7 @@ SimulationRunner::callMdrun(const CommandLine &callerRef)
        to this function. Passing a non-const reference might make it
        easier to write code that incorrectly re-uses callerRef after
        the call to this function. */
+
     CommandLine caller(callerRef);
     caller.addOption("-s", tprFileName_);
 
@@ -200,9 +217,22 @@ SimulationRunner::callMdrun(const CommandLine &callerRef)
         caller.addOption("-nsteps", nsteps_);
     }
 
+#ifdef GMX_MPI
+#  ifdef GMX_GPU
+#    ifdef GMX_THREAD_MPI
+    int         numGpusNeeded = g_numThreads;
+#    else   /* Must be real MPI */
+    int         numGpusNeeded = gmx_node_num();
+#    endif
+    std::string gpuIdString(numGpusNeeded, '0');
+    caller.addOption("-gpu_id", gpuIdString.c_str());
+#  endif
+#endif
+
 #ifdef GMX_THREAD_MPI
     caller.addOption("-nt", g_numThreads);
 #endif
+
 #ifdef GMX_OPENMP
     caller.addOption("-ntomp", g_numOpenMPThreads);
 #endif
index d021dc862180227bd8e4bb69ace3a993a0b70e14..6a5a56eae3fb4ceb2bbbfd4d9ef8915fdf6ad63b 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015, 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.
@@ -87,6 +87,8 @@ class SimulationRunner
         void useStringAsNdxFile(const char *ndxString);
         //! Use a standard .top and .gro file as input to grompp
         void useTopGroAndNdxFromDatabase(const char *name);
+        //! Use a standard .gro file as input to grompp
+        void useGroFromDatabase(const char *name);
         //! Calls grompp (on rank 0) to prepare for the mdrun test
         int callGrompp();
         //! Calls grompp (on this rank) to prepare for the mdrun test
diff --git a/src/programs/mdrun/tests/spc2_and_vacuum.gro b/src/programs/mdrun/tests/spc2_and_vacuum.gro
new file mode 100644 (file)
index 0000000..b2f66ca
--- /dev/null
@@ -0,0 +1,9 @@
+Two nearby SPC waters, and a bunch of empty vacuum
+ 6
+    1SOL     OW    1    .230    .628    .113
+    1SOL    HW1    2    .137    .626    .150
+    1SOL    HW2    3    .231    .589    .021
+    2SOL     OW    4    .225    .275   -.866
+    2SOL    HW1    5    .260    .258   -.774
+    2SOL    HW2    6    .137    .230   -.878
+   3.01000   3.01000   3.01000
index b6ef24f748f6797d26c9b15eaba8c9f85e17cffb..bdfa55782c13183e69b9357335a8e60d9f134067 100644 (file)
@@ -92,39 +92,53 @@ endfunction ()
 # that ctest can run the test binary over a range of numbers of MPI
 # ranks.
 function (gmx_register_mpi_integration_test NAME EXENAME NUMPROC)
-    if (GMX_BUILD_UNITTESTS AND BUILD_TESTING AND GMX_MPI)
-        foreach(VARNAME MPIEXEC MPIEXEC_NUMPROC_FLAG MPIEXEC_PREFLAGS MPIEXEC_POSTFLAGS)
-            # These variables need a valid value for the test to run
-            # and pass, but conceivably any of them might be valid
-            # with arbitrary (including empty) content. They can't be
-            # valid if they've been populated with the CMake
-            # find_package magic suffix/value "NOTFOUND", though.
-            if (${VARNAME} MATCHES ".*NOTFOUND")
-                message(STATUS "CMake variable ${VARNAME} was not detected to be a valid value. To test GROMACS correctly, check the advice in the install guide.")
+    if (GMX_BUILD_UNITTESTS AND BUILD_TESTING)
+        if (GMX_MPI)
+            foreach(VARNAME MPIEXEC MPIEXEC_NUMPROC_FLAG MPIEXEC_PREFLAGS MPIEXEC_POSTFLAGS)
+                # These variables need a valid value for the test to run
+                # and pass, but conceivably any of them might be valid
+                # with arbitrary (including empty) content. They can't be
+                # valid if they've been populated with the CMake
+                # find_package magic suffix/value "NOTFOUND", though.
+                if (${VARNAME} MATCHES ".*NOTFOUND")
+                    message(STATUS "CMake variable ${VARNAME} was not detected to be a valid value. To test GROMACS correctly, check the advice in the install guide.")
+                    set(_cannot_run_mpi_tests 1)
+                endif()
+                if (NOT VARNAME STREQUAL MPIEXEC AND ${VARNAME})
+                    set(_an_mpi_variable_had_content 1)
+                endif()
+            endforeach()
+            if(_an_mpi_variable_had_content AND NOT MPIEXEC)
+                message(STATUS "CMake variable MPIEXEC must have a valid value if one of the other related MPIEXEC variables does. To test GROMACS correctly, check the advice in the install guide.")
                 set(_cannot_run_mpi_tests 1)
             endif()
-            if (NOT VARNAME STREQUAL MPIEXEC AND ${VARNAME})
-                set(_an_mpi_variable_had_content 1)
+            if(NOT _cannot_run_mpi_tests)
+                add_test(NAME ${NAME}
+                    COMMAND
+                    ${MPIEXEC} ${MPIEXEC_NUMPROC_FLAG} ${NUMPROC}
+                    ${MPIEXEC_PREFLAGS} $<TARGET_FILE:${EXENAME}> ${MPIEXEC_POSTFLAGS}
+                    --gtest_output=xml:${CMAKE_BINARY_DIR}/Testing/Temporary/${NAME}.xml
+                    )
+                set_tests_properties(${testname} PROPERTIES LABELS "MpiIntegrationTest")
+                add_dependencies(tests ${EXENAME})
             endif()
-        endforeach()
-        if(_an_mpi_variable_had_content AND NOT MPIEXEC)
-            message(STATUS "CMake variable MPIEXEC must have a valid value if one of the other related MPIEXEC variables does. To test GROMACS correctly, check the advice in the install guide.")
-            set(_cannot_run_mpi_tests 1)
-        endif()
-        if(NOT _cannot_run_mpi_tests)
+
+            # GMX_EXTRA_LIBRARIES might be needed for mdrun integration tests at
+            # some point.
+            # target_link_libraries(${EXENAME} ${GMX_EXTRA_LIBRARIES})
+        elseif(GMX_THREAD_MPI)
             add_test(NAME ${NAME}
                 COMMAND
-                ${MPIEXEC} ${MPIEXEC_NUMPROC_FLAG} ${NUMPROC}
-                ${MPIEXEC_PREFLAGS} $<TARGET_FILE:${EXENAME}> ${MPIEXEC_POSTFLAGS}
+                $<TARGET_FILE:${EXENAME}> -nt ${NUMPROC}
                 --gtest_output=xml:${CMAKE_BINARY_DIR}/Testing/Temporary/${NAME}.xml
                 )
             set_tests_properties(${testname} PROPERTIES LABELS "MpiIntegrationTest")
             add_dependencies(tests ${EXENAME})
-        endif()
 
-        # GMX_EXTRA_LIBRARIES might be needed for mdrun integration tests at
-        # some point.
-        # target_link_libraries(${EXENAME} ${GMX_EXTRA_LIBRARIES})
+            # GMX_EXTRA_LIBRARIES might be needed for mdrun integration tests at
+            # some point.
+            # target_link_libraries(${EXENAME} ${GMX_EXTRA_LIBRARIES})
+        endif()
     endif()
 endfunction ()