Introduce function that checks if the computations can be done on a GPU
authorArtem Zhmurov <zhmurov@gmail.com>
Tue, 23 Jul 2019 13:54:53 +0000 (15:54 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 31 Jul 2019 12:02:00 +0000 (14:02 +0200)
The following conditions should be fulfilled in order to compute on a GPU:
1. The code should be compiled with GPU support.
2. There should be compatible GPUs installed in the system.
3. These GPUs shuold be detectible from the code.
4. Environmental variable GMX_DISABLE_GPU_DETECTION should not be set.

TODO: This function is used in tests, so it should be eventually moved
      to src/testutils.

Change-Id: I40ce381b68fcebcfab330db4b48a1bec2a2b7a57

src/gromacs/gpu_utils/CMakeLists.txt
src/gromacs/gpu_utils/gpu_testutils.cpp [new file with mode: 0644]
src/gromacs/gpu_utils/gpu_testutils.h [new file with mode: 0644]
src/gromacs/mdlib/tests/constr.cpp
src/gromacs/mdlib/tests/constr_impl.h
src/gromacs/mdlib/tests/leapfrog.cpp
src/gromacs/mdlib/tests/settle.cpp
src/programs/mdrun/tests/pmetest.cpp

index 20375436798e9137ae5a63855096ccf7b34c59cd..f9aa1da03fef062c9ed55c66d92f732d8190dd8d 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2015,2016,2017,2018, by the GROMACS development team, led by
+# Copyright (c) 2015,2016,2017,2018,2019, 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.
@@ -40,6 +40,7 @@ gmx_add_libgromacs_sources(
         clfftinitializer.cpp
         hostallocator.cpp
         gpu_utils.cpp
+        gpu_testutils.cpp
         )
 if(GMX_USE_OPENCL)
     gmx_add_libgromacs_sources(
diff --git a/src/gromacs/gpu_utils/gpu_testutils.cpp b/src/gromacs/gpu_utils/gpu_testutils.cpp
new file mode 100644 (file)
index 0000000..50be2b3
--- /dev/null
@@ -0,0 +1,58 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2019, 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 Function definitions for GPU detection, specific for tests.
+ *
+ *  \author Artem Zhmurov <zhmurov@gmail.com>
+ */
+#include "gmxpre.h"
+
+#include "gpu_testutils.h"
+
+#include "gromacs/gpu_utils/gpu_utils.h"
+#include "gromacs/hardware/gpu_hw_info.h"
+
+bool canComputeOnGpu()
+{
+    bool           canComputeOnGpu = false;
+    gmx_gpu_info_t gpuInfo {};
+    if (canPerformGpuDetection())
+    {
+        findGpus(&gpuInfo);
+        canComputeOnGpu = !getCompatibleGpus(gpuInfo).empty();
+    }
+    free_gpu_info(&gpuInfo);
+    return canComputeOnGpu;
+}
diff --git a/src/gromacs/gpu_utils/gpu_testutils.h b/src/gromacs/gpu_utils/gpu_testutils.h
new file mode 100644 (file)
index 0000000..1ea8227
--- /dev/null
@@ -0,0 +1,61 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2019, 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.
+ */
+/*! \libinternal \file
+ *  \brief Declare functions for detection of GPU devices, specific for tests.
+ *
+ *  \todo This should eventually go to src/testutils
+ *
+ *  \author Artem Zhmurov <zhmurov@gmail.com>
+ *
+ *  \inlibraryapi
+ */
+
+#ifndef GMX_GPU_UTILS_GPU_TESTUTILS_H
+#define GMX_GPU_UTILS_GPU_TESTUTILS_H
+
+/*! \brief Checks if there is a compatible GPU to run the computations on
+ *
+ * There are several reasons why code can not rune on the GPU:
+ * 1. The GPU can not be detected, because there is none in the system.
+ * 2. GPU detection is disabled by GMX_DISABLE_GPU_DETECTION environmental variable.
+ * 3. GPUs are detected, but none of them is compatible.
+ * This function checks all these conditions and returns true only if there at least
+ * one GPU that can be used for computations.
+ *
+ * \returns True, if there a GPU that can be used for computations
+ */
+bool canComputeOnGpu();
+
+#endif // GMX_GPU_UTILS_GPU_TESTUTILS_H
index 2b14d9e73d410a30aff6fe7f3a8e73e8a8cf3d02..4d74583aec48c3c9b5369b45d9e960c56fe86504 100644 (file)
@@ -312,8 +312,7 @@ std::vector<std::string> getAlgorithmsNames()
 {
     algorithmsNames.emplace_back("SHAKE");
     algorithmsNames.emplace_back("LINCS");
-    // TODO: Here we should check that at least 1 suitable GPU is available
-    if (GMX_GPU == GMX_GPU_CUDA && canPerformGpuDetection())
+    if (GMX_GPU == GMX_GPU_CUDA && canComputeOnGpu())
     {
         algorithmsNames.emplace_back("LINCS_CUDA");
     }
index c46d00e2f4f7a6ff4abb69f1d10193489ee6ad2f..875591496895fb12c0978eb4e65b6e88dfdd1a5a 100644 (file)
@@ -59,7 +59,7 @@
 #include "gromacs/fileio/gmxfio.h"
 #include "gromacs/gmxlib/nrnb.h"
 #include "gromacs/gmxlib/nonbonded/nonbonded.h"
-#include "gromacs/gpu_utils/gpu_utils.h"
+#include "gromacs/gpu_utils/gpu_testutils.h"
 #include "gromacs/math/paddedvector.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/math/vectypes.h"
index 319bb41efd9f5743f751a062526e1b0f01e05ef3..1e3215d5b934d944f4f9898206fba93600f3a801 100644 (file)
@@ -61,7 +61,7 @@
 
 #include <gtest/gtest.h>
 
-#include "gromacs/gpu_utils/gpu_utils.h"
+#include "gromacs/gpu_utils/gpu_testutils.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/mdlib/leapfrog_cuda.h"
 #include "gromacs/utility/stringutil.h"
@@ -166,15 +166,29 @@ class IntegratorTest : public ::testing::TestWithParam<IntegratorTestParameters>
             mdAtoms_.invmass = inverseMasses_.data();
         }
 
+        //! Store whether any compatible GPUs exist.
+        static bool s_hasCompatibleGpus;
+        //! Before any test is run, work out whether any compatible GPUs exist.
+        static void SetUpTestCase()
+        {
+            s_hasCompatibleGpus = canComputeOnGpu();
+        }
 };
 
+bool IntegratorTest::s_hasCompatibleGpus = false;
+
+// The test will run only if:
+// 1. The code was compiled with CUDA
+// 2. There is a CUDA-capable GPU in a system
+// 3. This GPU is detectable
+// 4. GPU detection was not disabled by GMX_DISABLE_GPU_DETECTION environment variable
 TEST_P(IntegratorTest, SimpleIntegration)
 {
-    // TODO: Here we should check that at least 1 suitable GPU is available
-    if (!canPerformGpuDetection())
+    if (!s_hasCompatibleGpus)
     {
         return;
     }
+
     int  numAtoms;    // 1. Number of atoms
     real timestep;    // 2. Timestep
     rvec v0;          // 3. Velocity
index 642ceb4117b68f4ca4e050270ae820603a2665ac..21e3a641b0e1dca20145957466a6e56b26808985 100644 (file)
@@ -43,7 +43,7 @@
 
 #include <gtest/gtest.h>
 
-#include "gromacs/gpu_utils/gpu_utils.h"
+#include "gromacs/gpu_utils/gpu_testutils.h"
 #include "gromacs/math/paddedvector.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/math/vectypes.h"
@@ -133,8 +133,18 @@ class SettleTest : public ::testing::TestWithParam<SettleTestParameters>
                 velocities_[i / DIM][i % DIM] = 0.0;
             }
         }
