Introduce GMX_USE_SIMD_KERNELS cmake option
authorMark Abraham <mark.j.abraham@gmail.com>
Thu, 16 Apr 2020 15:36:33 +0000 (15:36 +0000)
committerPaul Bauer <paul.bauer.q@gmail.com>
Thu, 16 Apr 2020 15:36:33 +0000 (15:36 +0000)
Most GROMACS development does not need to recompile the SIMD nbnxm
(and fep) kernels whenever their dependencies change. These
dependencies are large in number, and include frequently changed
files, including config.h and various utility and nbnxm module
headers. This flag permits people to efficiently recompile while
working on code that doesn't directly target changes to the SIMD
kernels.

It also means that CI builds not aimed at efficient mdrun execution
times can instead minimize compilation times and ccache db sizes. This
will also have the side effect of testing more of the reference NBNXM
kernels.

There's other SIMD code (particularly PME, bonded, LINCS, update)
which still compiles and runs in the usual way. Currently these are
less costly to compile and harder to disable. That could change in
future.

src/config.h.cmakein
src/gromacs/CMakeLists.txt
src/gromacs/gmxlib/nonbonded/nb_free_energy.cpp
src/gromacs/mdlib/forcerec.cpp
src/gromacs/nbnxm/CMakeLists.txt
src/gromacs/nbnxm/kernels_simd_2xmm/CMakeLists.txt [new file with mode: 0644]
src/gromacs/nbnxm/kernels_simd_4xm/CMakeLists.txt [new file with mode: 0644]
src/gromacs/nbnxm/nbnxm_simd.h
src/gromacs/nbnxm/prunekerneldispatch.cpp

index 5def31c6e3c90d268d0aa373359202666f207193..9be4bb51d1d6d311e1ae088c70b5429409c2449f 100644 (file)
 /* Enable code that requires AVX-512 instruction support, without GMX_SIMD=AVX_512 */
 #cmakedefine01 SIMD_AVX_512_CXX_SUPPORTED
 
+/* Whether NBNXM and other SIMD kernels should be compiled */
+#cmakedefine01 GMX_USE_SIMD_KERNELS
+
 /* Whether a double-precision configuration may target accuracy equivalent to single precision */
 #cmakedefine01 GMX_RELAXED_DOUBLE_PRECISION
 
index 66a547b647ba506d51e6100a454b15d39c29f6d3..ea7176bf6952771073d44e03090f5fd3cdaa3e54 100644 (file)
@@ -271,6 +271,14 @@ target_include_directories(libgromacs PUBLIC
                            $<INSTALL_INTERFACE:include>
                            )
 
+# Permit the configuration to disable compiling the many nbnxm kernels
+# and others involved in force calculations. Currently only
+# short-ranged and bonded kernels are disabled this way, but in future
+# others may be appropriate. Thus the cmake option is not specific to
+# nbnxm module.
+option(GMX_USE_SIMD_KERNELS "Whether to compile NBNXM and other SIMD kernels" ON)
+mark_as_advanced(GMX_USE_SIMD_KERNELS)
+
 if(SIMD_AVX_512_CXX_SUPPORTED AND NOT ("${GMX_SIMD_ACTIVE}" STREQUAL "AVX_512_KNL"))
     # Since we might be overriding -march=core-avx2, add a flag so we don't warn for this specific file.
     # On KNL this can cause illegal instruction because the compiler might use non KNL AVX instructions
