Move initialization of clFFT
authorMark Abraham <mark.j.abraham@gmail.com>
Mon, 2 Sep 2019 05:45:39 +0000 (07:45 +0200)
committerSzilárd Páll <pall.szilard@gmail.com>
Mon, 2 Sep 2019 21:44:57 +0000 (23:44 +0200)
Gave ClfftInitializer the responsibility for mutual exclusion, which
means the initialization is now convenient to do alongside other
PME-on-GPU initialization tasks. This simplifies the code.

Removed mention of lazy initialization, which was not implemented at
all.

Refs #2535

Change-Id: I429767b059ddc3b4f16f2f4a8830b54ed1f49fab

src/gromacs/ewald/pme_gpu_internal.cpp
src/gromacs/ewald/pme_gpu_types_host.h
src/gromacs/gpu_utils/clfftinitializer.cpp
src/gromacs/gpu_utils/clfftinitializer.h
src/gromacs/gpu_utils/tests/CMakeLists.txt
src/gromacs/gpu_utils/tests/clfftinitializer.cpp [new file with mode: 0644]
src/gromacs/mdrun/runner.cpp

index 7ca9cabf6d2ede478b022e4f6286c03cbc71665d..6e328ff24a14ea1f1cb9ac5b1eb4792bb291fe20 100644 (file)
@@ -823,6 +823,8 @@ static void pme_gpu_init(gmx_pme_t               *pme,
     GMX_ASSERT(pmeGpuProgram != nullptr, "GPU kernels must be already compiled");
     pmeGpu->programHandle_ = pmeGpuProgram;
 
+    pmeGpu->initializedClfftLibrary_ = std::make_unique<gmx::ClfftInitializer>();
+
     pme_gpu_init_internal(pmeGpu);
     pme_gpu_alloc_energy_virial(pmeGpu);
 
index eacc911dfb5d504b9facc142cf78f301c62caae7..f14da288d4db42093d60189ce664eee21c86b58b 100644 (file)
@@ -55,6 +55,7 @@
 
 #include "gromacs/ewald/pme.h"
 #include "gromacs/ewald/pme_gpu_program.h"
+#include "gromacs/gpu_utils/clfftinitializer.h"
 #include "gromacs/gpu_utils/gpu_utils.h"      // for GpuApiCallBehavior
 #include "gromacs/gpu_utils/hostallocator.h"
 #include "gromacs/math/vectypes.h"
@@ -189,6 +190,9 @@ struct PmeGpu
     //! A handle to the program created by buildPmeGpuProgram()
     PmeGpuProgramHandle programHandle_;
 
+    //! Handle that ensures the clFFT library has been initialized once per process.
+    std::unique_ptr<gmx::ClfftInitializer> initializedClfftLibrary_;
+
     /*! \brief The settings. */
     PmeGpuSettings settings;
 
index 693c32983a2ed019c9c16998deaa6e1b4c5fcfab..458afc963481c932de75654ba61d978e5458d363 100644 (file)
@@ -45,9 +45,9 @@
 
 #include "config.h"
 
-#include <memory>
-
+#include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/mutex.h"
 #include "gromacs/utility/stringutil.h"
 
 #if GMX_GPU == GMX_GPU_OPENCL
 namespace gmx
 {
 
+namespace
+{
+
+/*! \brief The clFFT library may only be initialized once per process,
+ * and this is orchestrated by this shared value and mutex.
+ *
+ * This ensures that thread-MPI and OpenMP builds can't accidentally
+ * initialize it more than once. */
+//! @{
+gmx_unused bool g_clfftInitialized = false;
+gmx::Mutex g_clfftMutex;
+//! @}
+
+}   // namespace
+
 ClfftInitializer::ClfftInitializer()
 {
 #if GMX_GPU == GMX_GPU_OPENCL
-    clfftSetupData fftSetup;
-    int            initErrorCode = clfftInitSetupData(&fftSetup);
+    gmx::lock_guard<gmx::Mutex> guard(g_clfftMutex);
+    clfftSetupData              fftSetup;
+    int                         initErrorCode = clfftInitSetupData(&fftSetup);
     if (initErrorCode != 0)
     {
         GMX_THROW(InternalError(formatString("Failed to initialize the clFFT library, error code %d", initErrorCode)));
@@ -71,20 +87,23 @@ ClfftInitializer::ClfftInitializer()
     {
         GMX_THROW(InternalError(formatString("Failed to initialize the clFFT library, error code %d", initErrorCode)));
     }
+    g_clfftInitialized = true;
+#else
+    GMX_UNUSED_VALUE(g_clfftInitialized);
 #endif
 }
 
 ClfftInitializer::~ClfftInitializer()
 {
 #if GMX_GPU == GMX_GPU_OPENCL
-    clfftTeardown();
-    // TODO: log non-0 return values (errors)
+    gmx::lock_guard<gmx::Mutex> guard(g_clfftMutex);
+    if (g_clfftInitialized)
+    {
+        clfftTeardown();
+        // TODO: log non-zero return values (errors)
+    }
+    g_clfftInitialized = false;
 #endif
 }
 
-std::unique_ptr<ClfftInitializer> initializeClfftLibrary()
-{
-    return std::make_unique<ClfftInitializer>();
-}
-
 }  // namespace gmx
index 18c63eb5823c5bac2d4eea669f8a93047011d699..4e7d07d4065fc88fef7acc4df85998bc4253a28c 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 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.
@@ -59,27 +59,16 @@ namespace gmx
 {
 
 /*! \libinternal
- * \brief Handle clFFT library init and tear down in RAII style. */
-class ClfftInitializer
-{
-    public:
-        ClfftInitializer();
-        ~ClfftInitializer();
-
-        GMX_DISALLOW_COPY_AND_ASSIGN(ClfftInitializer);
-};
-
-/*! \brief This routine should be called during setup for running
- * an FFT task on an OpenCL device.
- *
- * It should be called once per process with such a task.
+ * \brief Handle clFFT library init and tear down in RAII style also
+ * with mutual exclusion.
  *
- * It implements lazy initialization, so that we can have a lifetime
- * that begins when we know that PME on an OpenCL device will run an
- * FFT task on the device, and should continue until we know that the
- * last such task has completed. Any time required for this
- * initialization or tear down should not be accrued to per-MD-step
- * counters.
+ * Only one thread per process needs to attempt to set up and tear
+ * down the clFFT library, but this wrapper object ensures that it is
+ * safe to do so more than once, or from any thread. It is the
+ * responsibility of the caller to ensure that no use of the clFFT
+ * library is made before making an object of this type, or after the
+ * first such object is destroyed. If no more objects remain to be
+ * destroyed, then it is safe to create another and resume clFFT work.
  *
  * \todo Consider making a composite object that also handles
  * on-demand compilation, managing lifetime of PME FFT kernel programs
@@ -90,9 +79,23 @@ class ClfftInitializer
  * counter idiom so that both static initialization and
  * deinitialization can work in a fast, leak-free, and thread-safe way
  * without imposing constraints on the calling code.
- * See Redmine #2535.
- */
-std::unique_ptr<ClfftInitializer> initializeClfftLibrary();
+ * See Redmine #2535. */
+class ClfftInitializer
+{
+    public:
+        /*! \brief Constructor
+         *
+         * This initializers the clFFT library if there is a
+         * possibility of an FFT task on the device, and preserves it
+         * until destruction. Any time required for this
+         * initialization or tear down should not be accrued to
+         * per-MD-step counters. */
+        ClfftInitializer();
+        //! Destructor
+        ~ClfftInitializer();
+
+        GMX_DISALLOW_COPY_AND_ASSIGN(ClfftInitializer);
+};
 
 }  // namespace gmx
 
index 5a61beea2ad4a76e9795c154602be849b98cc6b2..c2fc3426bb09a97d7727b1b60ae05fde192c038d 100644 (file)
@@ -40,6 +40,7 @@
 
 # Always compiled as plain C++
 file(GLOB SOURCES_FROM_CXX
+    clfftinitializer.cpp
     hostallocator.cpp
     )
 
diff --git a/src/gromacs/gpu_utils/tests/clfftinitializer.cpp b/src/gromacs/gpu_utils/tests/clfftinitializer.cpp
new file mode 100644 (file)
index 0000000..d16a4d5
--- /dev/null
@@ -0,0 +1,66 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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.
+ *
+ * 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 for clFFT initialization logic
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ */
+#include "gmxpre.h"
+
+#include "gromacs/gpu_utils/clfftinitializer.h"
+
+#include <gtest/gtest.h>
+
+namespace gmx
+{
+
+namespace test
+{
+
+TEST(ClfftInitializer, SingleInitializationWorks)
+{
+    ClfftInitializer clfft;
+}
+
+TEST(ClfftInitializer, TwoInitializationsWork)
+{
+    // If the mutual exclusion was wrong, this would trigger any check
+    // that clFFT might have for double initialization.
+    ClfftInitializer first, second;
+}
+
+} // namespace test
+} // namespace gmx
index a35c67c4f7ba96948578ee3184a50babda676197..6333756b1a0cfaa4e64372f1ebbbee31a1ab63d5 100644 (file)
@@ -70,7 +70,6 @@
 #include "gromacs/fileio/tpxio.h"
 #include "gromacs/gmxlib/network.h"
 #include "gromacs/gmxlib/nrnb.h"
-#include "gromacs/gpu_utils/clfftinitializer.h"
 #include "gromacs/gpu_utils/gpu_utils.h"
 #include "gromacs/hardware/cpuinfo.h"
 #include "gromacs/hardware/detecthardware.h"
@@ -1198,8 +1197,6 @@ int Mdrunner::mdrunner()
         }
     }
 
-    std::unique_ptr<ClfftInitializer> initializedClfftLibrary;
-
     gmx_device_info_t                *pmeDeviceInfo = nullptr;
     // Later, this program could contain kernels that might be later
     // re-used as auto-tuning progresses, or subsequent simulations
@@ -1213,13 +1210,6 @@ int Mdrunner::mdrunner()
         pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_);
         init_gpu(pmeDeviceInfo);
         pmeGpuProgram = buildPmeGpuProgram(pmeDeviceInfo);
-        // TODO It would be nice to move this logic into the factory
-        // function. See Redmine #2535.
-        bool isMasterThread = !GMX_THREAD_MPI || MASTER(cr);
-        if (pmeRunMode == PmeRunMode::GPU && !initializedClfftLibrary && isMasterThread)
-        {
-            initializedClfftLibrary = initializeClfftLibrary();
-        }
     }
 
     /* getting number of PP/PME threads on this MPI / tMPI rank.