+
+        //! Store whether any compatible GPUs exist.
+        static bool s_hasCompatibleGpus;
+        //! Before any test is run, work out whether any compatible GPUs exist.
+        static void SetUpTestCase()
+        {
+            s_hasCompatibleGpus = canComputeOnGpu();
+        }
 };
 
+bool SettleTest::s_hasCompatibleGpus = false;
+
 TEST_P(SettleTest, SatisfiesConstraints)
 {
     int  numSettles;
@@ -284,12 +294,13 @@ TEST_P(SettleTest, SatisfiesConstraints)
         }
     }
 
-    // CUDA version will be tested only if
-    // 1. The code was compiled with cuda
+    // CUDA version will be tested only if:
+    // 1. The code was compiled with CUDA
     // 2. There is a CUDA-capable GPU in a system
+    // 3. This GPU is detectable
+    // 4. GPU detection was not disabled by GMX_DISABLE_GPU_DETECTION environment variable
 #if GMX_GPU == GMX_GPU_CUDA
-    // TODO: Here we should check that at least 1 suitable GPU is available
-    if (canPerformGpuDetection())
+    if (s_hasCompatibleGpus)
     {
         // Run the CUDA code and check if it gives identical results to CPU code
         t_idef idef;
index 1ab8fa444bee11d1430b8cec76ca61b56de3ddd0..89d4291f6d2fec8159a2e6ffdbe8c7116021dc95 100644 (file)
@@ -55,7 +55,7 @@
 #include <gtest/gtest-spi.h>
 
 #include "gromacs/ewald/pme.h"
-#include "gromacs/gpu_utils/gpu_utils.h"
+#include "gromacs/gpu_utils/gpu_testutils.h"
 #include "gromacs/hardware/detecthardware.h"
 #include "gromacs/hardware/gpu_hw_info.h"
 #include "gromacs/trajectory/energyframe.h"
@@ -98,18 +98,7 @@ bool PmeTest::s_hasCompatibleGpus = false;
 
 void PmeTest::SetUpTestCase()
 {
-    gmx_gpu_info_t gpuInfo {};
-    // It would be nicer to do this detection once and have mdrun
-    // re-use it, but this is OK. Note that this also caters for when
-    // there is no GPU support in the build.
-    //
-    // TODO report any error messages gracefully.
-    if (canPerformGpuDetection())
-    {
-        findGpus(&gpuInfo);
-        s_hasCompatibleGpus = (gpuInfo.n_dev_compatible > 0);
-    }
-    free_gpu_info(&gpuInfo);
+    s_hasCompatibleGpus = canComputeOnGpu();
 }
 
 void PmeTest::runTest(const RunModesList &runModes)