Pass ArrayRefs to dispatchNonbondedKernel
authorJoe Jordan <ejjordan12@gmail.com>
Wed, 5 May 2021 07:35:26 +0000 (07:35 +0000)
committerJoe Jordan <ejjordan12@gmail.com>
Wed, 5 May 2021 07:35:26 +0000 (07:35 +0000)
api/nblib/gmxcalculator.cpp
src/gromacs/mdlib/sim_util.cpp
src/gromacs/nbnxm/benchmark/bench_setup.cpp
src/gromacs/nbnxm/kerneldispatch.cpp
src/gromacs/nbnxm/nbnxm.h

index 12ab681d5ee61c29827ba50890026b349b83e1c8..751ed51518d63043fb7b3b06b145b306bd9bf4ca 100644 (file)
@@ -73,13 +73,16 @@ void GmxForceCalculator::compute(gmx::ArrayRef<const gmx::RVec> coordinateInput,
     // update the coordinates in the backend
     nbv_->convertCoordinates(gmx::AtomLocality::Local, coordinateInput);
 
-    nbv_->dispatchNonbondedKernel(gmx::InteractionLocality::Local,
-                                  *interactionConst_,
-                                  *stepWork_,
-                                  enbvClearFYes,
-                                  *forcerec_,
-                                  enerd_.get(),
-                                  nrnb_.get());
+    nbv_->dispatchNonbondedKernel(
+            gmx::InteractionLocality::Local,
+            *interactionConst_,
+            *stepWork_,
+            enbvClearFYes,
+            forcerec_->shift_vec,
+            enerd_->grpp.energyGroupPairTerms[forcerec_->haveBuckingham ? NonBondedEnergyTerms::BuckinghamSR
+                                                                        : NonBondedEnergyTerms::LJSR],
+            enerd_->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR],
+            nrnb_.get());
 
     nbv_->atomdata_add_nbat_f_to_f(gmx::AtomLocality::All, forceOutput);
 }
index 540bad7f0046adb39e9cef08b7f7e13873831eea..342c08bb625d865c9eb287baee7ec84e4279ae4c 100644 (file)
@@ -444,7 +444,16 @@ static void do_nb_verlet(t_forcerec*                fr,
         }
     }
 
-    nbv->dispatchNonbondedKernel(ilocality, *ic, stepWork, clearF, *fr, enerd, nrnb);
+    nbv->dispatchNonbondedKernel(
+            ilocality,
+            *ic,
+            stepWork,
+            clearF,
+            fr->shift_vec,
+            enerd->grpp.energyGroupPairTerms[fr->haveBuckingham ? NonBondedEnergyTerms::BuckinghamSR
+                                                                : NonBondedEnergyTerms::LJSR],
+            enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR],
+            nrnb);
 }
 
 static inline void clearRVecs(ArrayRef<RVec> v, const bool useOpenmpThreading)
index 642818194676e3233a85c3d5f4400b7ea6a4b121..5ac216df225f38e106e55bcdfb2c0ef5943e8a06 100644 (file)
@@ -49,6 +49,7 @@
 
 #include "gromacs/gmxlib/nrnb.h"
 #include "gromacs/math/units.h"
+#include "gromacs/math/vectypes.h"
 #include "gromacs/mdlib/dispersioncorrection.h"
 #include "gromacs/mdlib/force_flags.h"
 #include "gromacs/mdlib/forcerec.h"
@@ -358,7 +359,15 @@ static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
     for (int iter = 0; iter < options.numPreIterations; iter++)
     {
         nbv->dispatchNonbondedKernel(
-                gmx::InteractionLocality::Local, ic, stepWork, enbvClearFYes, system.forceRec, &enerd, &nrnb);
+                gmx::InteractionLocality::Local,
+                ic,
+                stepWork,
+                enbvClearFYes,
+                system.forceRec.shift_vec,
+                enerd.grpp.energyGroupPairTerms[system.forceRec.haveBuckingham ? NonBondedEnergyTerms::BuckinghamSR
+                                                                               : NonBondedEnergyTerms::LJSR],
+                enerd.grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR],
+                &nrnb);
     }
 
     const int numIterations = (doWarmup ? options.numWarmupIterations : options.numIterations);
@@ -369,7 +378,15 @@ static void setupAndRunInstance(const gmx::BenchmarkSystem& system,
     {
         // Run the kernel without force clearing
         nbv->dispatchNonbondedKernel(
-                gmx::InteractionLocality::Local, ic, stepWork, enbvClearFNo, system.forceRec, &enerd, &nrnb);
+                gmx::InteractionLocality::Local,
+                ic,
+                stepWork,
+                enbvClearFNo,
+                system.forceRec.shift_vec,
+                enerd.grpp.energyGroupPairTerms[system.forceRec.haveBuckingham ? NonBondedEnergyTerms::BuckinghamSR
+                                                                               : NonBondedEnergyTerms::LJSR],
+                enerd.grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR],
+                &nrnb);
     }
     cycles = gmx_cycles_read() - cycles;
     if (!doWarmup)
