From a41344a06a7e6706b209571da2777574011e3980 Mon Sep 17 00:00:00 2001 From: Aleksei Iupinov Date: Mon, 7 May 2018 17:33:25 +0200 Subject: [PATCH] Added the bundled clFFT into OpenCL builds Used an object library, since we have no need of a real library, to have or to install, whether shared or static. Checked for the availability of dynamic loading, and made it available portably to libgromacs. Clfft initialization class is added and used in mdrunner to initialize/tear down clFFT library resources in a thread-safe manner, and only on ranks that require such setup. Noted TODOs for future work. Noted a useful style for explicit listing of source files. Refs #2500 Refs #2515 Refs #2535 Change-Id: I62d7d66f65e147bde17929ccc30abad36e2373c6 --- src/external/clFFT/src/CMakeLists.txt | 30 +++++++-- src/gromacs/CMakeLists.txt | 19 +++++- src/gromacs/gpu_utils/CMakeLists.txt | 7 ++- src/gromacs/gpu_utils/clfftinitializer.cpp | 89 +++++++++++++++++++++++++++ src/gromacs/gpu_utils/clfftinitializer.h | 99 ++++++++++++++++++++++++++++++ src/gromacs/mdrun/runner.cpp | 12 +++- 6 files changed, 248 insertions(+), 8 deletions(-) create mode 100644 src/gromacs/gpu_utils/clfftinitializer.cpp create mode 100644 src/gromacs/gpu_utils/clfftinitializer.h diff --git a/src/external/clFFT/src/CMakeLists.txt b/src/external/clFFT/src/CMakeLists.txt index baad15159d..b5bcdd87b5 100644 --- a/src/external/clFFT/src/CMakeLists.txt +++ b/src/external/clFFT/src/CMakeLists.txt @@ -5,10 +5,30 @@ file(GLOB CLFFT_SOURCES statTimer/*.cpp library/*.cpp library/*.c) file(GLOB CLFFT_SOURCES_UNWANTED statTimer/dllmain.cpp library/dllmain.cpp) list(REMOVE_ITEM CLFFT_SOURCES ${CLFFT_SOURCES_UNWANTED}) +add_library(clFFT OBJECT ${CLFFT_SOURCES}) + # disable all warnings -set_source_files_properties(${CLFFT_SOURCES} PROPERTIES COMPILE_FLAGS "${COMPILE_FLAGS} -w") +target_compile_options(clFFT PRIVATE "-w") +if (BUILD_SHARED_LIBS) + set_target_properties(clFFT PROPERTIES + POSITION_INDEPENDENT_CODE ON + ) +endif() + +# Windows doesn't need anything special to load the dynamic libraries +# that the AMD clFFT implementation uses, but *nix and BSD do. +if (NOT WIN32) + include(CheckCXXSymbolExists) + check_cxx_symbol_exists(dlopen dlfcn.h HAVE_DLOPEN) + if (NOT HAVE_DLOPEN) + message(FATAL_ERROR "Cannot use clFFT for OpenCL support unless dlopen is available") + endif() +endif() -# link -add_library(clFFT STATIC ${CLFFT_SOURCES}) -target_include_directories(clFFT PUBLIC ${CMAKE_CURRENT_SOURCE_DIR}/include) -set_property(TARGET clFFT PROPERTY POSITION_INDEPENDENT_CODE ON) +target_include_directories(clFFT + PUBLIC + ${CMAKE_CURRENT_SOURCE_DIR}/include + PRIVATE + ${CMAKE_CURRENT_SOURCE_DIR}/library + ${CMAKE_CURRENT_SOURCE_DIR}/statTimer + ) diff --git a/src/gromacs/CMakeLists.txt b/src/gromacs/CMakeLists.txt index b84d783323..52bc8e1900 100644 --- a/src/gromacs/CMakeLists.txt +++ b/src/gromacs/CMakeLists.txt @@ -196,6 +196,22 @@ else() add_library(libgromacs ${LIBGROMACS_SOURCES}) endif() +# TODO: installed clFFT should be preferred, if possible; +# src/external/src/clFFT/FindclFFT.cmake should help with that. +# Redmine issue #2500 +if (GMX_USE_OPENCL) + set (_clFFT_dir ../external/clFFT/src) + add_subdirectory(${_clFFT_dir} clFFT-build) + target_sources(libgromacs PRIVATE + $ + ) + target_include_directories(libgromacs PRIVATE ${_clFFT_dir}/include) + # Use the magic variable for how to link any library needed for + # dlopen, etc. which is -ldl where needed, and empty otherwise + # (e.g. Windows, BSD, Mac). + target_link_libraries(libgromacs PRIVATE "${CMAKE_DL_LIBS}") +endif() + # Recent versions of gcc and clang give warnings on scanner.cpp, which # is a generated source file. These are awkward to suppress inline, so # we do it in the compilation command (after testing that the compiler @@ -221,7 +237,8 @@ target_link_libraries(libgromacs ${GMX_COMMON_LIBRARIES} ${FFT_LIBRARIES} ${LINEAR_ALGEBRA_LIBRARIES} ${LMFIT_LIBRARIES_TO_LINK} - ${THREAD_LIB} ${GMX_SHARED_LINKER_FLAGS} ${OpenCL_LIBRARIES} + ${THREAD_LIB} ${GMX_SHARED_LINKER_FLAGS} + ${OpenCL_LIBRARIES} ${GMX_STDLIB_LIBRARIES} PUBLIC ${GMX_PUBLIC_LIBRARIES} diff --git a/src/gromacs/gpu_utils/CMakeLists.txt b/src/gromacs/gpu_utils/CMakeLists.txt index ec822afe8a..1ee6f7d5ad 100644 --- a/src/gromacs/gpu_utils/CMakeLists.txt +++ b/src/gromacs/gpu_utils/CMakeLists.txt @@ -1,7 +1,7 @@ # # This file is part of the GROMACS molecular simulation package. # -# Copyright (c) 2015,2016,2017, by the GROMACS development team, led by +# Copyright (c) 2015,2016,2017,2018, 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. @@ -32,7 +32,12 @@ # To help us fund GROMACS development, we humbly ask that you cite # the research papers on the package. Check out http://www.gromacs.org. +# Note, source files should be listed in alphabetical order, to reduce +# the incidence of textual clashes when adding/moving files that +# otherwise make the end of the list a hotspot. + gmx_add_libgromacs_sources( + clfftinitializer.cpp hostallocator.cpp ) if(GMX_USE_OPENCL) diff --git a/src/gromacs/gpu_utils/clfftinitializer.cpp b/src/gromacs/gpu_utils/clfftinitializer.cpp new file mode 100644 index 0000000000..fab5285ab5 --- /dev/null +++ b/src/gromacs/gpu_utils/clfftinitializer.cpp @@ -0,0 +1,89 @@ +/* + * This file is part of the GROMACS molecular simulation package. + * + * Copyright (c) 2018, 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 + * Implements ClfftInitializer class. + * + * \author Aleksei Iupinov + */ + +#include "gmxpre.h" + +#include "clfftinitializer.h" + +#include "config.h" + +#include "gromacs/compat/make_unique.h" +#include "gromacs/utility/exceptions.h" +#include "gromacs/utility/stringutil.h" + +#if GMX_GPU == GMX_GPU_OPENCL +#include +#endif + +namespace gmx +{ + +ClfftInitializer::ClfftInitializer() +{ +#if GMX_GPU == GMX_GPU_OPENCL + clfftSetupData fftSetup; + int initErrorCode = clfftInitSetupData(&fftSetup); + if (initErrorCode != 0) + { + GMX_THROW(InternalError(formatString("Failed to initialize the clFFT library, error code %d", initErrorCode))); + } + initErrorCode = clfftSetup(&fftSetup); + if (initErrorCode != 0) + { + GMX_THROW(InternalError(formatString("Failed to initialize the clFFT library, error code %d", initErrorCode))); + } +#endif +} + +ClfftInitializer::~ClfftInitializer() +{ +#if GMX_GPU == GMX_GPU_OPENCL + clfftTeardown(); + // TODO: log non-0 return values (errors) +#endif +} + +std::unique_ptr initializeClfftLibrary() +{ + return compat::make_unique(); +} + +} diff --git a/src/gromacs/gpu_utils/clfftinitializer.h b/src/gromacs/gpu_utils/clfftinitializer.h new file mode 100644 index 0000000000..80c48f4bc9 --- /dev/null +++ b/src/gromacs/gpu_utils/clfftinitializer.h @@ -0,0 +1,99 @@ +/* + * This file is part of the GROMACS molecular simulation package. + * + * Copyright (c) 2018, 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 + * Declares ClfftInitializer class, which initializes and + * tears down the clFFT library resources in OpenCL builds, + * and does nothing in other builds, and a factory function for it. + * + * clFFT itself is used in the OpenCL implementation of PME + * for 3D R2C/C2R transforms. It is know to work with NVidia + * OpenCL, AMD fglrx and AMDGPU-PRO drivers, and to not work with + * AMD Rocm dev preview as of May 2018 (#2515). + * TODO: find out compatibility with Intel once the rest of PME + * gets there (#2516), or by building clFFT own tests. + * + * \inlibraryapi + * \author Aleksei Iupinov + */ +#ifndef GMX_GPU_UTILS_CLFFTINITIALIZER_H +#define GMX_GPU_UTILS_CLFFTINITIALIZER_H + +#include + +#include "gromacs/utility/classhelpers.h" + +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. + * + * 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. + * + * \todo Consider making a composite object that also handles + * on-demand compilation, managing lifetime of PME FFT kernel programs + * to avoid exhausting resources and/or recompiling kernels previously + * used. See Redmine #2535. + * + * \todo Consider implementing an appropriate flavor of the nifty + * 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 initializeClfftLibrary(); + +} + +#endif diff --git a/src/gromacs/mdrun/runner.cpp b/src/gromacs/mdrun/runner.cpp index 3a77c00351..995a69680a 100644 --- a/src/gromacs/mdrun/runner.cpp +++ b/src/gromacs/mdrun/runner.cpp @@ -65,6 +65,7 @@ #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" @@ -1036,7 +1037,9 @@ int Mdrunner::mdrunner() } } - gmx_device_info_t *pmeDeviceInfo = nullptr; + std::unique_ptr 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 // are invoked. @@ -1048,6 +1051,13 @@ int Mdrunner::mdrunner() pmeDeviceInfo = getDeviceInfo(hwinfo->gpu_info, pmeGpuTaskMapping->deviceId_); init_gpu(mdlog, 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 -- 2.11.4.GIT