Move the copy of the velocities for the output closer to the output
authorArtem Zhmurov <zhmurov@gmail.com>
Tue, 29 Oct 2019 20:18:08 +0000 (21:18 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 20 Nov 2019 18:42:25 +0000 (19:42 +0100)
The D2H copy of velocities is split into two calls, one of which is
moved closer to the output, where it is consumed. This is in preparation
to creation a unified method to copy all data for the outuput when
it is needed.

Change-Id: I4b81aa95e2c0763bb3ebb4a7d8b18611765c54de

src/gromacs/mdrun/md.cpp

index 4756c64ffe78f149d19f32f018e8f045b02d8ac3..b0b006b8bc785224c54251d7605f2184d1cb0e61 100644 (file)
@@ -799,25 +799,16 @@ void gmx::LegacySimulator::do_md()
         do_verbose = mdrunOptions.verbose
                      && (step % mdrunOptions.verboseStepPrintInterval == 0 || bFirstStep || bLastStep);
 
-        if (useGpuForUpdate && !bFirstStep)
+        if (useGpuForUpdate && !bFirstStep && bNS)
         {
-            // Copy velocities from the GPU when needed:
-            // - On search steps to keep copy on host (device buffers are reinitialized).
-            // - When needed for the output.
-            if (bNS || do_per_step(step, ir->nstvout))
-            {
-                stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
-                stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
-            }
-
+            // Copy velocities from the GPU on search steps to keep a copy on host (device buffers are reinitialized).
+            stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
+            stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
             // Copy coordinate from the GPU when needed at the search step.
             // NOTE: The cases when coordinates needed on CPU for force evaluation are handled in sim_utils.
             // NOTE: If the coordinates are to be written into output file they are also copied separately before the output.
-            if (bNS)
-            {
-                stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
-                stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
-            }
+            stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
+            stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
         }
 
         if (bNS && !(bFirstStep && ir->bContinuation))
@@ -1124,7 +1115,13 @@ void gmx::LegacySimulator::do_md()
             stateGpu->copyCoordinatesFromGpu(state->x, AtomLocality::Local);
             stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
         }
-
+        // Copy velocities if needed for the output.
+        // NOTE: Copy on the search steps is done at the beginning of the step.
+        if (useGpuForUpdate && !bNS && do_per_step(step, ir->nstvout))
+        {
+            stateGpu->copyVelocitiesFromGpu(state->v, AtomLocality::Local);
+            stateGpu->waitVelocitiesReadyOnHost(AtomLocality::Local);
+        }
         // Copy forces for the output if the forces were reduced on the GPU (not the case on virial steps)
         // and update is offloaded hence forces are kept on the GPU for update and have not been
         // already transferred in do_force().