index e84bdaf037d2b7f1fb32a71c7113ebed0229385b..8155e5beaed6b4253c7644c1fd56a7afbf4aa90d 100644 (file)
@@ -833,7 +833,7 @@ static KernelFunction dispatchKernelOnUseSimd(const bool useSimd)
 {
     if (useSimd)
     {
-#if GMX_SIMD_HAVE_REAL && GMX_SIMD_HAVE_INT32_ARITHMETICS
+#if GMX_SIMD_HAVE_REAL && GMX_SIMD_HAVE_INT32_ARITHMETICS && GMX_USE_SIMD_KERNELS
         /* FIXME: Here SimdDataTypes should be used to enable SIMD. So far, the code in nb_free_energy_kernel is not adapted to SIMD */
         return (nb_free_energy_kernel<ScalarDataTypes, useSoftCore, scLambdasOrAlphasDiffer, vdwInteractionTypeIsEwald,
                                       elecInteractionTypeIsEwald, vdwModifierIsPotSwitch>);
index 7f78c142baf6ea5679f1c60044e5ff29ad79bbb9..2f771943f6b04566c5469b4b764633cfd5080855 100644 (file)
@@ -38,6 +38,8 @@
 
 #include "forcerec.h"
 
+#include "config.h"
+
 #include <cassert>
 #include <cmath>
 #include <cstdlib>
@@ -926,8 +928,8 @@ void init_forcerec(FILE*                            fp,
                    gmx::ArrayRef<const std::string> tabbfnm,
                    real                             print_force)
 {
-    /* By default we turn SIMD kernels on, but it might be turned off further down... */
-    fr->use_simd_kernels = TRUE;
+    /* The CMake default turns SIMD kernels on, but it might be turned off further down... */
+    fr->use_simd_kernels = GMX_USE_SIMD_KERNELS;
 
     if (check_box(ir->pbcType, box))
     {
index 087c2696a4af2541594d537c23eaae5f9ccb2176..ecfb452e4f4edc4e8afef978324c43b2ec4467ab 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2019, by the GROMACS development team, led by
+# Copyright (c) 2019,2020, 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.
 # To help us fund GROMACS development, we humbly ask that you cite
 # the research papers on the package. Check out http://www.gromacs.org.
 
-file(GLOB NBNXM_SOURCES kernels_reference/*.cpp kernels_simd_4xm/*.cpp kernels_simd_2xmm/*.cpp *.cpp benchmark/*.cpp)
+add_subdirectory(kernels_simd_4xm)
+add_subdirectory(kernels_simd_2xmm)
+
+file (GLOB NBNXM_SOURCES
+    # Source files
+    atomdata.cpp
+    grid.cpp
+    gridset.cpp
+    kernel_common.cpp
+    kerneldispatch.cpp
+    nbnxm.cpp
+    nbnxm_geometry.cpp
+    nbnxm_setup.cpp
+    pairlist.cpp
+    pairlistparams.cpp
+    pairlistset.cpp
+    pairlist_tuning.cpp
+    pairsearch.cpp
+    prunekerneldispatch.cpp
+    # Reference kernel source files
+    kernels_reference/kernel_gpu_ref.cpp
+    kernels_reference/kernel_ref.cpp
+    kernels_reference/kernel_ref_prune.cpp
+    # Benchmark source files
+    # TODO these should not be in libgromacs
+    benchmark/bench_setup.cpp
+    benchmark/bench_system.cpp
+    )
 
 if(GMX_USE_CUDA)
     add_subdirectory(cuda)
diff --git a/src/gromacs/nbnxm/kernels_simd_2xmm/CMakeLists.txt b/src/gromacs/nbnxm/kernels_simd_2xmm/CMakeLists.txt
new file mode 100644 (file)
index 0000000..1896bfc
--- /dev/null
@@ -0,0 +1,131 @@
+#
+# This file is part of the GROMACS molecular simulation package.
+#
+# Copyright (c) 2020, 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.
+
+if (GMX_USE_SIMD_KERNELS)
+    file(GLOB KERNEL_SOURCES
+        kernel_ElecEwTwinCut_VdwLJCombGeom_F.cpp
+        kernel_ElecEwTwinCut_VdwLJCombGeom_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJCombLB_F.cpp
+        kernel_ElecEwTwinCut_VdwLJCombLB_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJEwCombGeom_F.cpp
+        kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJ_F.cpp
+        kernel_ElecEwTwinCut_VdwLJFSw_F.cpp
+        kernel_ElecEwTwinCut_VdwLJFSw_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJFSw_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJPSw_F.cpp
+        kernel_ElecEwTwinCut_VdwLJPSw_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJPSw_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJ_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJ_VgrpF.cpp
+        kernel_ElecEw_VdwLJCombGeom_F.cpp
+        kernel_ElecEw_VdwLJCombGeom_VF.cpp
+        kernel_ElecEw_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecEw_VdwLJCombLB_F.cpp
+        kernel_ElecEw_VdwLJCombLB_VF.cpp
+        kernel_ElecEw_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecEw_VdwLJEwCombGeom_F.cpp
+        kernel_ElecEw_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecEw_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecEw_VdwLJ_F.cpp
+        kernel_ElecEw_VdwLJFSw_F.cpp
+        kernel_ElecEw_VdwLJFSw_VF.cpp
+        kernel_ElecEw_VdwLJFSw_VgrpF.cpp
+        kernel_ElecEw_VdwLJPSw_F.cpp
+        kernel_ElecEw_VdwLJPSw_VF.cpp
+        kernel_ElecEw_VdwLJPSw_VgrpF.cpp
+        kernel_ElecEw_VdwLJ_VF.cpp
+        kernel_ElecEw_VdwLJ_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombGeom_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombGeom_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombLB_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombLB_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJ_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJFSw_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJFSw_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJFSw_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJPSw_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJPSw_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJPSw_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJ_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJ_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJCombGeom_F.cpp
+        kernel_ElecQSTab_VdwLJCombGeom_VF.cpp
+        kernel_ElecQSTab_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJCombLB_F.cpp
+        kernel_ElecQSTab_VdwLJCombLB_VF.cpp
+        kernel_ElecQSTab_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJEwCombGeom_F.cpp
+        kernel_ElecQSTab_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecQSTab_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJ_F.cpp
+        kernel_ElecQSTab_VdwLJFSw_F.cpp
+        kernel_ElecQSTab_VdwLJFSw_VF.cpp
+        kernel_ElecQSTab_VdwLJFSw_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJPSw_F.cpp
+        kernel_ElecQSTab_VdwLJPSw_VF.cpp
+        kernel_ElecQSTab_VdwLJPSw_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJ_VF.cpp
+        kernel_ElecQSTab_VdwLJ_VgrpF.cpp
+        kernel_ElecRF_VdwLJCombGeom_F.cpp
+        kernel_ElecRF_VdwLJCombGeom_VF.cpp
+        kernel_ElecRF_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecRF_VdwLJCombLB_F.cpp
+        kernel_ElecRF_VdwLJCombLB_VF.cpp
+        kernel_ElecRF_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecRF_VdwLJEwCombGeom_F.cpp
+        kernel_ElecRF_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecRF_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecRF_VdwLJ_F.cpp
+        kernel_ElecRF_VdwLJFSw_F.cpp
+        kernel_ElecRF_VdwLJFSw_VF.cpp
+        kernel_ElecRF_VdwLJFSw_VgrpF.cpp
+        kernel_ElecRF_VdwLJPSw_F.cpp
+        kernel_ElecRF_VdwLJPSw_VF.cpp
+        kernel_ElecRF_VdwLJPSw_VgrpF.cpp
+        kernel_ElecRF_VdwLJ_VF.cpp
+        kernel_ElecRF_VdwLJ_VgrpF.cpp
+        kernel_prune.cpp
+        )
+endif()
+
+set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${KERNEL_SOURCES} PARENT_SCOPE)
diff --git a/src/gromacs/nbnxm/kernels_simd_4xm/CMakeLists.txt b/src/gromacs/nbnxm/kernels_simd_4xm/CMakeLists.txt
new file mode 100644 (file)
index 0000000..1896bfc
--- /dev/null
@@ -0,0 +1,131 @@
+#
+# This file is part of the GROMACS molecular simulation package.
+#
+# Copyright (c) 2020, 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.
+
+if (GMX_USE_SIMD_KERNELS)
+    file(GLOB KERNEL_SOURCES
+        kernel_ElecEwTwinCut_VdwLJCombGeom_F.cpp
+        kernel_ElecEwTwinCut_VdwLJCombGeom_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJCombLB_F.cpp
+        kernel_ElecEwTwinCut_VdwLJCombLB_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJEwCombGeom_F.cpp
+        kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJ_F.cpp
+        kernel_ElecEwTwinCut_VdwLJFSw_F.cpp
+        kernel_ElecEwTwinCut_VdwLJFSw_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJFSw_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJPSw_F.cpp
+        kernel_ElecEwTwinCut_VdwLJPSw_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJPSw_VgrpF.cpp
+        kernel_ElecEwTwinCut_VdwLJ_VF.cpp
+        kernel_ElecEwTwinCut_VdwLJ_VgrpF.cpp
+        kernel_ElecEw_VdwLJCombGeom_F.cpp
+        kernel_ElecEw_VdwLJCombGeom_VF.cpp
+        kernel_ElecEw_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecEw_VdwLJCombLB_F.cpp
+        kernel_ElecEw_VdwLJCombLB_VF.cpp
+        kernel_ElecEw_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecEw_VdwLJEwCombGeom_F.cpp
+        kernel_ElecEw_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecEw_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecEw_VdwLJ_F.cpp
+        kernel_ElecEw_VdwLJFSw_F.cpp
+        kernel_ElecEw_VdwLJFSw_VF.cpp
+        kernel_ElecEw_VdwLJFSw_VgrpF.cpp
+        kernel_ElecEw_VdwLJPSw_F.cpp
+        kernel_ElecEw_VdwLJPSw_VF.cpp
+        kernel_ElecEw_VdwLJPSw_VgrpF.cpp
+        kernel_ElecEw_VdwLJ_VF.cpp
+        kernel_ElecEw_VdwLJ_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombGeom_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombGeom_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombLB_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombLB_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJ_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJFSw_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJFSw_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJFSw_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJPSw_F.cpp
+        kernel_ElecQSTabTwinCut_VdwLJPSw_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJPSw_VgrpF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJ_VF.cpp
+        kernel_ElecQSTabTwinCut_VdwLJ_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJCombGeom_F.cpp
+        kernel_ElecQSTab_VdwLJCombGeom_VF.cpp
+        kernel_ElecQSTab_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJCombLB_F.cpp
+        kernel_ElecQSTab_VdwLJCombLB_VF.cpp
+        kernel_ElecQSTab_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJEwCombGeom_F.cpp
+        kernel_ElecQSTab_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecQSTab_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJ_F.cpp
+        kernel_ElecQSTab_VdwLJFSw_F.cpp
+        kernel_ElecQSTab_VdwLJFSw_VF.cpp
+        kernel_ElecQSTab_VdwLJFSw_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJPSw_F.cpp
+        kernel_ElecQSTab_VdwLJPSw_VF.cpp
+        kernel_ElecQSTab_VdwLJPSw_VgrpF.cpp
+        kernel_ElecQSTab_VdwLJ_VF.cpp
+        kernel_ElecQSTab_VdwLJ_VgrpF.cpp
+        kernel_ElecRF_VdwLJCombGeom_F.cpp
+        kernel_ElecRF_VdwLJCombGeom_VF.cpp
+        kernel_ElecRF_VdwLJCombGeom_VgrpF.cpp
+        kernel_ElecRF_VdwLJCombLB_F.cpp
+        kernel_ElecRF_VdwLJCombLB_VF.cpp
+        kernel_ElecRF_VdwLJCombLB_VgrpF.cpp
+        kernel_ElecRF_VdwLJEwCombGeom_F.cpp
+        kernel_ElecRF_VdwLJEwCombGeom_VF.cpp
+        kernel_ElecRF_VdwLJEwCombGeom_VgrpF.cpp
+        kernel_ElecRF_VdwLJ_F.cpp
+        kernel_ElecRF_VdwLJFSw_F.cpp
+        kernel_ElecRF_VdwLJFSw_VF.cpp
+        kernel_ElecRF_VdwLJFSw_VgrpF.cpp
+        kernel_ElecRF_VdwLJPSw_F.cpp
+        kernel_ElecRF_VdwLJPSw_VF.cpp
+        kernel_ElecRF_VdwLJPSw_VgrpF.cpp
+        kernel_ElecRF_VdwLJ_VF.cpp
+        kernel_ElecRF_VdwLJ_VgrpF.cpp
+        kernel_prune.cpp
+        )
+endif()
+
+set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${KERNEL_SOURCES} PARENT_SCOPE)
index a81c52fdd29be68092b1ad07560877ba6fe0eb31..53fa3b9606992c29f94909be8345fea14518436b 100644 (file)
@@ -48,7 +48,7 @@
 #include "gromacs/simd/simd.h"
 #include "gromacs/utility/real.h"
 
-#if GMX_SIMD
+#if GMX_SIMD && GMX_USE_SIMD_KERNELS
 /*! \brief The nbnxn SIMD 4xN and 2x(N+N) kernels can be added independently.
  * Currently the 2xNN SIMD kernels only make sense with:
  *  8-way SIMD: 4x4 setup, works with AVX-256 in single precision
@@ -65,6 +65,6 @@
 #        error "No SIMD kernel type defined"
 #    endif
 
-#endif // GMX_SIMD
+#endif // GMX_SIMD && GMX_USE_SIMD_KERNELS
 
 #endif
index 574d9c143ae851755d994a7b3461f23b57a10e11..7b55de41379986d5b36adca6585330da4f536940 100644 (file)
@@ -42,6 +42,7 @@
 
 #include "clusterdistancekerneltype.h"
 #include "nbnxm_gpu.h"
+#include "nbnxm_simd.h"
 #include "pairlistset.h"
 #include "pairlistsets.h"
 #include "kernels_reference/kernel_ref_prune.h"
@@ -72,12 +73,16 @@ void PairlistSet::dispatchPruneKernel(const nbnxn_atomdata_t* nbat, const rvec*
 
         switch (getClusterDistanceKernelType(params_.pairlistType, *nbat))
         {
+#ifdef GMX_NBNXN_SIMD_4XN
             case ClusterDistanceKernelType::CpuSimd_4xM:
                 nbnxn_kernel_prune_4xn(nbl, nbat, shift_vec, rlistInner);
                 break;
+#endif
+#ifdef GMX_NBNXN_SIMD_2XNN
             case ClusterDistanceKernelType::CpuSimd_2xMM:
                 nbnxn_kernel_prune_2xnn(nbl, nbat, shift_vec, rlistInner);
                 break;
+#endif
             case ClusterDistanceKernelType::CpuPlainC:
                 nbnxn_kernel_prune_ref(nbl, nbat, shift_vec, rlistInner);
                 break;