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 Leap-Frog using CUDA
39 * This file contains implementation of basic Leap-Frog integrator
40 * using CUDA, including class initialization, data-structures management
43 * \author Artem Zhmurov <zhmurov@gmail.com>
45 * \ingroup module_mdlib
49 #include "leapfrog_gpu.h"
58 #include "gromacs/gpu_utils/cudautils.cuh"
59 #include "gromacs/gpu_utils/devicebuffer.h"
60 #include "gromacs/gpu_utils/vectype_ops.cuh"
61 #include "gromacs/math/vec.h"
62 #include "gromacs/mdtypes/group.h"
63 #include "gromacs/pbcutil/pbc.h"
64 #include "gromacs/pbcutil/pbc_aiuc_cuda.cuh"
65 #include "gromacs/utility/arrayref.h"
70 /*!\brief Number of CUDA threads in a block
72 * \todo Check if using smaller block size will lead to better performance.
74 constexpr static int c_threadsPerBlock = 256;
75 //! Maximum number of threads in a block (for __launch_bounds__)
76 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
78 /*! \brief Main kernel for Leap-Frog integrator.
80 * The coordinates and velocities are updated on the GPU. Also saves the intermediate values of the coordinates for
81 * further use in constraints.
83 * Each GPU thread works with a single particle. Empty declaration is needed to
84 * avoid "no previous prototype for function" clang warning.
86 * \todo Check if the force should be set to zero here.
87 * \todo This kernel can also accumulate incidental temperatures for each atom.
89 * \tparam numTempScaleValues The number of different T-couple values.
90 * \tparam velocityScaling Type of the Parrinello-Rahman velocity rescaling.
91 * \param[in] numAtoms Total number of atoms.
92 * \param[in,out] gm_x Coordinates to update upon integration.
93 * \param[out] gm_xp A copy of the coordinates before the integration (for constraints).
94 * \param[in,out] gm_v Velocities to update.
95 * \param[in] gm_f Atomic forces.
96 * \param[in] gm_inverseMasses Reciprocal masses.
97 * \param[in] dt Timestep.
98 * \param[in] gm_lambdas Temperature scaling factors (one per group)
99 * \param[in] gm_tempScaleGroups Mapping of atoms into groups.
100 * \param[in] prVelocityScalingMatrixDiagonal Diagonal elements of Parrinello-Rahman velocity scaling matrix
102 template<NumTempScaleValues numTempScaleValues, VelocityScalingType velocityScaling>
103 __launch_bounds__(c_maxThreadsPerBlock) __global__
104 void leapfrog_kernel(const int numAtoms,
105 float3* __restrict__ gm_x,
106 float3* __restrict__ gm_xp,
107 float3* __restrict__ gm_v,
108 const float3* __restrict__ gm_f,
109 const float* __restrict__ gm_inverseMasses,
111 const float* __restrict__ gm_lambdas,
112 const unsigned short* __restrict__ gm_tempScaleGroups,
113 const float3 prVelocityScalingMatrixDiagonal);
115 template<NumTempScaleValues numTempScaleValues, VelocityScalingType velocityScaling>
116 __launch_bounds__(c_maxThreadsPerBlock) __global__
117 void leapfrog_kernel(const int numAtoms,
118 float3* __restrict__ gm_x,
119 float3* __restrict__ gm_xp,
120 float3* __restrict__ gm_v,
121 const float3* __restrict__ gm_f,
122 const float* __restrict__ gm_inverseMasses,
124 const float* __restrict__ gm_lambdas,
125 const unsigned short* __restrict__ gm_tempScaleGroups,
126 const float3 prVelocityScalingMatrixDiagonal)
128 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
129 if (threadIndex < numAtoms)
131 float3 x = gm_x[threadIndex];
132 float3 v = gm_v[threadIndex];
133 float3 f = gm_f[threadIndex];
134 float im = gm_inverseMasses[threadIndex];
135 float imdt = im * dt;
137 // Swapping places for xp and x so that the x will contain the updated coordinates and xp - the
138 // coordinates before update. This should be taken into account when (if) constraints are applied
139 // after the update: x and xp have to be passed to constraints in the 'wrong' order.
141 gm_xp[threadIndex] = x;
143 if (numTempScaleValues != NumTempScaleValues::None || velocityScaling != VelocityScalingType::None)
147 if (numTempScaleValues != NumTempScaleValues::None)
150 if (numTempScaleValues == NumTempScaleValues::Single)
152 lambda = gm_lambdas[0];
154 else if (numTempScaleValues == NumTempScaleValues::Multiple)
156 int tempScaleGroup = gm_tempScaleGroups[threadIndex];
157 lambda = gm_lambdas[tempScaleGroup];
162 if (velocityScaling == VelocityScalingType::Diagonal)
164 vp.x -= prVelocityScalingMatrixDiagonal.x * v.x;
165 vp.y -= prVelocityScalingMatrixDiagonal.y * v.y;
166 vp.z -= prVelocityScalingMatrixDiagonal.z * v.z;
175 gm_v[threadIndex] = v;
176 gm_x[threadIndex] = x;
181 /*! \brief Select templated kernel.
183 * Returns pointer to a CUDA kernel based on the number of temperature coupling groups and
184 * whether or not the temperature and(or) pressure coupling is enabled.
186 * \param[in] doTemperatureScaling If the kernel with temperature coupling velocity scaling
187 * should be selected.
188 * \param[in] numTempScaleValues Number of temperature coupling groups in the system.
189 * \param[in] prVelocityScalingType Type of the Parrinello-Rahman velocity scaling.
191 * \retrun Pointer to CUDA kernel
193 inline auto selectLeapFrogKernelPtr(bool doTemperatureScaling,
194 int numTempScaleValues,
195 VelocityScalingType prVelocityScalingType)
197 // Check input for consistency: if there is temperature coupling, at least one coupling group should be defined.
198 GMX_ASSERT(!doTemperatureScaling || (numTempScaleValues > 0),
199 "Temperature coupling was requested with no temperature coupling groups.");
200 auto kernelPtr = leapfrog_kernel<NumTempScaleValues::None, VelocityScalingType::None>;
202 if (prVelocityScalingType == VelocityScalingType::None)
204 if (!doTemperatureScaling)
206 kernelPtr = leapfrog_kernel<NumTempScaleValues::None, VelocityScalingType::None>;
208 else if (numTempScaleValues == 1)
210 kernelPtr = leapfrog_kernel<NumTempScaleValues::Single, VelocityScalingType::None>;
212 else if (numTempScaleValues > 1)
214 kernelPtr = leapfrog_kernel<NumTempScaleValues::Multiple, VelocityScalingType::None>;
217 else if (prVelocityScalingType == VelocityScalingType::Diagonal)
219 if (!doTemperatureScaling)
221 kernelPtr = leapfrog_kernel<NumTempScaleValues::None, VelocityScalingType::Diagonal>;
223 else if (numTempScaleValues == 1)
225 kernelPtr = leapfrog_kernel<NumTempScaleValues::Single, VelocityScalingType::Diagonal>;
227 else if (numTempScaleValues > 1)
229 kernelPtr = leapfrog_kernel<NumTempScaleValues::Multiple, VelocityScalingType::Diagonal>;
234 GMX_RELEASE_ASSERT(false,
235 "Only isotropic Parrinello-Rahman pressure coupling is supported.");
240 void LeapFrogGpu::integrate(const DeviceBuffer<float3> d_x,
241 DeviceBuffer<float3> d_xp,
242 DeviceBuffer<float3> d_v,
243 const DeviceBuffer<float3> d_f,
245 const bool doTemperatureScaling,
246 gmx::ArrayRef<const t_grp_tcstat> tcstat,
247 const bool doParrinelloRahman,
248 const float dtPressureCouple,
249 const matrix prVelocityScalingMatrix)
252 ensureNoPendingCudaError("In CUDA version of Leap-Frog integrator");
254 auto kernelPtr = leapfrog_kernel<NumTempScaleValues::None, VelocityScalingType::None>;
255 if (doTemperatureScaling || doParrinelloRahman)
257 if (doTemperatureScaling)
259 GMX_ASSERT(numTempScaleValues_ == ssize(h_lambdas_),
260 "Number of temperature scaling factors changed since it was set for the "
262 for (int i = 0; i < numTempScaleValues_; i++)
264 h_lambdas_[i] = tcstat[i].lambda;
266 copyToDeviceBuffer(&d_lambdas_, h_lambdas_.data(), 0, numTempScaleValues_,
267 deviceStream_, GpuApiCallBehavior::Async, nullptr);
269 VelocityScalingType prVelocityScalingType = VelocityScalingType::None;
270 if (doParrinelloRahman)
272 prVelocityScalingType = VelocityScalingType::Diagonal;
273 GMX_ASSERT(prVelocityScalingMatrix[YY][XX] == 0 && prVelocityScalingMatrix[ZZ][XX] == 0
274 && prVelocityScalingMatrix[ZZ][YY] == 0
275 && prVelocityScalingMatrix[XX][YY] == 0
276 && prVelocityScalingMatrix[XX][ZZ] == 0
277 && prVelocityScalingMatrix[YY][ZZ] == 0,
278 "Fully anisotropic Parrinello-Rahman pressure coupling is not yet supported "
279 "in GPU version of Leap-Frog integrator.");
280 prVelocityScalingMatrixDiagonal_ =
281 make_float3(dtPressureCouple * prVelocityScalingMatrix[XX][XX],
282 dtPressureCouple * prVelocityScalingMatrix[YY][YY],
283 dtPressureCouple * prVelocityScalingMatrix[ZZ][ZZ]);
285 kernelPtr = selectLeapFrogKernelPtr(doTemperatureScaling, numTempScaleValues_, prVelocityScalingType);
288 const auto kernelArgs = prepareGpuKernelArguments(
289 kernelPtr, kernelLaunchConfig_, &numAtoms_, &d_x, &d_xp, &d_v, &d_f, &d_inverseMasses_,
290 &dt, &d_lambdas_, &d_tempScaleGroups_, &prVelocityScalingMatrixDiagonal_);
291 launchGpuKernel(kernelPtr, kernelLaunchConfig_, deviceStream_, nullptr, "leapfrog_kernel", kernelArgs);
296 LeapFrogGpu::LeapFrogGpu(const DeviceContext& deviceContext, const DeviceStream& deviceStream) :
297 deviceContext_(deviceContext),
298 deviceStream_(deviceStream)
302 changePinningPolicy(&h_lambdas_, gmx::PinningPolicy::PinnedIfSupported);
304 kernelLaunchConfig_.blockSize[0] = c_threadsPerBlock;
305 kernelLaunchConfig_.blockSize[1] = 1;
306 kernelLaunchConfig_.blockSize[2] = 1;
307 kernelLaunchConfig_.sharedMemorySize = 0;
310 LeapFrogGpu::~LeapFrogGpu()
312 freeDeviceBuffer(&d_inverseMasses_);
315 void LeapFrogGpu::set(const int numAtoms,
316 const real* inverseMasses,
317 const int numTempScaleValues,
318 const unsigned short* tempScaleGroups)
320 numAtoms_ = numAtoms;
321 kernelLaunchConfig_.gridSize[0] = (numAtoms_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
323 numTempScaleValues_ = numTempScaleValues;
325 reallocateDeviceBuffer(&d_inverseMasses_, numAtoms_, &numInverseMasses_,
326 &numInverseMassesAlloc_, deviceContext_);
327 copyToDeviceBuffer(&d_inverseMasses_, (float*)inverseMasses, 0, numAtoms_, deviceStream_,
328 GpuApiCallBehavior::Sync, nullptr);
330 // Temperature scale group map only used if there are more then one group
331 if (numTempScaleValues > 1)
333 reallocateDeviceBuffer(&d_tempScaleGroups_, numAtoms_, &numTempScaleGroups_,
334 &numTempScaleGroupsAlloc_, deviceContext_);
335 copyToDeviceBuffer(&d_tempScaleGroups_, tempScaleGroups, 0, numAtoms_, deviceStream_,
336 GpuApiCallBehavior::Sync, nullptr);
339 // If the temperature coupling is enabled, we need to make space for scaling factors
340 if (numTempScaleValues_ > 0)
342 h_lambdas_.resize(numTempScaleValues);
343 reallocateDeviceBuffer(&d_lambdas_, numTempScaleValues_, &numLambdas_, &numLambdasAlloc_,