Move ScalingMatrix to common GPU update code.
authorAndrey Alekseenko <al42and@gmail.com>
Fri, 9 Apr 2021 10:54:46 +0000 (13:54 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Sun, 11 Apr 2021 05:42:58 +0000 (05:42 +0000)
Also, add a constructor to initialize it from an old-style matrix.

Currently, this struct is used only in CUDA, but would be used in SYCL
soon.

Refs #3932.

src/gromacs/mdlib/update_constrain_gpu_impl.cu
src/gromacs/mdlib/update_constrain_gpu_impl.h

index 3b428a183b591af953d6100f31ec8b7975dd3b1a..5fc4fb86092938a64cfcd6687904792916728b6e 100644 (file)
@@ -80,15 +80,6 @@ constexpr static int c_threadsPerBlock = 256;
 //! Maximum number of threads in a block (for __launch_bounds__)
 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
 
-/*! \brief Scaling matrix struct.
- *
- * \todo Should be generalized.
- */
-struct ScalingMatrix
-{
-    float xx, yy, zz, yx, zx, zy;
-};
-
 __launch_bounds__(c_maxThreadsPerBlock) __global__
         static void scaleCoordinates_kernel(const int numAtoms,
                                             float3* __restrict__ gm_x,
@@ -161,13 +152,7 @@ void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
 
-    ScalingMatrix mu;
-    mu.xx = scalingMatrix[XX][XX];
-    mu.yy = scalingMatrix[YY][YY];
-    mu.zz = scalingMatrix[ZZ][ZZ];
-    mu.yx = scalingMatrix[YY][XX];
-    mu.zx = scalingMatrix[ZZ][XX];
-    mu.zy = scalingMatrix[ZZ][YY];
+    ScalingMatrix mu(scalingMatrix);
 
     const auto kernelArgs = prepareGpuKernelArguments(
             scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_x_, &mu);
@@ -190,13 +175,7 @@ void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
 
-    ScalingMatrix mu;
-    mu.xx = scalingMatrix[XX][XX];
-    mu.yy = scalingMatrix[YY][YY];
-    mu.zz = scalingMatrix[ZZ][ZZ];
-    mu.yx = scalingMatrix[YY][XX];
-    mu.zx = scalingMatrix[ZZ][XX];
-    mu.zy = scalingMatrix[ZZ][YY];
+    ScalingMatrix mu(scalingMatrix);
 
     const auto kernelArgs = prepareGpuKernelArguments(
             scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_v_, &mu);
index 1d1a8fd0de1a8453e6a20a022fedcf5da60c54bb..ece9859374200fe10bec26895f028d7a0c3fe80f 100644 (file)
@@ -228,6 +228,24 @@ private:
     gmx_wallcycle* wcycle_ = nullptr;
 };
 
+/*! \brief Scaling matrix struct.
+ *
+ * \todo Should be generalized.
+ */
+struct ScalingMatrix
+{
+    ScalingMatrix(const matrix m) :
+        xx(m[XX][XX]),
+        yy(m[YY][YY]),
+        zz(m[ZZ][ZZ]),
+        yx(m[YY][XX]),
+        zx(m[ZZ][XX]),
+        zy(m[ZZ][YY])
+    {
+    }
+    float xx, yy, zz, yx, zx, zy;
+};
+
 } // namespace gmx
 
 #endif // GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H