Use unique_ptr<ListedForcesGpu> in forcerec
authorJoe Jordan <ejjordan12@gmail.com>
Wed, 16 Jun 2021 09:42:22 +0000 (09:42 +0000)
committerJoe Jordan <ejjordan12@gmail.com>
Wed, 16 Jun 2021 09:42:22 +0000 (09:42 +0000)
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdrun/runner.cpp
src/gromacs/mdtypes/forcerec.h

index 8795ec6188643a4d001b362b55928739c0831ed7..5b3f72dd871f7a64a4efb47cf984141695ce9544 100644 (file)
@@ -1470,7 +1470,7 @@ void do_force(FILE*                               fplog,
         /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
         nbv->constructPairlist(InteractionLocality::Local, top->excls, step, nrnb);
 
-        nbv->setupGpuShortRangeWork(fr->listedForcesGpu, InteractionLocality::Local);
+        nbv->setupGpuShortRangeWork(fr->listedForcesGpu.get(), InteractionLocality::Local);
 
         wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchLocal);
         wallcycle_stop(wcycle, WallCycleCounter::NS);
@@ -1559,7 +1559,7 @@ void do_force(FILE*                               fplog,
             /* Note that with a GPU the launch overhead of the list transfer is not timed separately */
             nbv->constructPairlist(InteractionLocality::NonLocal, top->excls, step, nrnb);
 
-            nbv->setupGpuShortRangeWork(fr->listedForcesGpu, InteractionLocality::NonLocal);
+            nbv->setupGpuShortRangeWork(fr->listedForcesGpu.get(), InteractionLocality::NonLocal);
             wallcycle_sub_stop(wcycle, WallCycleSubCounter::NBSSearchNonLocal);
             wallcycle_stop(wcycle, WallCycleCounter::NS);
             // TODO refactor this GPU halo exchange re-initialisation
@@ -2285,7 +2285,8 @@ void do_force(FILE*                               fplog,
         }
     }
 
-    launchGpuEndOfStepTasks(nbv, fr->listedForcesGpu, fr->pmedata, enerd, *runScheduleWork, step, wcycle);
+    launchGpuEndOfStepTasks(
+            nbv, fr->listedForcesGpu.get(), fr->pmedata, enerd, *runScheduleWork, step, wcycle);
 
     if (DOMAINDECOMP(cr))
     {
index a1cc4729582b04c96ce03259a9b69e80d90ec91c..45c569630c2b778683cfc54e3808fbd41f8bd8a0 100644 (file)
@@ -1611,7 +1611,6 @@ int Mdrunner::mdrunner()
     const bool               thisRankHasPmeGpuTask = gpuTaskAssignments.thisRankHasPmeGpuTask();
     std::unique_ptr<MDAtoms> mdAtoms;
     std::unique_ptr<VirtualSitesHandler> vsite;
-    std::unique_ptr<ListedForcesGpu>     listedForcesGpu;
 
     t_nrnb nrnb;
     if (thisRankHasDuty(cr, DUTY_PP))
@@ -1677,13 +1676,12 @@ int Mdrunner::mdrunner()
             GMX_RELEASE_ASSERT(deviceStreamManager != nullptr,
                                "GPU device stream manager should be valid in order to use GPU "
                                "version of bonded forces.");
-            listedForcesGpu = std::make_unique<ListedForcesGpu>(
+            fr->listedForcesGpu = std::make_unique<ListedForcesGpu>(
                     mtop.ffparams,
                     fr->ic->epsfac * fr->fudgeQQ,
                     deviceStreamManager->context(),
                     deviceStreamManager->bondedStream(havePPDomainDecomposition(cr)),
                     wcycle.get());
-            fr->listedForcesGpu = listedForcesGpu.get();
         }
 
         /* Initialize the mdAtoms structure.
@@ -2073,8 +2071,7 @@ int Mdrunner::mdrunner()
     mdAtoms.reset(nullptr);
     globalState.reset(nullptr);
     mdModules_.reset(nullptr); // destruct force providers here as they might also use the GPU
-    listedForcesGpu.reset(nullptr);
-    fr.reset(nullptr); // destruct forcerec before gpu
+    fr.reset(nullptr);         // destruct forcerec before gpu
     // TODO convert to C++ so we can get rid of these frees
     sfree(disresdata);
 
index 8abdadfa9f068c4718962e64d6ca1b8c8d01233a..d4bcb4c350add584f95969cebec9eca4bf1c791f 100644 (file)
@@ -257,8 +257,8 @@ struct t_forcerec
     // The listed forces calculation data, 1 entry or multiple entries with multiple time stepping
     std::vector<ListedForces> listedForces;
 
-    /* TODO: Replace the pointer by an object once we got rid of C */
-    gmx::ListedForcesGpu* listedForcesGpu = nullptr;
+    // The listed forces calculation data for GPU
+    std::unique_ptr<gmx::ListedForcesGpu> listedForcesGpu;
 
     /* Ewald correction thread local virial and energy data */
     int                              nthread_ewc = 0;