Change default Ewald kernel to tabulated on NVIDIA CC 7.0
authorSzilárd Páll <pall.szilard@gmail.com>
Tue, 26 Jan 2021 16:33:16 +0000 (16:33 +0000)
committerPaul Bauer <paul.bauer.q@gmail.com>
Tue, 26 Jan 2021 16:33:16 +0000 (16:33 +0000)
Partially addressed #3845

docs/release-notes/2021/major/performance.rst
src/gromacs/nbnxm/cuda/nbnxm_cuda_data_mgmt.cu
src/gromacs/nbnxm/nbnxm_gpu_data_mgmt.cpp
src/gromacs/nbnxm/nbnxm_gpu_data_mgmt.h
src/gromacs/nbnxm/opencl/nbnxm_ocl_data_mgmt.cpp

index 9c499fde539986941fdef5dc8753941eea374fba..06888c99c9b90e743915afdaa7e0bd861085d034 100644 (file)
@@ -46,3 +46,10 @@ Allow offloading GPU update and constraints without direct GPU communication
 Allow domain-decomposition and separate PME rank parallel runs to offload update and
 constraints to a GPU with CUDA without requiring the (experimental) direct GPU
 communication features to be also enabled.
+
+Tune CUDA short-range nonbonded kernel parameters on NVIDIA Volta and Ampere A100
+"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+Recent compilers allowed re-tuning the nonbonded kernel defaults on NVIDIA Volta and
+Ampere A100GPUs which improves performance of the Ewald kernels, especially those that
+also compute energies.
index b1d6774a26eb329362fc11cee409be22644c11a8..bedabd85c8bb7eef62d8596c6ccd25c6a526fca7 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
- * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2017,2018,2019,2020,2021, 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.
@@ -193,7 +193,7 @@ static void init_nbparam(NBParamGpu*                     nbp,
     }
     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
     {
-        nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
+        nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic, deviceContext.deviceInfo());
     }
     else
     {
index 0c3586f7d68e510a7d4611e7820a8a8399470095..0af6e94e6189a623768d0eb334890e2c52834b67 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
- * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2017,2018,2019,2020,2021, 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.
@@ -58,6 +58,7 @@
 
 #include "nbnxm_gpu_data_mgmt.h"
 
+#include "gromacs/hardware/device_information.h"
 #include "gromacs/nbnxm/gpu_data_mgmt.h"
 #include "gromacs/timing/gpu_timing.h"
 #include "gromacs/utility/cstringutil.h"
@@ -95,7 +96,8 @@ void inline printEnvironmentVariableDeprecationMessage(bool               isEnvi
     }
 }
 
-int nbnxn_gpu_pick_ewald_kernel_type(const interaction_const_t& ic)
+int nbnxn_gpu_pick_ewald_kernel_type(const interaction_const_t& ic,
+                                     const DeviceInformation gmx_unused& deviceInfo)
 {
     bool bTwinCut = (ic.rcoulomb != ic.rvdw);
     int  kernel_type;
@@ -129,13 +131,19 @@ int nbnxn_gpu_pick_ewald_kernel_type(const interaction_const_t& ic)
                 "requested through environment variables.");
     }
 
-    /* By default, use analytical Ewald
-     * TODO: tabulated does not work in OpenCL, it needs fixing, see init_nbparam() in nbnxn_ocl_data_mgmt.cpp
-     *
+    /* By default, use analytical Ewald except with CUDA on NVIDIA CC 7.0 and 8.0.
      */
-    bool bUseAnalyticalEwald = true;
+    const bool c_useTabulatedEwaldDefault =
+#if GMX_GPU_CUDA
+            (deviceInfo.prop.major == 7 && deviceInfo.prop.minor == 0)
+            || (deviceInfo.prop.major == 8 && deviceInfo.prop.minor == 0);
+#else
+            false;
+#endif
+    bool bUseAnalyticalEwald = !c_useTabulatedEwaldDefault;
     if (forceAnalyticalEwald)
     {
+        bUseAnalyticalEwald = true;
         if (debug)
         {
             fprintf(debug, "Using analytical Ewald GPU kernels\n");
@@ -198,7 +206,7 @@ void gpu_pme_loadbal_update_param(const nonbonded_verlet_t* nbv, const interacti
 
     set_cutoff_parameters(nbp, ic, nbv->pairlistSets().params());
 
-    nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
+    nbp->eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic, nb->deviceContext_->deviceInfo());
 
     GMX_RELEASE_ASSERT(ic->coulombEwaldTables, "Need valid Coulomb Ewald correction tables");
     init_ewald_coulomb_force_table(*ic->coulombEwaldTables, nbp, *nb->deviceContext_);
index 761737ddf0c50a2f19a1c335796bb4e0434ac85c..b417566db4fd96e10356b33e7a3877e42fc10d4f 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,2020,2021, 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.
@@ -68,7 +68,8 @@ void init_ewald_coulomb_force_table(const EwaldCorrectionTables& tables,
                                     const DeviceContext&         deviceContext);
 
 /*! \brief Selects the Ewald kernel type, analytical or tabulated, single or twin cut-off. */
-int nbnxn_gpu_pick_ewald_kernel_type(const interaction_const_t gmx_unused& ic);
+int nbnxn_gpu_pick_ewald_kernel_type(const interaction_const_t gmx_unused& ic,
+                                     const DeviceInformation&              deviceInfo);
 
 /*! \brief Copies all parameters related to the cut-off from ic to nbp
  */
index 29989c8095b5bbf0a0117b39902e223c45739db8..ac99cf926d7a94b25427712ac54888df86e4bcda 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
- * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2017,2018,2019,2020,2021, 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.
@@ -132,7 +132,8 @@ static void init_atomdata_first(cl_atomdata_t* ad, int ntypes, const DeviceConte
 static void map_interaction_types_to_gpu_kernel_flavors(const interaction_const_t* ic,
                                                         int                        combRule,
                                                         int*                       gpu_eeltype,
-                                                        int*                       gpu_vdwtype)
+                                                        int*                       gpu_vdwtype,
+                                                        const DeviceContext&       deviceContext)
 {
     if (ic->vdwtype == evdwCUT)
     {
@@ -185,7 +186,7 @@ static void map_interaction_types_to_gpu_kernel_flavors(const interaction_const_
     }
     else if ((EEL_PME(ic->eeltype) || ic->eeltype == eelEWALD))
     {
-        *gpu_eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic);
+        *gpu_eeltype = nbnxn_gpu_pick_ewald_kernel_type(*ic, deviceContext.deviceInfo());
     }
     else
     {
@@ -206,7 +207,8 @@ static void init_nbparam(NBParamGpu*                     nbp,
 {
     set_cutoff_parameters(nbp, ic, listParams);
 
-    map_interaction_types_to_gpu_kernel_flavors(ic, nbatParams.comb_rule, &(nbp->eeltype), &(nbp->vdwtype));
+    map_interaction_types_to_gpu_kernel_flavors(ic, nbatParams.comb_rule, &(nbp->eeltype),
+                                                &(nbp->vdwtype), deviceContext);
 
     if (ic->vdwtype == evdwPME)
     {