index 52667776a0b05993c801c9408d282b33ce74b236..f84ddf3a56f68b3b4a4efe9afed979c7eefb5e5b 100644 (file)
@@ -441,13 +441,14 @@ static void accountFlops(t_nrnb*                    nrnb,
     }
 }
 
-void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality   iLocality,
-                                                 const interaction_const_t& ic,
-                                                 const gmx::StepWorkload&   stepWork,
-                                                 int                        clearF,
-                                                 const t_forcerec&          fr,
-                                                 gmx_enerdata_t*            enerd,
-                                                 t_nrnb*                    nrnb)
+void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality       iLocality,
+                                                 const interaction_const_t&     ic,
+                                                 const gmx::StepWorkload&       stepWork,
+                                                 int                            clearF,
+                                                 gmx::ArrayRef<const gmx::RVec> shiftvec,
+                                                 gmx::ArrayRef<real> repulsionDispersionSR,
+                                                 gmx::ArrayRef<real> CoulombSR,
+                                                 t_nrnb*             nrnb) const
 {
     const PairlistSet& pairlistSet = pairlistSets().pairlistSet(iLocality);
 
@@ -456,19 +457,16 @@ void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality   iLoc
         case Nbnxm::KernelType::Cpu4x4_PlainC:
         case Nbnxm::KernelType::Cpu4xN_Simd_4xN:
         case Nbnxm::KernelType::Cpu4xN_Simd_2xNN:
-            nbnxn_kernel_cpu(
-                    pairlistSet,
-                    kernelSetup(),
-                    nbat.get(),
-                    ic,
-                    fr.shift_vec,
-                    stepWork,
-                    clearF,
-                    enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
-                    fr.haveBuckingham
-                            ? enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR].data()
-                            : enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data(),
-                    wcycle_);
+            nbnxn_kernel_cpu(pairlistSet,
+                             kernelSetup(),
+                             nbat.get(),
+                             ic,
+                             shiftvec,
+                             stepWork,
+                             clearF,
+                             CoulombSR.data(),
+                             repulsionDispersionSR.data(),
+                             wcycle_);
             break;
 
         case Nbnxm::KernelType::Gpu8x8x8:
@@ -476,19 +474,16 @@ void nonbonded_verlet_t::dispatchNonbondedKernel(gmx::InteractionLocality   iLoc
             break;
 
         case Nbnxm::KernelType::Cpu8x8x8_PlainC:
-            nbnxn_kernel_gpu_ref(
-                    pairlistSet.gpuList(),
-                    nbat.get(),
-                    &ic,
-                    fr.shift_vec,
-                    stepWork,
-                    clearF,
-                    nbat->out[0].f,
-                    nbat->out[0].fshift.data(),
-                    enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::CoulombSR].data(),
-                    fr.haveBuckingham
-                            ? enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::BuckinghamSR].data()
-                            : enerd->grpp.energyGroupPairTerms[NonBondedEnergyTerms::LJSR].data());
+            nbnxn_kernel_gpu_ref(pairlistSet.gpuList(),
+                                 nbat.get(),
+                                 &ic,
+                                 shiftvec,
+                                 stepWork,
+                                 clearF,
+                                 nbat->out[0].f,
+                                 nbat->out[0].fshift.data(),
+                                 CoulombSR.data(),
+                                 repulsionDispersionSR.data());
             break;
 
         default: GMX_RELEASE_ASSERT(false, "Invalid nonbonded kernel type passed!");
index b44fcaa8087b1b904f79104ebd0eaa0b38f8642f..6d5cdccd694427c814b57850cca351b399b3a3a3 100644 (file)
@@ -358,13 +358,14 @@ public:
     void dispatchPruneKernelGpu(int64_t step);
 
     //! \brief Executes the non-bonded kernel of the GPU or launches it on the GPU
-    void dispatchNonbondedKernel(gmx::InteractionLocality   iLocality,
-                                 const interaction_const_t& ic,
-                                 const gmx::StepWorkload&   stepWork,
-                                 int                        clearF,
-                                 const t_forcerec&          fr,
-                                 gmx_enerdata_t*            enerd,
-                                 t_nrnb*                    nrnb);
+    void dispatchNonbondedKernel(gmx::InteractionLocality       iLocality,
+                                 const interaction_const_t&     ic,
+                                 const gmx::StepWorkload&       stepWork,
+                                 int                            clearF,
+                                 gmx::ArrayRef<const gmx::RVec> shiftvec,
+                                 gmx::ArrayRef<real>            repulsionDispersionSR,
+                                 gmx::ArrayRef<real>            CoulombSR,
+                                 t_nrnb*                        nrnb) const;
 
     //! Executes the non-bonded free-energy kernel, always runs on the CPU
     void dispatchFreeEnergyKernel(gmx::InteractionLocality       iLocality,