Don't allow multiple energy groups for GPU runs
authorErik Lindahl <erik@kth.se>
Sun, 31 Dec 2017 16:24:23 +0000 (17:24 +0100)
committerErik Lindahl <erik.lindahl@gmail.com>
Wed, 3 Jan 2018 12:26:05 +0000 (13:26 +0100)
Exit with a fatal error instead of only warning, since the
latter leads to writing data for energy groups that
is incorrect to the energy file.

Fixes #1822.

Change-Id: I34ccb10bba6d6e1350283e34ebc908c6f830baab

src/gromacs/mdlib/nbnxn_atomdata.cpp
src/gromacs/taskassignment/decidegpuusage.cpp
src/programs/mdrun/runner.cpp

index d6c6be8e23e6e76c2ad4154b695e7d73da5361d5..42a5a244142315bb9837a82d77809f1f33d7da58 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,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.
@@ -701,12 +701,8 @@ void nbnxn_atomdata_init(const gmx::MDLogger &mdlog,
     nbat->nenergrp = n_energygroups;
     if (!simple)
     {
-        /* Energy groups not supported yet for super-sub lists */
-        if (n_energygroups > 1)
-        {
-            GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: With GPUs, reporting energy group contributions is not supported");
-        }
-        nbat->nenergrp = 1;
+        // We now check for energy groups already when starting mdrun
+        GMX_RELEASE_ASSERT(n_energygroups == 1, "GPU kernels do not support energy groups");
     }
     /* Temporary storage goes as #grp^3*simd_width^2/2, so limit to 64 */
     if (nbat->nenergrp > 64)
index d5132add212e7c986de4c39332978081643d1468..1e397586a27addbf15ff57d6325cf77baebdcb83 100644 (file)
@@ -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.
@@ -288,8 +288,8 @@ bool decideWhetherToUseGpusForNonbonded(const TaskTarget           nonbondedTarg
         if (nonbondedTarget == TaskTarget::Gpu)
         {
             GMX_THROW(InconsistentInputError
-                          ("Nonbonded interactions on the GPU were required, but this would not be "
-                          "useful. Probably you should not require using GPUs."));
+                          ("Nonbonded interactions on the GPU were required, but not supported for these "
+                          "simulation settings. Change your settings, or do not require using GPUs."));
         }
 
         return false;
index ab416babb52082e15ca0503d23348bfae8501441..b96cf345fdd24a6f819f6c17a3edc9f05fcd58ae 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2011,2012,2013,2014,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.
@@ -316,27 +316,23 @@ static void override_nsteps_cmdline(const gmx::MDLogger &mdlog,
 namespace gmx
 {
 
-/*! \brief Return whether GPU acceleration of nonbondeds is useful with the given settings.
+/*! \brief Return whether GPU acceleration of nonbondeds is supported with the given settings.
  *
- * If not, logs a message about falling back to CPU code. */
-static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger   &mdlog,
-                                               const t_inputrec *ir,
-                                               bool              doRerun)
+ * If not, logs a message about falling back to CPU code.
+ */
+static bool gpuAccelerationOfNonbondedIsUseful(const MDLogger     &mdlog,
+                                               const t_inputrec *  ir)
 {
-    if (doRerun && ir->opts.ngener > 1)
-    {
-        /* Rerun execution time is dominated by I/O and pair search,
-         * so GPUs are not very useful, plus they do not support more
-         * than one energy group. If the user requested GPUs
-         * explicitly, a fatal error is given later.  With non-reruns,
-         * we fall back to a single whole-of system energy group
-         * (which runs much faster than a multiple-energy-groups
-         * implementation would), and issue a note in the .log
-         * file. Users can re-run if they want the information. */
-        GMX_LOG(mdlog.warning).asParagraph().appendText("Multiple energy groups is not implemented for GPUs, so is not useful for this rerun, so falling back to the CPU");
+    if (ir->opts.ngener > 1)
+    {
+        /* The GPU code does not support more than one energy group.
+         * If the user requested GPUs explicitly, a fatal error is given later.
+         */
+        GMX_LOG(mdlog.warning).asParagraph().appendText("Multiple energy groups is not implemented for GPUs, falling back to the CPU. "
+                                                        "For better performance, run on the GPU without energy groups and then do "
+                                                        "gmx mdrun -rerun option on the trajectory with an energy group .tpr file.");
         return false;
     }
-
     return true;
 }
 
@@ -594,15 +590,17 @@ int Mdrunner::mdrunner()
             useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
                     (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
                     inputrec->cutoff_scheme == ecutsVERLET,
-                    gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, doRerun),
+                    gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec),
                     hw_opt.nthreads_tmpi);
             auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
             auto canUseGpuForPme   = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
                     (useGpuForNonbonded, pmeTarget, gpuIdsToUse, userGpuTaskAssignment,
                     canUseGpuForPme, hw_opt.nthreads_tmpi, domdecOptions.numPmeRanks);
+
         }
         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
+
         /* Determine how many thread-MPI ranks to start.
          *
          * TODO Over-writing the user-supplied value here does
@@ -654,13 +652,14 @@ int Mdrunner::mdrunner()
         bool gpusWereDetected = hwinfo->ngpu_compatible_tot > 0;
         useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
                                                                 emulateGpuNonbonded, inputrec->cutoff_scheme == ecutsVERLET,
-                                                                gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, doRerun),
+                                                                gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec),
                                                                 gpusWereDetected);
         auto inputSystemHasPme = EEL_PME(inputrec->coulombtype) || EVDW_PME(inputrec->vdwtype);
         auto canUseGpuForPme   = inputSystemHasPme && pme_gpu_supports_input(inputrec, nullptr);
         useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
                                                     canUseGpuForPme, cr->nnodes, domdecOptions.numPmeRanks,
                                                     gpusWereDetected);
+
         pmeRunMode   = (useGpuForPme ? PmeRunMode::GPU : PmeRunMode::CPU);
         if (pmeRunMode == PmeRunMode::GPU)
         {