2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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/gpueventsynchronizer.cuh"
64 #include "gromacs/gpu_utils/gputraits.cuh"
65 #include "gromacs/gpu_utils/vectype_ops.cuh"
66 #include "gromacs/mdlib/leapfrog_gpu.h"
67 #include "gromacs/mdlib/lincs_gpu.cuh"
68 #include "gromacs/mdlib/settle_gpu.cuh"
69 #include "gromacs/mdlib/update_constrain_gpu.h"
70 #include "gromacs/mdtypes/mdatom.h"
71 #include "gromacs/timing/wallcycle.h"
75 /*!\brief Number of CUDA threads in a block
77 * \todo Check if using smaller block size will lead to better performance.
79 constexpr static int c_threadsPerBlock = 256;
80 //! Maximum number of threads in a block (for __launch_bounds__)
81 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
83 __launch_bounds__(c_maxThreadsPerBlock) __global__
84 static void scaleCoordinates_kernel(const int numAtoms,
85 float3* __restrict__ gm_x,
86 const ScalingMatrix scalingMatrix)
88 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
89 if (threadIndex < numAtoms)
91 float3 x = gm_x[threadIndex];
93 x.x = scalingMatrix.xx * x.x + scalingMatrix.yx * x.y + scalingMatrix.zx * x.z;
94 x.y = scalingMatrix.yy * x.y + scalingMatrix.zy * x.z;
95 x.z = scalingMatrix.zz * x.z;
97 gm_x[threadIndex] = x;
101 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer* fReadyOnDevice,
103 const bool updateVelocities,
104 const bool computeVirial,
106 const bool doTemperatureScaling,
107 gmx::ArrayRef<const t_grp_tcstat> tcstat,
108 const bool doParrinelloRahman,
109 const float dtPressureCouple,
110 const matrix prVelocityScalingMatrix)
112 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
113 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
115 // Clearing virial matrix
116 // TODO There is no point in having separate virial matrix for constraints
119 // Make sure that the forces are ready on device before proceeding with the update.
120 fReadyOnDevice->enqueueWaitEvent(deviceStream_);
122 // The integrate should save a copy of the current coordinates in d_xp_ and write updated
123 // once into d_x_. The d_xp_ is only needed by constraints.
124 integrator_->integrate(
125 d_x_, d_xp_, d_v_, d_f_, dt, doTemperatureScaling, tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
126 // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
127 // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
128 // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
129 lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
130 settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
132 // scaledVirial -> virial (methods above returns scaled values)
133 float scaleFactor = 0.5f / (dt * dt);
134 for (int i = 0; i < DIM; i++)
136 for (int j = 0; j < DIM; j++)
138 virial[i][j] = scaleFactor * virial[i][j];
142 coordinatesReady_->markEvent(deviceStream_);
144 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
145 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
150 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
152 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
153 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
155 ScalingMatrix mu(scalingMatrix);
157 const auto kernelArgs = prepareGpuKernelArguments(
158 scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_x_, &mu);
159 launchGpuKernel(scaleCoordinates_kernel,
160 coordinateScalingKernelLaunchConfig_,
163 "scaleCoordinates_kernel",
165 // TODO: Although this only happens on the pressure coupling steps, this synchronization
166 // can affect the performance if nstpcouple is small.
167 deviceStream_.synchronize();
169 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
170 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
173 void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
175 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
176 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
178 ScalingMatrix mu(scalingMatrix);
180 const auto kernelArgs = prepareGpuKernelArguments(
181 scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_v_, &mu);
182 launchGpuKernel(scaleCoordinates_kernel,
183 coordinateScalingKernelLaunchConfig_,
186 "scaleCoordinates_kernel",
188 // TODO: Although this only happens on the pressure coupling steps, this synchronization
189 // can affect the performance if nstpcouple is small.
190 deviceStream_.synchronize();
192 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
193 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
196 UpdateConstrainGpu::Impl::Impl(const t_inputrec& ir,
197 const gmx_mtop_t& mtop,
198 const int numTempScaleValues,
199 const DeviceContext& deviceContext,
200 const DeviceStream& deviceStream,
201 GpuEventSynchronizer* xUpdatedOnDevice,
202 gmx_wallcycle* wcycle) :
203 deviceContext_(deviceContext),
204 deviceStream_(deviceStream),
205 coordinatesReady_(xUpdatedOnDevice),
208 GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
211 integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_, numTempScaleValues);
212 lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
213 settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
215 coordinateScalingKernelLaunchConfig_.blockSize[0] = c_threadsPerBlock;
216 coordinateScalingKernelLaunchConfig_.blockSize[1] = 1;
217 coordinateScalingKernelLaunchConfig_.blockSize[2] = 1;
218 coordinateScalingKernelLaunchConfig_.sharedMemorySize = 0;
221 UpdateConstrainGpu::Impl::~Impl() {}
223 void UpdateConstrainGpu::Impl::set(DeviceBuffer<Float3> d_x,
224 DeviceBuffer<Float3> d_v,
225 const DeviceBuffer<Float3> d_f,
226 const InteractionDefinitions& idef,
230 wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
231 wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
233 GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
234 GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
235 GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
243 reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
245 reallocateDeviceBuffer(
246 &d_inverseMasses_, numAtoms_, &numInverseMasses_, &numInverseMassesAlloc_, deviceContext_);
248 // Integrator should also update something, but it does not even have a method yet
249 integrator_->set(numAtoms_, md.invmass, md.cTC);
250 lincsGpu_->set(idef, numAtoms_, md.invmass);
251 settleGpu_->set(idef);
253 coordinateScalingKernelLaunchConfig_.gridSize[0] =
254 (numAtoms_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
256 wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
257 wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
260 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
263 setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
266 GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
268 return coordinatesReady_;
271 UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec& ir,
272 const gmx_mtop_t& mtop,
273 const int numTempScaleValues,
274 const DeviceContext& deviceContext,
275 const DeviceStream& deviceStream,
276 GpuEventSynchronizer* xUpdatedOnDevice,
277 gmx_wallcycle* wcycle) :
278 impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
282 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
284 void UpdateConstrainGpu::integrate(GpuEventSynchronizer* fReadyOnDevice,
286 const bool updateVelocities,
287 const bool computeVirial,
289 const bool doTemperatureScaling,
290 gmx::ArrayRef<const t_grp_tcstat> tcstat,
291 const bool doParrinelloRahman,
292 const float dtPressureCouple,
293 const matrix prVelocityScalingMatrix)
295 impl_->integrate(fReadyOnDevice,
300 doTemperatureScaling,
304 prVelocityScalingMatrix);
307 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
309 impl_->scaleCoordinates(scalingMatrix);
312 void UpdateConstrainGpu::scaleVelocities(const matrix scalingMatrix)
314 impl_->scaleVelocities(scalingMatrix);
317 void UpdateConstrainGpu::set(DeviceBuffer<Float3> d_x,
318 DeviceBuffer<Float3> d_v,
319 const DeviceBuffer<Float3> d_f,
320 const InteractionDefinitions& idef,
323 impl_->set(d_x, d_v, d_f, idef, md);
326 void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
328 impl_->setPbc(pbcType, box);
331 GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
333 return impl_->getCoordinatesReadySync();
336 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
338 return LincsGpu::isNumCoupledConstraintsSupported(mtop);