Remove GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE legacy flag in do_force()
authorSzilárd Páll <pall.szilard@gmail.com>
Tue, 11 May 2021 10:19:42 +0000 (10:19 +0000)
committerArtem Zhmurov <zhmurov@gmail.com>
Tue, 11 May 2021 10:19:42 +0000 (10:19 +0000)
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdrun/md.cpp
src/gromacs/mdtypes/simulation_workload.h

index 2515cc148ca384d759b2dbc24b7267311c509353..a93b67770d1fe46ffadc4c9c590bc52fc3443b30 100644 (file)
@@ -959,14 +959,15 @@ static StepWorkload setupStepWorkload(const int                     legacyFlags,
     const bool computeSlowForces = (mtsLevels.empty() || step % mtsLevels[1].stepFactor == 0);
 
     StepWorkload flags;
-    flags.stateChanged        = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
-    flags.haveDynamicBox      = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
-    flags.doNeighborSearch    = ((legacyFlags & GMX_FORCE_NS) != 0);
-    flags.computeSlowForces   = computeSlowForces;
-    flags.computeVirial       = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
-    flags.computeEnergy       = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
-    flags.computeForces       = ((legacyFlags & GMX_FORCE_FORCES) != 0);
-    flags.computeListedForces = ((legacyFlags & GMX_FORCE_LISTED) != 0);
+    flags.stateChanged                  = ((legacyFlags & GMX_FORCE_STATECHANGED) != 0);
+    flags.haveDynamicBox                = ((legacyFlags & GMX_FORCE_DYNAMICBOX) != 0);
+    flags.doNeighborSearch              = ((legacyFlags & GMX_FORCE_NS) != 0);
+    flags.computeSlowForces             = computeSlowForces;
+    flags.computeVirial                 = ((legacyFlags & GMX_FORCE_VIRIAL) != 0);
+    flags.computeEnergy                 = ((legacyFlags & GMX_FORCE_ENERGY) != 0);
+    flags.computeForces                 = ((legacyFlags & GMX_FORCE_FORCES) != 0);
+    flags.useOnlyMtsCombinedForceBuffer = ((legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0);
+    flags.computeListedForces           = ((legacyFlags & GMX_FORCE_LISTED) != 0);
     flags.computeNonbondedForces =
             ((legacyFlags & GMX_FORCE_NONBONDED) != 0) && simulationWork.computeNonbonded
             && !(simulationWork.computeNonbondedAtMtsLevel1 && !computeSlowForces);
@@ -2052,8 +2053,7 @@ void do_force(FILE*                               fplog,
      * avoids an extra halo exchange (when DD is used) and post-processing step.
      */
     const bool combineMtsForcesBeforeHaloExchange =
-            (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces
-             && (legacyFlags & GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE) != 0
+            (stepWork.computeForces && fr->useMts && stepWork.computeSlowForces && stepWork.useOnlyMtsCombinedForceBuffer
              && !(stepWork.computeVirial || simulationWork.useGpuNonbonded || stepWork.haveGpuPmeOnThisRank));
     if (combineMtsForcesBeforeHaloExchange)
     {
index b82675c014d4f50a3762f7625a23483fe0640d4e..68428e84ef0690f840481d09fb36a7b3f8bcea24 100644 (file)
@@ -1127,6 +1127,7 @@ void gmx::LegacySimulator::do_md()
                        | (bCalcEner ? GMX_FORCE_ENERGY : 0) | (bDoFEP ? GMX_FORCE_DHDL : 0));
         if (fr->useMts && !do_per_step(step, ir->nstfout))
         {
+            // TODO: merge this with stepWork.useOnlyMtsCombinedForceBuffer
             force_flags |= GMX_FORCE_DO_NOT_NEED_NORMAL_FORCE;
         }
 
index 281c1fc6eb7cd4f23cbdce5b8c39f75a7d5fdafd..4b7580f3ad6365ff83e85c1ce70b4b46522accb8 100644 (file)
@@ -77,6 +77,8 @@ public:
     bool computeEnergy = false;
     //! Whether (any) forces need to be computed this step, not only energies
     bool computeForces = false;
+    //! Whether only the MTS combined force buffers are needed and not the separate normal force buffer.
+    bool useOnlyMtsCombinedForceBuffer = false;
     //! Whether nonbonded forces need to be computed this step
     bool computeNonbondedForces = false;
     //! Whether listed forces need to be computed this step