2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 * \brief Implements update and constraints class using CUDA.
39 * The class combines Leap-Frog integrator with LINCS and SETTLE constraints.
41 * \todo The computational procedures in members should be integrated to improve
42 * computational performance.
44 * \author Artem Zhmurov <zhmurov@gmail.com>
46 * \ingroup module_mdlib
50 #include "update_constrain_gpu_impl.h"
59 #include "gromacs/gpu_utils/cudautils.cuh"
60 #include "gromacs/gpu_utils/device_context.h"
61 #include "gromacs/gpu_utils/device_stream.h"
62 #include "gromacs/gpu_utils/devicebuffer.h"
63 #include "gromacs/gpu_utils/gputraits.cuh"
64 #include "gromacs/gpu_utils/vectype_ops.cuh"
65 #include "gromacs/mdlib/leapfrog_gpu.h"
66 #include "gromacs/mdlib/lincs_gpu.cuh"
67 #include "gromacs/mdlib/settle_gpu.cuh"
68 #include "gromacs/mdlib/update_constrain_gpu.h"
69 #include "gromacs/mdtypes/mdatom.h"
70 #include "gromacs/timing/wallcycle.h"
74 /*!\brief Number of CUDA threads in a block
76 * \todo Check if using smaller block size will lead to better performance.
78 constexpr static int c_threadsPerBlock = 256;
79 //! Maximum number of threads in a block (for __launch_bounds__)
80 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
82 /*! \brief Scaling matrix struct.
84 * \todo Should be generalized.
88 float xx, yy, zz, yx, zx, zy;
91 __launch_bounds__(c_maxThreadsPerBlock) __global__
92 static void scaleCoordinates_kernel(const int numAtoms,
93 float3* __restrict__ gm_x,
94 const ScalingMatrix scalingMatrix)
96 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
97 if (threadIndex < numAtoms)
99 float3 x = gm_x[threadIndex];
101 x.x = scalingMatrix.xx * x.x + scalingMatrix.yx * x.y + scalingMatrix.zx * x.z;
102 x.y = scalingMatrix.yy * x.y + scalingMatrix.zy * x.z;
103 x.z = scalingMatrix.zz * x.z;
105 gm_x[threadIndex] = x;
109 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer* fReadyOnDevice,
111 const bool updateVelocities,
112 const bool computeVirial,
114 const bool doTemperatureScaling,
115 gmx::ArrayRef<const t_grp_tcstat> tcstat,
116 const bool doParrinelloRahman,
117 const float dtPressureCouple,
118 const matrix prVelocityScalingMatrix)
120 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
121 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
123 // Clearing virial matrix
124 // TODO There is no point in having separate virial matrix for constraints
127 // Make sure that the forces are ready on device before proceeding with the update.
128 fReadyOnDevice->enqueueWaitEvent(deviceStream_);
130 // The integrate should save a copy of the current coordinates in d_xp_ and write updated
131 // once into d_x_. The d_xp_ is only needed by constraints.
132 integrator_->integrate(d_x_, d_xp_, d_v_, d_f_, dt, doTemperatureScaling, tcstat,
133 doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
134 // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
135 // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
136 // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
137 lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
138 settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
140 // scaledVirial -> virial (methods above returns scaled values)
141 float scaleFactor = 0.5f / (dt * dt);
142 for (int i = 0; i < DIM; i++)
144 for (int j = 0; j < DIM; j++)
146 virial[i][j] = scaleFactor * virial[i][j];
150 coordinatesReady_->markEvent(deviceStream_);
152 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
153 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
158 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
160 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
161 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
164 mu.xx = scalingMatrix[XX][XX];
165 mu.yy = scalingMatrix[YY][YY];
166 mu.zz = scalingMatrix[ZZ][ZZ];
167 mu.yx = scalingMatrix[YY][XX];
168 mu.zx = scalingMatrix[ZZ][XX];
169 mu.zy = scalingMatrix[ZZ][YY];
171 const auto kernelArgs = prepareGpuKernelArguments(
172 scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_x_, &mu);
173 launchGpuKernel(scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, deviceStream_,
174 nullptr, "scaleCoordinates_kernel", kernelArgs);
175 // TODO: Although this only happens on the pressure coupling steps, this synchronization
176 // can affect the performance if nstpcouple is small.
177 deviceStream_.synchronize();
179 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
180 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
183 void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
185 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
186 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
189 mu.xx = scalingMatrix[XX][XX];
190 mu.yy = scalingMatrix[YY][YY];
191 mu.zz = scalingMatrix[ZZ][ZZ];
192 mu.yx = scalingMatrix[YY][XX];
193 mu.zx = scalingMatrix[ZZ][XX];
194 mu.zy = scalingMatrix[ZZ][YY];
196 const auto kernelArgs = prepareGpuKernelArguments(
197 scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_v_, &mu);
198 launchGpuKernel(scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, deviceStream_,
199 nullptr, "scaleCoordinates_kernel", kernelArgs);
200 // TODO: Although this only happens on the pressure coupling steps, this synchronization
201 // can affect the performance if nstpcouple is small.
202 deviceStream_.synchronize();
204 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
205 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
208 UpdateConstrainGpu::Impl::Impl(const t_inputrec& ir,
209 const gmx_mtop_t& mtop,
210 const DeviceContext& deviceContext,
211 const DeviceStream& deviceStream,
212 GpuEventSynchronizer* xUpdatedOnDevice,
213 gmx_wallcycle* wcycle) :
214 deviceContext_(deviceContext),
215 deviceStream_(deviceStream),
216 coordinatesReady_(xUpdatedOnDevice),
219 GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
222 integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_);
223 lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
224 settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
226 coordinateScalingKernelLaunchConfig_.blockSize[0] = c_threadsPerBlock;
227 coordinateScalingKernelLaunchConfig_.blockSize[1] = 1;
228 coordinateScalingKernelLaunchConfig_.blockSize[2] = 1;
229 coordinateScalingKernelLaunchConfig_.sharedMemorySize = 0;
232 UpdateConstrainGpu::Impl::~Impl() {}
234 void UpdateConstrainGpu::Impl::set(DeviceBuffer<RVec> d_x,
235 DeviceBuffer<RVec> d_v,
236 const DeviceBuffer<RVec> d_f,
237 const InteractionDefinitions& idef,
239 const int numTempScaleValues)
242 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
243 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
245 GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
246 GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
247 GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
249 d_x_ = reinterpret_cast<float3*>(d_x);
250 d_v_ = reinterpret_cast<float3*>(d_v);
251 d_f_ = reinterpret_cast<float3*>(d_f);
255 reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
257 reallocateDeviceBuffer(&d_inverseMasses_, numAtoms_, &numInverseMasses_,
258 &numInverseMassesAlloc_, deviceContext_);
260 // Integrator should also update something, but it does not even have a method yet
261 integrator_->set(numAtoms_, md.invmass, numTempScaleValues, md.cTC);
262 lincsGpu_->set(idef, numAtoms_, md.invmass);
263 settleGpu_->set(idef);
265 coordinateScalingKernelLaunchConfig_.gridSize[0] =
266 (numAtoms_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
268 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
269 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
272 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
275 setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
278 GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
280 return coordinatesReady_;
283 UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec& ir,
284 const gmx_mtop_t& mtop,
285 const DeviceContext& deviceContext,
286 const DeviceStream& deviceStream,
287 GpuEventSynchronizer* xUpdatedOnDevice,
288 gmx_wallcycle* wcycle) :
289 impl_(new Impl(ir, mtop, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
293 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
295 void UpdateConstrainGpu::integrate(GpuEventSynchronizer* fReadyOnDevice,
297 const bool updateVelocities,
298 const bool computeVirial,
300 const bool doTemperatureScaling,
301 gmx::ArrayRef<const t_grp_tcstat> tcstat,
302 const bool doParrinelloRahman,
303 const float dtPressureCouple,
304 const matrix prVelocityScalingMatrix)
306 impl_->integrate(fReadyOnDevice, dt, updateVelocities, computeVirial, virialScaled, doTemperatureScaling,
307 tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
310 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
312 impl_->scaleCoordinates(scalingMatrix);
315 void UpdateConstrainGpu::scaleVelocities(const matrix scalingMatrix)
317 impl_->scaleVelocities(scalingMatrix);
320 void UpdateConstrainGpu::set(DeviceBuffer<RVec> d_x,
321 DeviceBuffer<RVec> d_v,
322 const DeviceBuffer<RVec> d_f,
323 const InteractionDefinitions& idef,
325 const int numTempScaleValues)
327 impl_->set(d_x, d_v, d_f, idef, md, numTempScaleValues);
330 void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
332 impl_->setPbc(pbcType, box);
335 GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
337 return impl_->getCoordinatesReadySync();
340 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
342 return LincsGpu::isNumCoupledConstraintsSupported(mtop);