Properly check for frozen atoms when disabling GPU update
authorAndrey Alekseenko <al42and@gmail.com>
Thu, 11 Mar 2021 05:33:39 +0000 (05:33 +0000)
committerArtem Zhmurov <zhmurov@gmail.com>
Thu, 11 Mar 2021 05:33:39 +0000 (05:33 +0000)
Follow-up of !1241 (9550c3e8564068c86a0f34baee0f578e8ce6d5a9)

Initially, the check was overly relaxed, and was triggered even when
there were no frozen atoms.

Thanks to @alangray3 for noticing!

src/gromacs/mdlib/mdatoms.cpp
src/gromacs/mdrun/runner.cpp
src/gromacs/mdtypes/inputrec.cpp
src/gromacs/mdtypes/inputrec.h
src/gromacs/modularsimulator/modularsimulator.cpp

index 66e412cee39e361828fc4294aeb43d0df88cd961..1fc7f919e8cbd86eaceacbe6c5634039e4677e37 100644 (file)
@@ -288,8 +288,7 @@ void atoms2md(const gmx_mtop_t*  mtop,
             /* We always copy cTC with domain decomposition */
         }
         srenew(md->cENER, md->nalloc);
-        if (opts->nFreeze
-            && (opts->ngfrz > 1 || opts->nFreeze[0][XX] || opts->nFreeze[0][YY] || opts->nFreeze[0][ZZ]))
+        if (inputrecFrozenAtoms(ir))
         {
             srenew(md->cFREEZE, md->nalloc);
         }
index 05a7239325948736ff713491b597bd6675ad5eea..8202effc13b78408d5340807e8a4d20905d3f252 100644 (file)
@@ -1361,6 +1361,7 @@ int Mdrunner::mdrunner()
     try
     {
         const bool useUpdateGroups = cr->dd ? ddUsesUpdateGroups(*cr->dd) : false;
+        const bool haveFrozenAtoms = inputrecFrozenAtoms(inputrec.get());
 
         useGpuForUpdate = decideWhetherToUseGpuForUpdate(useDomainDecomposition,
                                                          useUpdateGroups,
@@ -1374,7 +1375,7 @@ int Mdrunner::mdrunner()
                                                          doEssentialDynamics,
                                                          gmx_mtop_ftype_count(mtop, F_ORIRES) > 0,
                                                          replExParams.exchangeInterval > 0,
-                                                         inputrec->opts.nFreeze != nullptr,
+                                                         haveFrozenAtoms,
                                                          doRerun,
                                                          devFlags,
                                                          mdlog);
index ab292ec5231c3c768cae3b427614cdc12dec386e..0300c4e5d891c58625de3ed393a355a960093655 100644 (file)
@@ -1584,6 +1584,13 @@ bool inputrecPbcXY2Walls(const t_inputrec* ir)
     return (ir->pbcType == PbcType::XY && ir->nwall == 2);
 }
 
+bool inputrecFrozenAtoms(const t_inputrec* ir)
+{
+    return ((ir->opts.nFreeze != nullptr)
+            && (ir->opts.ngfrz > 1 || ir->opts.nFreeze[0][XX] != 0 || ir->opts.nFreeze[0][YY] != 0
+                || ir->opts.nFreeze[0][ZZ] != 0));
+}
+
 bool integratorHasConservedEnergyQuantity(const t_inputrec* ir)
 {
     if (!EI_MD(ir->eI))
index 7b50eb0b9b07efbc03c3ca6453c99743a9e8fe98..5f97a2e4edbc198be0f9085298ed73439e833e0a 100644 (file)
@@ -643,6 +643,9 @@ gmx_bool inputrecNphTrotter(const t_inputrec* ir);
 /*! \brief Return true if the simulation is 2D periodic with two walls. */
 bool inputrecPbcXY2Walls(const t_inputrec* ir);
 
+//! \brief Return true if the simulation has frozen atoms (non-trivial freeze groups).
+bool inputrecFrozenAtoms(const t_inputrec* ir);
+
 /*! \brief Returns true for MD integator with T and/or P-coupling that supports
  * calculating a conserved energy quantity.
  *
index fbb962be3dd49771f9fee1e43147ca0ce0b7c126..b285fc054c5819f8825e5aea0f2b74cdbc248d4b 100644 (file)
@@ -262,9 +262,7 @@ bool ModularSimulator::isInputCompatible(bool                             exitOn
                                  "Acceleration is not supported by the modular simulator.");
     isInputCompatible =
             isInputCompatible
-            && conditionalAssert(inputrec->opts.ngfrz == 1 && inputrec->opts.nFreeze[0][XX] == 0
-                                         && inputrec->opts.nFreeze[0][YY] == 0
-                                         && inputrec->opts.nFreeze[0][ZZ] == 0,
+            && conditionalAssert(!inputrecFrozenAtoms(inputrec),
                                  "Freeze groups are not supported by the modular simulator.");
     isInputCompatible =
             isInputCompatible