Reorder PBC setup in do_force_lowlevel
authorSzilárd Páll <pall.szilard@gmail.com>
Tue, 12 Feb 2019 17:40:19 +0000 (18:40 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 21 Mar 2019 13:03:47 +0000 (14:03 +0100)
Eliminated redundant set_pbc() call so now there is always a single call
to the DD-optimized set_pbc_dd() called only when needed with clarified
conditionals for the different use cases (listed forces and Group
scheme).

Also improve the t_forcerec.bMolPBC documentation.

Change-Id: I8a0782bf81e05f1c7221ee5f0da42bf5e4d24ac9

src/gromacs/mdlib/force.cpp
src/gromacs/mdtypes/forcerec.h

index 3ca8b6de2925eee675f4824decc4a34bd1084543..fb8359ea45c45b43017aa7984ddd86748edd7543 100644 (file)
@@ -165,15 +165,12 @@ void do_force_lowlevel(t_forcerec                   *fr,
     int         i, j;
     int         donb_flags;
     int         pme_flags;
-    t_pbc       pbc;
     real        dvdl_dum[efptNR], dvdl_nb[efptNR];
 
 #if GMX_MPI
     double  t0 = 0.0, t1, t2, t3; /* time measurement for coarse load balancing */
 #endif
 
-    set_pbc(&pbc, fr->ePBC, box);
-
     /* reset free energy components */
     for (i = 0; i < efptNR; i++)
     {
@@ -320,19 +317,24 @@ void do_force_lowlevel(t_forcerec                   *fr,
             inc_nrnb(nrnb, eNR_SHIFTX, graph->nnodes);
         }
     }
-    /* Check whether we need to do listed interactions or correct for exclusions */
-    if (fr->bMolPBC &&
-        ((flags & GMX_FORCE_LISTED)
-         || EEL_RF(fr->ic->eeltype) || EEL_FULL(fr->ic->eeltype) || EVDW_PME(fr->ic->vdwtype)))
+
+    t_pbc      pbc;
+    const bool useRfWithGroupScheme = (fr->cutoff_scheme == ecutsGROUP) && EEL_RF(fr->ic->eeltype);
+    /* Check whether we need to take into account PBC in the following force tasks:
+     * listed interactions or when correcting for exclusions (for the Group scheme with RF) */
+    if (fr->bMolPBC)
     {
-        /* TODO There are no electrostatics methods that require this
-           transformation, when using the Verlet scheme, so update the
-           above conditional. */
-        /* Since all atoms are in the rectangular or triclinic unit-cell,
-         * only single box vector shifts (2 in x) are required.
-         */
-        set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
-                   TRUE, box);
+        // TODO: if all listed forces are offloaded to the GPU, both PBC calls below could be skipped.
+        // TODO: after rebase add && haveCpuListedForces(*fr, idef, *fcd) to the condition
+        const auto needPbcForListedForces = bool(flags & GMX_FORCE_LISTED);
+        if (needPbcForListedForces || useRfWithGroupScheme)
+        {
+            /* Since all atoms are in the rectangular or triclinic unit-cell,
+             * only single box vector shifts (2 in x) are required.
+             */
+            set_pbc_dd(&pbc, fr->ePBC, DOMAINDECOMP(cr) ? cr->dd->nc : nullptr,
+                       TRUE, box);
+        }
     }
 
     do_force_listed(wcycle, box, ir->fepvals, cr, ms,
@@ -550,7 +552,7 @@ void do_force_lowlevel(t_forcerec                   *fr,
          * With the Verlet scheme, exclusion forces are calculated
          * in the non-bonded kernel.
          */
-        if (ir->cutoff_scheme != ecutsVERLET && EEL_RF(fr->ic->eeltype))
+        if (useRfWithGroupScheme)
         {
             real dvdl_rf_excl      = 0;
             enerd->term[F_RF_EXCL] =
index e90470a62a8145907ce1b1b00a0995db9340c8a4..2d3e5498b1a9207b958791084a5f610fed7a39c7 100644 (file)
@@ -130,7 +130,9 @@ struct t_forcerec { // NOLINT (clang-analyzer-optin.performance.Padding)
 
     /* PBC stuff */
     int                         ePBC = 0;
-    //! Whether PBC must be considered for bonded interactions.
+    //! Tells whether atoms inside a molecule can be in different periodic images,
+    //  i.e. whether we need to take into account PBC when computing distances inside molecules.
+    //  This determines whether PBC must be considered for e.g. bonded interactions.
     gmx_bool                    bMolPBC     = FALSE;
     int                         rc_scaling  = 0;
     rvec                        posres_com  = { 0 };