Ensure coordinates are copied for dipole moment calculation
authorAndrey Alekseenko <al42and@gmail.com>
Thu, 1 Apr 2021 13:40:30 +0000 (15:40 +0200)
committerArtem Zhmurov <zhmurov@gmail.com>
Fri, 2 Apr 2021 13:53:50 +0000 (13:53 +0000)
With GPU update, we should wait for the coordinates to be copied to host
before using them for calculating the dipole moment.

This did not seem to cause any issues with CUDA but was causing a unit
test failure with SYCL (with !1329):

    SYCL_BE=PI_OPENCL GMX_USE_GPU_BUFFER_OPS=1 GMX_FORCE_UPDATE_DEFAULT_GPU=1 ./bin/mdrun-test --gtest_filter=EwaldSurfaceTerm/EwaldSurfaceTermTest.WithinTolerances/0

Refs #3930, #3932

src/gromacs/mdlib/sim_util.cpp

index d72fed314630230e086b2f7094fce9f190dcee6e..9f6711386b86fdb16f1bf263aff148ba0b77a65d 100644 (file)
@@ -1291,7 +1291,7 @@ void do_force(FILE*                               fplog,
     bool gmx_used_in_debug haveCopiedXFromGpu = false;
     if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch
         && (runScheduleWork->domainWork.haveCpuLocalForceWork || stepWork.computeVirial
-            || haveHostPmePpComms || haveHostHaloExchangeComms))
+            || haveHostPmePpComms || haveHostHaloExchangeComms || simulationWork.computeMuTot))
     {
         stateGpu->copyCoordinatesFromGpu(x.unpaddedArrayRef(), AtomLocality::Local);
         haveCopiedXFromGpu = true;
@@ -1641,6 +1641,13 @@ void do_force(FILE*                               fplog,
     {
         const int start = 0;
 
+        if (simulationWork.useGpuUpdate && !stepWork.doNeighborSearch)
+        {
+            GMX_ASSERT(haveCopiedXFromGpu,
+                       "a wait should only be triggered if copy has been scheduled");
+            stateGpu->waitCoordinatesReadyOnHost(AtomLocality::Local);
+        }
+
         /* Calculate total (local) dipole moment in a temporary common array.
          * This makes it possible to sum them over nodes faster.
          */