Remove group scheme checks from runner
authorKevin Boyd <kevin.boyd@uconn.edu>
Sat, 12 Jan 2019 14:47:54 +0000 (09:47 -0500)
committerMark Abraham <mark.j.abraham@gmail.com>
Sun, 30 Jun 2019 21:13:05 +0000 (23:13 +0200)
refs #1852

Change-Id: I906ffc7c063694fbbd17128ee0d1e62259040306

src/gromacs/mdlib/gmx_omp_nthreads.cpp
src/gromacs/mdlib/gmx_omp_nthreads.h
src/gromacs/mdrun/runner.cpp
src/gromacs/taskassignment/decidegpuusage.cpp
src/gromacs/taskassignment/decidegpuusage.h
src/gromacs/taskassignment/resourcedivision.cpp
src/gromacs/taskassignment/resourcedivision.h

index 6c9fa21577e5428679093adf43107a3692a09803..587917990798673f6309c959de1db5c72e64223f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,2015,2016,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.
@@ -102,13 +102,9 @@ static omp_module_nthreads_t modth = { 0, 0, {0, 0, 0, 0, 0, 0, 0, 0, 0}, FALSE}
  *
  *  Each number of threads per module takes the default value unless
  *  GMX_*_NUM_THERADS env var is set, case in which its value overrides
