Make forcerec const in dispatchNonbondedKernel()
authorBerk Hess <hess@kth.se>
Tue, 18 Jun 2019 08:28:19 +0000 (10:28 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 18 Jun 2019 19:51:57 +0000 (21:51 +0200)
Changed forcrecin Nbnxm::dispatchNonbondedKernel() from a pointer
to a const reference. The GPU emulation kernel now writes the shift
forces to nbnxn_atomdata_t.

Change-Id: I170165529ba447db9b5d1824eea10675ae1bef7d

src/gromacs/mdlib/sim_util.cpp
src/gromacs/nbnxm/kerneldispatch.cpp
src/gromacs/nbnxm/nbnxm.h

index d53d02d6945ac5ed76aae10e1bc09bfcba50e9e0..84736d802d2d5abd9571b009fe75db5fe428237b 100644 (file)
@@ -350,7 +350,7 @@ static void do_nb_verlet(t_forcerec                       *fr,
         wallcycle_sub_start(wcycle, ewcsNONBONDED);
     }
 
-    nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, fr, enerd, nrnb);
+    nbv->dispatchNonbondedKernel(ilocality, *ic, flags, clearF, *fr, enerd, nrnb);
 
     if (!nbv->useGpu())
     {
@@ -1354,6 +1354,12 @@ void do_force(FILE                                     *fplog,
 
             nbv->atomdata_add_nbat_f_to_f(Nbnxm::AtomLocality::NonLocal,
                                           forceOut.f, wcycle);
+
+            if (fr->nbv->emulateGpu() && (flags & GMX_FORCE_VIRIAL))
+            {
+                nbnxn_atomdata_add_nbat_fshift_to_fshift(nbv->nbat.get(),
+                                                         fr->fshift);
+            }
         }
     }
 
index 721ef867f174ceb9f9a794b2052f7b93c10d9b8e..eb8905525c8fd9f0113b08d9586b96160f0bb644 100644 (file)
@@ -454,7 +454,7 @@ nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality
                                             const interaction_const_t &ic,
                                             int                        forceFlags,
                                             int                        clearF,
-                                            t_forcerec                *fr,
+                                            const t_forcerec          &fr,
                                             gmx_enerdata_t            *enerd,
                                             t_nrnb                    *nrnb)
 {
@@ -469,11 +469,11 @@ nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality
                              kernelSetup(),
                              nbat.get(),
                              ic,
-                             fr->shift_vec,
+                             fr.shift_vec,
                              forceFlags,
                              clearF,
                              enerd->grpp.ener[egCOULSR].data(),
-                             fr->bBHAM ?
+                             fr.bBHAM ?
                              enerd->grpp.ener[egBHAMSR].data() :
                              enerd->grpp.ener[egLJSR].data());
             break;
@@ -485,13 +485,13 @@ nonbonded_verlet_t::dispatchNonbondedKernel(Nbnxm::InteractionLocality iLocality
         case Nbnxm::KernelType::Cpu8x8x8_PlainC:
             nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
                                  nbat.get(), &ic,
-                                 fr->shift_vec,
+                                 fr.shift_vec,
                                  forceFlags,
                                  clearF,
                                  nbat->out[0].f,
-                                 fr->fshift[0],
+                                 nbat->out[0].fshift.data(),
                                  enerd->grpp.ener[egCOULSR].data(),
-                                 fr->bBHAM ?
+                                 fr.bBHAM ?
                                  enerd->grpp.ener[egBHAMSR].data() :
                                  enerd->grpp.ener[egLJSR].data());
             break;
index ef5bda9b8789cca668784e3d7cea62bb42c5c7b5..32ae108d1bb627c48b1c75ffe9a407d8cde55330 100644 (file)
@@ -276,7 +276,7 @@ struct nonbonded_verlet_t
                                      const interaction_const_t  &ic,
                                      int                         forceFlags,
                                      int                         clearF,
-                                     t_forcerec                 *fr,
+                                     const t_forcerec           &fr,
                                      gmx_enerdata_t             *enerd,
                                      t_nrnb                     *nrnb);