Activate GPU update support in SYCL build
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_gpu_impl.cpp
index df6420d378f1de257e550a3f5a92ce410223f805..9b8d102d2533889d8a305248b7d939f08ec58954 100644 (file)
 #include "gromacs/gpu_utils/device_context.h"
 #include "gromacs/gpu_utils/device_stream.h"
 #include "gromacs/gpu_utils/devicebuffer.h"
-#if GMX_GPU_CUDA
-#    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
-#elif GMX_GPU_SYCL
-#    include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
-#endif
+#include "gromacs/gpu_utils/gpueventsynchronizer.h"
 #include "gromacs/gpu_utils/gputraits.h"
 #include "gromacs/mdlib/leapfrog_gpu.h"
 #include "gromacs/mdlib/update_constrain_gpu.h"
 #include "gromacs/mdlib/update_constrain_gpu_internal.h"
 #include "gromacs/mdtypes/mdatom.h"
 #include "gromacs/timing/wallcycle.h"
+#include "gromacs/topology/mtop_util.h"
+
+static constexpr bool sc_haveGpuConstraintSupport = GMX_GPU_CUDA || GMX_GPU_SYCL;
 
 namespace gmx
 {
@@ -102,8 +101,11 @@ void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer*             fRead
     // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
     // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
     // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
-    lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
-    settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
+    if (sc_haveGpuConstraintSupport)
+    {
+        lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
+        settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
+    }
 
     // scaledVirial -> virial (methods above returns scaled values)
     float scaleFactor = 0.5F / (dt * dt);
@@ -115,12 +117,10 @@ void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer*             fRead
         }
     }
 
-    coordinatesReady_->markEvent(deviceStream_);
+    xUpdatedOnDeviceEvent_.markEvent(deviceStream_);
 
     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
-
-    return;
 }
 
 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
@@ -149,23 +149,20 @@ void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
 }
 
-UpdateConstrainGpu::Impl::Impl(const t_inputrec&     ir,
-                               const gmx_mtop_t&     mtop,
-                               const int             numTempScaleValues,
-                               const DeviceContext&  deviceContext,
-                               const DeviceStream&   deviceStream,
-                               GpuEventSynchronizer* xUpdatedOnDevice,
-                               gmx_wallcycle*        wcycle) :
-    deviceContext_(deviceContext),
-    deviceStream_(deviceStream),
-    coordinatesReady_(xUpdatedOnDevice),
-    wcycle_(wcycle)
+UpdateConstrainGpu::Impl::Impl(const t_inputrec&    ir,
+                               const gmx_mtop_t&    mtop,
+                               const int            numTempScaleValues,
+                               const DeviceContext& deviceContext,
+                               const DeviceStream&  deviceStream,
+                               gmx_wallcycle*       wcycle) :
+    deviceContext_(deviceContext), deviceStream_(deviceStream), wcycle_(wcycle)
 {
-    GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
-
     integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_, numTempScaleValues);
-    lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
-    settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
+    if (sc_haveGpuConstraintSupport)
+    {
+        lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
+        settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
+    }
 }
 
 UpdateConstrainGpu::Impl::~Impl() {}
@@ -196,8 +193,16 @@ void UpdateConstrainGpu::Impl::set(DeviceBuffer<Float3>          d_x,
 
     // Integrator should also update something, but it does not even have a method yet
     integrator_->set(numAtoms_, md.invmass, md.cTC);
-    lincsGpu_->set(idef, numAtoms_, md.invmass);
-    settleGpu_->set(idef);
+    if (sc_haveGpuConstraintSupport)
+    {
+        lincsGpu_->set(idef, numAtoms_, md.invmass);
+        settleGpu_->set(idef);
+    }
+    else
+    {
+        GMX_ASSERT(idef.il[F_SETTLE].empty(), "SETTLE not supported");
+        GMX_ASSERT(idef.il[F_CONSTR].empty(), "LINCS not supported");
+    }
 
     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
@@ -209,19 +214,18 @@ void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
     setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
 }
 
-GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
+GpuEventSynchronizer* UpdateConstrainGpu::Impl::xUpdatedOnDeviceEvent()
 {
-    return coordinatesReady_;
+    return &xUpdatedOnDeviceEvent_;
 }
 
-UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec&     ir,
-                                       const gmx_mtop_t&     mtop,
-                                       const int             numTempScaleValues,
-                                       const DeviceContext&  deviceContext,
-                                       const DeviceStream&   deviceStream,
-                                       GpuEventSynchronizer* xUpdatedOnDevice,
-                                       gmx_wallcycle*        wcycle) :
-    impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
+UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec&    ir,
+                                       const gmx_mtop_t&    mtop,
+                                       const int            numTempScaleValues,
+                                       const DeviceContext& deviceContext,
+                                       const DeviceStream&  deviceStream,
+                                       gmx_wallcycle*       wcycle) :
+    impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, wcycle))
 {
 }
 
@@ -274,9 +278,9 @@ void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
     impl_->setPbc(pbcType, box);
 }
 
-GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
+GpuEventSynchronizer* UpdateConstrainGpu::xUpdatedOnDeviceEvent()
 {
-    return impl_->getCoordinatesReadySync();
+    return impl_->xUpdatedOnDeviceEvent();
 }
 
 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
@@ -284,4 +288,9 @@ bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop
     return LincsGpu::isNumCoupledConstraintsSupported(mtop);
 }
 
+bool UpdateConstrainGpu::areConstraintsSupported()
+{
+    return sc_haveGpuConstraintSupport;
+}
+
 } // namespace gmx