- *  the deafult.
- *
- *  The "group" scheme supports OpenMP only in PME and in thise case all but
- *  the PME nthread values default to 1.
+ *  the default.
  */
 static void pick_module_nthreads(const gmx::MDLogger &mdlog, int m,
-                                 gmx_bool bFullOmpSupport,
                                  gmx_bool bSepPME)
 {
     char      *env;
@@ -137,25 +133,13 @@ static void pick_module_nthreads(const gmx::MDLogger &mdlog, int m,
 
         /* with the verlet codepath, when any GMX_*_NUM_THREADS env var is set,
          * OMP_NUM_THREADS also has to be set */
-        if (bFullOmpSupport && getenv("OMP_NUM_THREADS") == nullptr)
+        if (getenv("OMP_NUM_THREADS") == nullptr)
         {
             gmx_warning("%s=%d is set, the default number of threads also "
                         "needs to be set with OMP_NUM_THREADS!",
                         modth_env_var[m], nth);
         }
 
-        /* with the group scheme warn if any env var except PME is set */
-        if (!bFullOmpSupport)
-        {
-            if (m != emntPME)
-            {
-                gmx_warning("%s=%d is set, but OpenMP multithreading is not "
-                            "supported in %s!",
-                            modth_env_var[m], nth, mod_name[m]);
-                nth = 1;
-            }
-        }
-
         /* only babble if we are really overriding with a different value */
         if ((bSepPME && m == emntPME && nth != modth.gnth_pme) || (nth != modth.gnth))
         {
@@ -231,7 +215,6 @@ static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
                                             int                  omp_nthreads_req,
                                             int                  omp_nthreads_pme_req,
                                             gmx_bool gmx_unused  bThisNodePMEOnly,
-                                            gmx_bool             bFullOmpSupport,
                                             int                  numRanksOnThisNode,
                                             gmx_bool             bSepPME)
 {
@@ -271,8 +254,6 @@ static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
      *   the respective module and it has to be used in conjunction with
      *   OMP_NUM_THREADS.
      *
-     * With the group scheme OpenMP multithreading is only supported in PME,
-     * for all other modules nthreads is set to 1.
      * The number of PME threads is equal to:
      * - 1 if not compiled with OpenMP or
      * - GMX_PME_NUM_THREADS if defined, otherwise
@@ -296,7 +277,7 @@ static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
     {
         nth = omp_nthreads_req;
     }
-    else if (bFullOmpSupport && bOMP)
+    else if (bOMP)
     {
         /* max available threads per node */
         nth = nthreads_hw_avail;
@@ -313,10 +294,10 @@ static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
     }
 
     /* now we have the global values, set them:
-     * - 1 if not compiled with OpenMP and for the group scheme
+     * - 1 if not compiled with OpenMP
      * - nth for the verlet scheme when compiled with OpenMP
      */
-    if (bFullOmpSupport && bOMP)
+    if (bOMP)
     {
         modth.gnth = nth;
     }
@@ -343,15 +324,15 @@ static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
 
     /* now set the per-module values */
     modth.nth[emntDefault] = modth.gnth;
-    pick_module_nthreads(mdlog, emntDomdec, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntPairsearch, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntNonbonded, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntBonded, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntPME, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntUpdate, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntVSITE, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntLINCS, bFullOmpSupport, bSepPME);
-    pick_module_nthreads(mdlog, emntSETTLE, bFullOmpSupport, bSepPME);
+    pick_module_nthreads(mdlog, emntDomdec, bSepPME);
+    pick_module_nthreads(mdlog, emntPairsearch, bSepPME);
+    pick_module_nthreads(mdlog, emntNonbonded, bSepPME);
+    pick_module_nthreads(mdlog, emntBonded, bSepPME);
+    pick_module_nthreads(mdlog, emntPME, bSepPME);
+    pick_module_nthreads(mdlog, emntUpdate, bSepPME);
+    pick_module_nthreads(mdlog, emntVSITE, bSepPME);
+    pick_module_nthreads(mdlog, emntLINCS, bSepPME);
+    pick_module_nthreads(mdlog, emntSETTLE, bSepPME);
 
     /* set the number of threads globally */
     if (bOMP)
@@ -364,14 +345,7 @@ static void manage_number_of_openmp_threads(const gmx::MDLogger &mdlog,
         else
 #endif      /* GMX_THREAD_MPI */
         {
-            if (bFullOmpSupport)
-            {
-                gmx_omp_set_num_threads(nth);
-            }
-            else
-            {
-                gmx_omp_set_num_threads(1);
-            }
+            gmx_omp_set_num_threads(nth);
         }
     }
 
@@ -383,7 +357,6 @@ static void
 reportOpenmpSettings(const gmx::MDLogger &mdlog,
                      const t_commrec     *cr,
                      gmx_bool             bOMP,
-                     gmx_bool             bFullOmpSupport,
                      gmx_bool             bSepPME)
 {
 #if GMX_THREAD_MPI
@@ -426,23 +399,21 @@ reportOpenmpSettings(const gmx::MDLogger &mdlog,
         nth_pme_max = modth.gnth_pme;
     }
 
-    /* for group scheme we print PME threads info only */
-    if (bFullOmpSupport)
+
+    if (nth_max == nth_min)
     {
-        if (nth_max == nth_min)
-        {
-            GMX_LOG(mdlog.warning).appendTextFormatted(
-                    "Using %d OpenMP thread%s %s",
-                    nth_min, nth_min > 1 ? "s" : "",
-                    cr->nnodes > 1 ? mpi_str : "");
-        }
-        else
-        {
-            GMX_LOG(mdlog.warning).appendTextFormatted(
-                    "Using %d - %d OpenMP threads %s",
-                    nth_min, nth_max, mpi_str);
-        }
+        GMX_LOG(mdlog.warning).appendTextFormatted(
+                "Using %d OpenMP thread%s %s",
+                nth_min, nth_min > 1 ? "s" : "",
+                cr->nnodes > 1 ? mpi_str : "");
     }
+    else
+    {
+        GMX_LOG(mdlog.warning).appendTextFormatted(
+                "Using %d - %d OpenMP threads %s",
+                nth_min, nth_max, mpi_str);
+    }
+
     if (bSepPME && (nth_pme_min != nth_min || nth_pme_max != nth_max))
     {
         if (nth_pme_max == nth_pme_min)
@@ -467,8 +438,7 @@ void gmx_omp_nthreads_init(const gmx::MDLogger &mdlog, t_commrec *cr,
                            int numRanksOnThisNode,
                            int omp_nthreads_req,
                            int omp_nthreads_pme_req,
-                           gmx_bool bThisNodePMEOnly,
-                           gmx_bool bFullOmpSupport)
+                           gmx_bool bThisNodePMEOnly)
 {
     gmx_bool   bSepPME;
 
@@ -479,8 +449,7 @@ void gmx_omp_nthreads_init(const gmx::MDLogger &mdlog, t_commrec *cr,
     manage_number_of_openmp_threads(mdlog, cr, bOMP,
                                     nthreads_hw_avail,
                                     omp_nthreads_req, omp_nthreads_pme_req,
-                                    bThisNodePMEOnly, bFullOmpSupport,
-                                    numRanksOnThisNode, bSepPME);
+                                    bThisNodePMEOnly, numRanksOnThisNode, bSepPME);
 #if GMX_THREAD_MPI
     /* Non-master threads have to wait for the OpenMP management to be
      * done, so that code elsewhere that uses OpenMP can be certain
@@ -491,7 +460,7 @@ void gmx_omp_nthreads_init(const gmx::MDLogger &mdlog, t_commrec *cr,
     }
 #endif
 
-    reportOpenmpSettings(mdlog, cr, bOMP, bFullOmpSupport, bSepPME);
+    reportOpenmpSettings(mdlog, cr, bOMP, bSepPME);
 }
 
 int gmx_omp_nthreads_get(int mod)
index e7dd9907b5256e44d1291f4b6926f7e4ff80ffac..5470b463baf2102717f3c06d737a882c65e37518 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,2015,2016,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.
@@ -68,8 +68,7 @@ void gmx_omp_nthreads_init(const gmx::MDLogger &fplog, t_commrec *cr,
                            int numRanksOnThisNode,
                            int omp_nthreads_req,
                            int omp_nthreads_pme_req,
-                           gmx_bool bCurrNodePMEOnly,
-                           gmx_bool bFullOmpSupport);
+                           gmx_bool bCurrNodePMEOnly);
 
 /*! \brief
  * Returns the number of threads to be used in the given module \p mod. */
index b5419b0ac040914b618a5da76aa9d6c4d12f415e..e6835bfc0c24e49c024ff3afec5e667c02775af4 100644 (file)
@@ -641,8 +641,8 @@ int Mdrunner::mdrunner()
         read_tpx_state(ftp2fn(efTPR, filenames.size(), filenames.data()), inputrec, globalState.get(), &mtop);
     }
 
-    // Check and update the hardware options for internal consistency
-    check_and_update_hw_opt_1(mdlog, &hw_opt, SIMMASTER(cr), domdecOptions.numPmeRanks);
+    /* Check and update the hardware options for internal consistency */
+    checkAndUpdateHardwareOptions(mdlog, &hw_opt, SIMMASTER(cr), domdecOptions.numPmeRanks);
 
     if (GMX_THREAD_MPI && SIMMASTER(cr))
     {
@@ -657,7 +657,6 @@ int Mdrunner::mdrunner()
             useGpuForNonbonded = decideWhetherToUseGpusForNonbondedWithThreadMpi
                     (nonbondedTarget, gpuIdsToUse, userGpuTaskAssignment, emulateGpuNonbonded,
                     canUseGpuForNonbonded,
-                    inputrec->cutoff_scheme == ecutsVERLET,
                     gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, GMX_THREAD_MPI),
                     hw_opt.nthreads_tmpi);
             useGpuForPme = decideWhetherToUseGpusForPmeWithThreadMpi
@@ -722,12 +721,10 @@ int Mdrunner::mdrunner()
         // handle. If unsuitable, we will notice that during task
         // assignment.
         bool gpusWereDetected      = hwinfo->ngpu_compatible_tot > 0;
-        bool usingVerletScheme     = inputrec->cutoff_scheme == ecutsVERLET;
         auto canUseGpuForNonbonded = buildSupportsNonbondedOnGpu(nullptr);
         useGpuForNonbonded = decideWhetherToUseGpusForNonbonded(nonbondedTarget, userGpuTaskAssignment,
                                                                 emulateGpuNonbonded,
                                                                 canUseGpuForNonbonded,
-                                                                usingVerletScheme,
                                                                 gpuAccelerationOfNonbondedIsUseful(mdlog, inputrec, !GMX_THREAD_MPI),
                                                                 gpusWereDetected);
         useGpuForPme = decideWhetherToUseGpusForPme(useGpuForNonbonded, pmeTarget, userGpuTaskAssignment,
@@ -736,7 +733,7 @@ int Mdrunner::mdrunner()
                                                     gpusWereDetected);
         auto canUseGpuForBonded = buildSupportsGpuBondeds(nullptr) && inputSupportsGpuBondeds(*inputrec, mtop, nullptr);
         useGpuForBonded =
-            decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme, usingVerletScheme,
+            decideWhetherToUseGpusForBonded(useGpuForNonbonded, useGpuForPme,
                                             bondedTarget, canUseGpuForBonded,
                                             EVDW_PME(inputrec->vdwtype),
                                             EEL_PME_EWALD(inputrec->coulombtype),
@@ -937,11 +934,8 @@ int Mdrunner::mdrunner()
         gmx_fatal(FARGS, "This group-scheme .tpr file can no longer be run by mdrun. Please update to the Verlet scheme, or use an earlier version of GROMACS if necessary.");
     }
     /* Update rlist and nstlist. */
-    if (inputrec->cutoff_scheme == ecutsVERLET)
-    {
-        prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
-                              useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
-    }
+    prepare_verlet_scheme(fplog, cr, inputrec, nstlist_cmdline, &mtop, box,
+                          useGpuForNonbonded || (emulateGpuNonbonded == EmulateGpuNonbonded::Yes), *hwinfo->cpuInfo);
 
     LocalAtomSetManager atomSets;
 
@@ -1011,17 +1005,16 @@ int Mdrunner::mdrunner()
                           physicalNodeComm.size_,
                           hw_opt.nthreads_omp,
                           hw_opt.nthreads_omp_pme,
-                          !thisRankHasDuty(cr, DUTY_PP),
-                          inputrec->cutoff_scheme == ecutsVERLET);
+                          !thisRankHasDuty(cr, DUTY_PP));
 
-    // Enable FP exception detection for the Verlet scheme, but not in
+    // Enable FP exception detection, but not in
     // Release mode and not for compilers with known buggy FP
     // exception support (clang with any optimization) or suspected
     // buggy FP exception support (gcc 7.* with optimization).
 #if !defined NDEBUG && \
     !((defined __clang__ || (defined(__GNUC__) && !defined(__ICC) && __GNUC__ == 7)) \
     && defined __OPTIMIZE__)
-    const bool bEnableFPE = inputrec->cutoff_scheme == ecutsVERLET;
+    const bool bEnableFPE = true;
 #else
     const bool bEnableFPE = false;
 #endif
index fbb3683047f8f068016547e577522bff90fef390..d29cdaadc5857759767c0a02598e8c5f73a55fba 100644 (file)
@@ -106,14 +106,12 @@ decideWhetherToUseGpusForNonbondedWithThreadMpi(const TaskTarget          nonbon
                                                 const std::vector<int>   &userGpuTaskAssignment,
                                                 const EmulateGpuNonbonded emulateGpuNonbonded,
                                                 const bool                buildSupportsNonbondedOnGpu,
-                                                const bool                usingVerletScheme,
                                                 const bool                nonbondedOnGpuIsUseful,
                                                 const int                 numRanksPerSimulation)
 {
     // First, exclude all cases where we can't run NB on GPUs.
     if (nonbondedTarget == TaskTarget::Cpu ||
         emulateGpuNonbonded == EmulateGpuNonbonded::Yes ||
-        !usingVerletScheme ||
         !nonbondedOnGpuIsUseful ||
         !buildSupportsNonbondedOnGpu)
     {
@@ -242,7 +240,6 @@ bool decideWhetherToUseGpusForNonbonded(const TaskTarget           nonbondedTarg
                                         const std::vector<int>    &userGpuTaskAssignment,
                                         const EmulateGpuNonbonded  emulateGpuNonbonded,
                                         const bool                 buildSupportsNonbondedOnGpu,
-                                        const bool                 usingVerletScheme,
                                         const bool                 nonbondedOnGpuIsUseful,
                                         const bool                 gpusWereDetected)
 {
@@ -287,18 +284,6 @@ bool decideWhetherToUseGpusForNonbonded(const TaskTarget           nonbondedTarg
         return false;
     }
 
-    if (!usingVerletScheme)
-    {
-        if (nonbondedTarget == TaskTarget::Gpu)
-        {
-            GMX_THROW(InconsistentInputError
-                          ("Nonbonded interactions on the GPU were required, which requires using "
-                          "the Verlet scheme. Either use the Verlet scheme, or do not require using GPUs."));
-        }
-
-        return false;
-    }
-
     if (!nonbondedOnGpuIsUseful)
     {
         if (nonbondedTarget == TaskTarget::Gpu)
@@ -446,7 +431,6 @@ bool decideWhetherToUseGpusForPme(const bool              useGpuForNonbonded,
 
 bool decideWhetherToUseGpusForBonded(const bool       useGpuForNonbonded,
                                      const bool       useGpuForPme,
-                                     const bool       usingVerletScheme,
                                      const TaskTarget bondedTarget,
                                      const bool       canUseGpuForBonded,
                                      const bool       usingLJPme,
@@ -459,18 +443,6 @@ bool decideWhetherToUseGpusForBonded(const bool       useGpuForNonbonded,
         return false;
     }
 
-    if (!usingVerletScheme)
-    {
-        if (bondedTarget == TaskTarget::Gpu)
-        {
-            GMX_THROW(InconsistentInputError
-                          ("Bonded interactions on the GPU were required, which requires using "
-                          "the Verlet scheme. Either use the Verlet scheme, or do not require using GPUs."));
-        }
-
-        return false;
-    }
-
     if (!canUseGpuForBonded)
     {
         if (bondedTarget == TaskTarget::Gpu)
index df0a50d2317eb66a9db68ff7324435e6fe72e20f..818841b722fc82b5777276bd651c4fa1c776c2fe 100644 (file)
@@ -82,7 +82,6 @@ enum class EmulateGpuNonbonded : bool
  * \param[in]  userGpuTaskAssignment       The user-specified assignment of GPU tasks to device IDs.
  * \param[in]  emulateGpuNonbonded         Whether we will emulate GPU calculation of nonbonded interactions.
  * \param[in]  buildSupportsNonbondedOnGpu Whether GROMACS was built with GPU support.
- * \param[in]  usingVerletScheme           Whether the nonbondeds are using the Verlet scheme.
  * \param[in]  nonbondedOnGpuIsUseful    Whether computing nonbonded interactions on a GPU is useful for this calculation.
  * \param[in]  numRanksPerSimulation     The number of ranks in each simulation.
  *
@@ -95,7 +94,6 @@ bool decideWhetherToUseGpusForNonbondedWithThreadMpi(TaskTarget              non
                                                      const std::vector<int> &userGpuTaskAssignment,
                                                      EmulateGpuNonbonded     emulateGpuNonbonded,
                                                      bool                    buildSupportsNonbondedOnGpu,
-                                                     bool                    usingVerletScheme,
                                                      bool                    nonbondedOnGpuIsUseful,
                                                      int                     numRanksPerSimulation);
 
@@ -151,7 +149,6 @@ bool decideWhetherToUseGpusForPmeWithThreadMpi(bool                    useGpuFor
  * \param[in]  userGpuTaskAssignment       The user-specified assignment of GPU tasks to device IDs.
  * \param[in]  emulateGpuNonbonded         Whether we will emulate GPU calculation of nonbonded interactions.
  * \param[in]  buildSupportsNonbondedOnGpu Whether GROMACS was build with GPU support.
- * \param[in]  usingVerletScheme           Whether the nonbondeds are using the Verlet scheme.
  * \param[in]  nonbondedOnGpuIsUseful      Whether computing nonbonded interactions on a GPU is useful for this calculation.
  * \param[in]  gpusWereDetected            Whether compatible GPUs were detected on any node.
  *
@@ -163,7 +160,6 @@ bool decideWhetherToUseGpusForNonbonded(TaskTarget              nonbondedTarget,
                                         const std::vector<int> &userGpuTaskAssignment,
                                         EmulateGpuNonbonded     emulateGpuNonbonded,
                                         bool                    buildSupportsNonbondedOnGpu,
-                                        bool                    usingVerletScheme,
                                         bool                    nonbondedOnGpuIsUseful,
                                         bool                    gpusWereDetected);
 
@@ -211,7 +207,6 @@ bool decideWhetherToUseGpusForPme(bool                    useGpuForNonbonded,
  *
  * \param[in]  useGpuForNonbonded        Whether GPUs will be used for nonbonded interactions.
  * \param[in]  useGpuForPme              Whether GPUs will be used for PME interactions.
- * \param[in]  usingVerletScheme         Whether the nonbondeds are using the Verlet scheme.
  * \param[in]  bondedTarget              The user's choice for mdrun -bonded for where to assign tasks.
  * \param[in]  canUseGpuForBonded        Whether the bonded interactions can run on a GPU
  * \param[in]  usingLJPme                Whether Vdw interactions use LJ-PME.
@@ -225,7 +220,6 @@ bool decideWhetherToUseGpusForPme(bool                    useGpuForNonbonded,
  *             InconsistentInputError  If the user requirements are inconsistent. */
 bool decideWhetherToUseGpusForBonded(bool       useGpuForNonbonded,
                                      bool       useGpuForPme,
-                                     bool       usingVerletScheme,
                                      TaskTarget bondedTarget,
                                      bool       canUseGpuForBonded,
                                      bool       usingLJPme,
index 379cf67ad169fd3d393376f080ecd5de57d472b9..90a40414c783387b29de85a576ecd0cae4906175 100644 (file)
@@ -426,12 +426,6 @@ int get_nthreads_mpi(const gmx_hw_info_t    *hwinfo,
      * the group scheme, or is a rerun with energy groups. */
     ngpu = (nonbondedOnGpu ? gmx::ssize(gpuIdsToUse) : 0);
 
-    if (inputrec->cutoff_scheme == ecutsGROUP)
-    {
-        /* We checked this before, but it doesn't hurt to do it once more */
-        GMX_RELEASE_ASSERT(hw_opt->nthreads_omp == 1, "The group scheme only supports one OpenMP thread per rank");
-    }
-
     nrank =
         get_tmpi_omp_thread_division(hwinfo, *hw_opt, nthreads_tot_max, ngpu);
 
@@ -658,10 +652,10 @@ static void print_hw_opt(FILE *fp, const gmx_hw_opt_t *hw_opt)
             hw_opt->userGpuTaskAssignment.c_str());
 }
 
-void check_and_update_hw_opt_1(const gmx::MDLogger &mdlog,
-                               gmx_hw_opt_t        *hw_opt,
-                               const bool           isSimulationMasterRank,
-                               int                  nPmeRanks)
+void checkAndUpdateHardwareOptions(const gmx::MDLogger &mdlog,
+                                   gmx_hw_opt_t        *hw_opt,
+                                   const bool           isSimulationMasterRank,
+                                   const int            nPmeRanks)
 {
     /* Currently hw_opt only contains default settings or settings supplied
      * by the user on the command line.
index 0acdca32a5ac71e8c834a51763cfa15f76fbb4fd..6d8f1c235465098936741402880385adbe0d473d 100644 (file)
@@ -97,11 +97,18 @@ void check_resource_division_efficiency(const gmx_hw_info_t *hwinfo,
                                         t_commrec           *cr,
                                         const gmx::MDLogger &mdlog);
 
-/*! \brief Checks we can do when we don't (yet) know the cut-off scheme */
-void check_and_update_hw_opt_1(const gmx::MDLogger &mdlog,
-                               gmx_hw_opt_t        *hw_opt,
-                               bool                 isSimulationMasterRank,
-                               int                  nPmeRanks);
+/*! \brief Checks what our hardware options are based on how Gromacs was compiled
+ *  and user-set options
+ *
+ *  \param[in]      mdlog                     Logger
+ *  \param[in, out] hw_opt                    Hardware-related and threading options
+ *  \param[in]      isSimulationMasterRank
+ *  \param[in]      nPmeRanks                 Number of PME ranks
+ *  */
+void checkAndUpdateHardwareOptions(const gmx::MDLogger &mdlog,
+                                   gmx_hw_opt_t        *hw_opt,
+                                   bool                 isSimulationMasterRank,
+                                   int                  nPmeRanks);
 
 /*! \brief Check, and if necessary update, the number of OpenMP threads requested
  *