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 SETTLE using CUDA
39 * This file contains implementation of SETTLE constraints algorithm
40 * using CUDA, including class initialization, data-structures management
44 * \author Artem Zhmurov <zhmurov@gmail.com>
46 * \ingroup module_mdlib
50 #include "settle_gpu.cuh"
59 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
60 #include "gromacs/gpu_utils/cudautils.cuh"
61 #include "gromacs/gpu_utils/devicebuffer.h"
62 #include "gromacs/gpu_utils/gputraits.cuh"
63 #include "gromacs/gpu_utils/vectype_ops.cuh"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/pbcutil/pbc.h"
66 #include "gromacs/pbcutil/pbc_aiuc_cuda.cuh"
71 //! Number of CUDA threads in a block
72 constexpr static int c_threadsPerBlock = 256;
73 //! Maximum number of threads in a block (for __launch_bounds__)
74 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
76 /*! \brief SETTLE constraints kernel
78 * Each thread corresponds to a single constraints triangle (i.e. single water molecule).
80 * See original CPU version in settle.cpp
82 * \param [in] numSettles Number of constraints triangles (water molecules).
83 * \param [in] gm_settles Indexes of three atoms in the constraints triangle. The field .x of int3
84 * data type corresponds to Oxygen, fields .y and .z are two hydrogen atoms.
85 * \param [in] pars Parameters for the algorithm (i.e. masses, target distances, etc.).
86 * \param [in] gm_x Coordinates of atoms before the timestep.
87 * \param [in,out] gm_x Coordinates of atoms after the timestep (constrained coordinates will be
89 * \param [in] invdt Reciprocal timestep.
90 * \param [in] gm_v Velocities of the particles.
91 * \param [in] gm_virialScaled Virial tensor.
92 * \param [in] pbcAiuc Periodic boundary conditions data.
94 template<bool updateVelocities, bool computeVirial>
95 __launch_bounds__(c_maxThreadsPerBlock) __global__
96 void settle_kernel(const int numSettles,
97 const int3* __restrict__ gm_settles,
98 const SettleParameters pars,
99 const float3* __restrict__ gm_x,
100 float3* __restrict__ gm_xprime,
102 float3* __restrict__ gm_v,
103 float* __restrict__ gm_virialScaled,
104 const PbcAiuc pbcAiuc)
106 /* ******************************************************************* */
108 /* Original code by Shuichi Miyamoto, last update Oct. 1, 1992 ** */
110 /* Algorithm changes by Berk Hess: ** */
111 /* 2004-07-15 Convert COM to double precision to avoid drift ** */
112 /* 2006-10-16 Changed velocity update to use differences ** */
113 /* 2012-09-24 Use oxygen as reference instead of COM ** */
114 /* 2016-02 Complete rewrite of the code for SIMD ** */
116 /* Reference for the SETTLE algorithm ** */
117 /* S. Miyamoto et al., J. Comp. Chem., 13, 952 (1992). ** */
119 /* ******************************************************************* */
121 constexpr float almost_zero = real(1e-12);
123 extern __shared__ float sm_threadVirial[];
125 int tid = static_cast<int>(blockIdx.x * blockDim.x + threadIdx.x);
127 if (tid < numSettles)
129 // These are the indexes of three atoms in a single 'water' molecule.
130 // TODO Can be reduced to one integer if atoms are consecutive in memory.
131 int3 indices = gm_settles[tid];
133 float3 x_ow1 = gm_x[indices.x];
134 float3 x_hw2 = gm_x[indices.y];
135 float3 x_hw3 = gm_x[indices.z];
137 float3 xprime_ow1 = gm_xprime[indices.x];
138 float3 xprime_hw2 = gm_xprime[indices.y];
139 float3 xprime_hw3 = gm_xprime[indices.z];
141 float3 dist21 = pbcDxAiuc(pbcAiuc, x_hw2, x_ow1);
142 float3 dist31 = pbcDxAiuc(pbcAiuc, x_hw3, x_ow1);
143 float3 doh2 = pbcDxAiuc(pbcAiuc, xprime_hw2, xprime_ow1);
145 float3 sh_hw2 = xprime_hw2 - (xprime_ow1 + doh2);
146 xprime_hw2 = xprime_hw2 - sh_hw2;
148 float3 doh3 = pbcDxAiuc(pbcAiuc, xprime_hw3, xprime_ow1);
150 float3 sh_hw3 = xprime_hw3 - (xprime_ow1 + doh3);
151 xprime_hw3 = xprime_hw3 - sh_hw3;
153 float3 a1 = (-doh2 - doh3) * pars.wh;
154 float3 com = xprime_ow1 - a1;
156 float3 b1 = xprime_hw2 - com;
158 float3 c1 = xprime_hw3 - com;
160 float xakszd = dist21.y * dist31.z - dist21.z * dist31.y;
161 float yakszd = dist21.z * dist31.x - dist21.x * dist31.z;
162 float zakszd = dist21.x * dist31.y - dist21.y * dist31.x;
164 float xaksxd = a1.y * zakszd - a1.z * yakszd;
165 float yaksxd = a1.z * xakszd - a1.x * zakszd;
166 float zaksxd = a1.x * yakszd - a1.y * xakszd;
168 float xaksyd = yakszd * zaksxd - zakszd * yaksxd;
169 float yaksyd = zakszd * xaksxd - xakszd * zaksxd;
170 float zaksyd = xakszd * yaksxd - yakszd * xaksxd;
172 float axlng = rsqrt(xaksxd * xaksxd + yaksxd * yaksxd + zaksxd * zaksxd);
173 float aylng = rsqrt(xaksyd * xaksyd + yaksyd * yaksyd + zaksyd * zaksyd);
174 float azlng = rsqrt(xakszd * xakszd + yakszd * yakszd + zakszd * zakszd);
176 // TODO {1,2,3} indexes should be swapped with {.x, .y, .z} components.
177 // This way, we will be able to use vector ops more.
178 float3 trns1, trns2, trns3;
180 trns1.x = xaksxd * axlng;
181 trns2.x = yaksxd * axlng;
182 trns3.x = zaksxd * axlng;
184 trns1.y = xaksyd * aylng;
185 trns2.y = yaksyd * aylng;
186 trns3.y = zaksyd * aylng;
188 trns1.z = xakszd * azlng;
189 trns2.z = yakszd * azlng;
190 trns3.z = zakszd * azlng;
195 b0d.x = trns1.x * dist21.x + trns2.x * dist21.y + trns3.x * dist21.z;
196 b0d.y = trns1.y * dist21.x + trns2.y * dist21.y + trns3.y * dist21.z;
198 c0d.x = trns1.x * dist31.x + trns2.x * dist31.y + trns3.x * dist31.z;
199 c0d.y = trns1.y * dist31.x + trns2.y * dist31.y + trns3.y * dist31.z;
203 float a1d_z = trns1.z * a1.x + trns2.z * a1.y + trns3.z * a1.z;
205 b1d.x = trns1.x * b1.x + trns2.x * b1.y + trns3.x * b1.z;
206 b1d.y = trns1.y * b1.x + trns2.y * b1.y + trns3.y * b1.z;
207 b1d.z = trns1.z * b1.x + trns2.z * b1.y + trns3.z * b1.z;
209 c1d.x = trns1.x * c1.x + trns2.x * c1.y + trns3.x * c1.z;
210 c1d.y = trns1.y * c1.x + trns2.y * c1.y + trns3.y * c1.z;
211 c1d.z = trns1.z * c1.x + trns2.z * c1.y + trns3.z * c1.z;
214 float sinphi = a1d_z * rsqrt(pars.ra * pars.ra);
215 float tmp2 = 1.0f - sinphi * sinphi;
217 if (almost_zero > tmp2)
222 float tmp = rsqrt(tmp2);
223 float cosphi = tmp2 * tmp;
224 float sinpsi = (b1d.z - c1d.z) * pars.irc2 * tmp;
225 tmp2 = 1.0f - sinpsi * sinpsi;
227 float cospsi = tmp2 * rsqrt(tmp2);
229 float a2d_y = pars.ra * cosphi;
230 float b2d_x = -pars.rc * cospsi;
231 float t1 = -pars.rb * cosphi;
232 float t2 = pars.rc * sinpsi * sinphi;
233 float b2d_y = t1 - t2;
234 float c2d_y = t1 + t2;
236 /* --- Step3 al,be,ga --- */
237 float alpha = b2d_x * (b0d.x - c0d.x) + b0d.y * b2d_y + c0d.y * c2d_y;
238 float beta = b2d_x * (c0d.y - b0d.y) + b0d.x * b2d_y + c0d.x * c2d_y;
239 float gamma = b0d.x * b1d.y - b1d.x * b0d.y + c0d.x * c1d.y - c1d.x * c0d.y;
240 float al2be2 = alpha * alpha + beta * beta;
241 tmp2 = (al2be2 - gamma * gamma);
242 float sinthe = (alpha * gamma - beta * tmp2 * rsqrt(tmp2)) * rsqrt(al2be2 * al2be2);
244 /* --- Step4 A3' --- */
245 tmp2 = 1.0f - sinthe * sinthe;
246 float costhe = tmp2 * rsqrt(tmp2);
248 float3 a3d, b3d, c3d;
250 a3d.x = -a2d_y * sinthe;
251 a3d.y = a2d_y * costhe;
253 b3d.x = b2d_x * costhe - b2d_y * sinthe;
254 b3d.y = b2d_x * sinthe + b2d_y * costhe;
256 c3d.x = -b2d_x * costhe - c2d_y * sinthe;
257 c3d.y = -b2d_x * sinthe + c2d_y * costhe;
260 /* --- Step5 A3 --- */
263 a3.x = trns1.x * a3d.x + trns1.y * a3d.y + trns1.z * a3d.z;
264 a3.y = trns2.x * a3d.x + trns2.y * a3d.y + trns2.z * a3d.z;
265 a3.z = trns3.x * a3d.x + trns3.y * a3d.y + trns3.z * a3d.z;
267 b3.x = trns1.x * b3d.x + trns1.y * b3d.y + trns1.z * b3d.z;
268 b3.y = trns2.x * b3d.x + trns2.y * b3d.y + trns2.z * b3d.z;
269 b3.z = trns3.x * b3d.x + trns3.y * b3d.y + trns3.z * b3d.z;
271 c3.x = trns1.x * c3d.x + trns1.y * c3d.y + trns1.z * c3d.z;
272 c3.y = trns2.x * c3d.x + trns2.y * c3d.y + trns2.z * c3d.z;
273 c3.z = trns3.x * c3d.x + trns3.y * c3d.y + trns3.z * c3d.z;
276 /* Compute and store the corrected new coordinate */
277 xprime_ow1 = com + a3;
278 xprime_hw2 = com + b3 + sh_hw2;
279 xprime_hw3 = com + c3 + sh_hw3;
281 gm_xprime[indices.x] = xprime_ow1;
282 gm_xprime[indices.y] = xprime_hw2;
283 gm_xprime[indices.z] = xprime_hw3;
286 if (updateVelocities || computeVirial)
293 if (updateVelocities)
296 float3 v_ow1 = gm_v[indices.x];
297 float3 v_hw2 = gm_v[indices.y];
298 float3 v_hw3 = gm_v[indices.z];
300 /* Add the position correction divided by dt to the velocity */
301 v_ow1 = da * invdt + v_ow1;
302 v_hw2 = db * invdt + v_hw2;
303 v_hw3 = dc * invdt + v_hw3;
305 gm_v[indices.x] = v_ow1;
306 gm_v[indices.y] = v_hw2;
307 gm_v[indices.z] = v_hw3;
313 float3 mdb = pars.mH * db;
314 float3 mdc = pars.mH * dc;
315 float3 mdo = pars.mO * da + mdb + mdc;
317 sm_threadVirial[0 * blockDim.x + threadIdx.x] =
318 -(x_ow1.x * mdo.x + dist21.x * mdb.x + dist31.x * mdc.x);
319 sm_threadVirial[1 * blockDim.x + threadIdx.x] =
320 -(x_ow1.x * mdo.y + dist21.x * mdb.y + dist31.x * mdc.y);
321 sm_threadVirial[2 * blockDim.x + threadIdx.x] =
322 -(x_ow1.x * mdo.z + dist21.x * mdb.z + dist31.x * mdc.z);
323 sm_threadVirial[3 * blockDim.x + threadIdx.x] =
324 -(x_ow1.y * mdo.y + dist21.y * mdb.y + dist31.y * mdc.y);
325 sm_threadVirial[4 * blockDim.x + threadIdx.x] =
326 -(x_ow1.y * mdo.z + dist21.y * mdb.z + dist31.y * mdc.z);
327 sm_threadVirial[5 * blockDim.x + threadIdx.x] =
328 -(x_ow1.z * mdo.z + dist21.z * mdb.z + dist31.z * mdc.z);
334 // Filling data for dummy threads with zeroes
337 for (int d = 0; d < 6; d++)
339 sm_threadVirial[d * blockDim.x + threadIdx.x] = 0.0f;
343 // Basic reduction for the values inside single thread block
344 // TODO what follows should be separated out as a standard virial reduction subroutine
347 // This is to ensure that all threads saved the data before reduction starts
349 // This casts unsigned into signed integers to avoid clang warnings
350 int tib = static_cast<int>(threadIdx.x);
351 int blockSize = static_cast<int>(blockDim.x);
352 // Reduce up to one virial per thread block
353 // All blocks are divided by half, the first half of threads sums
354 // two virials. Then the first half is divided by two and the first half
355 // of it sums two values... The procedure continues until only one thread left.
356 // Only works if the threads per blocks is a power of two.
357 for (int divideBy = 2; divideBy <= blockSize; divideBy *= 2)
359 int dividedAt = blockSize / divideBy;
362 for (int d = 0; d < 6; d++)
364 sm_threadVirial[d * blockSize + tib] +=
365 sm_threadVirial[d * blockSize + (tib + dividedAt)];
368 if (dividedAt > warpSize / 2)
373 // First 6 threads in the block add the 6 components of virial to the global memory address
376 atomicAdd(&(gm_virialScaled[tib]), sm_threadVirial[tib * blockSize]);
383 /*! \brief Select templated kernel.
385 * Returns pointer to a CUDA kernel based on provided booleans.
387 * \param[in] updateVelocities If the velocities should be constrained.
388 * \param[in] bCalcVir If virial should be updated.
390 * \retrun Pointer to CUDA kernel
392 inline auto getSettleKernelPtr(const bool updateVelocities, const bool computeVirial)
395 auto kernelPtr = settle_kernel<true, true>;
396 if (updateVelocities && computeVirial)
398 kernelPtr = settle_kernel<true, true>;
400 else if (updateVelocities && !computeVirial)
402 kernelPtr = settle_kernel<true, false>;
404 else if (!updateVelocities && computeVirial)
406 kernelPtr = settle_kernel<false, true>;
408 else if (!updateVelocities && !computeVirial)
410 kernelPtr = settle_kernel<false, false>;
415 void SettleGpu::apply(const float3* d_x,
417 const bool updateVelocities,
420 const bool computeVirial,
422 const PbcAiuc pbcAiuc)
425 ensureNoPendingCudaError("In CUDA version SETTLE");
427 // Early exit if no settles
428 if (numSettles_ == 0)
435 // Fill with zeros so the values can be reduced to it
436 // Only 6 values are needed because virial is symmetrical
437 clearDeviceBufferAsync(&d_virialScaled_, 0, 6, commandStream_);
440 auto kernelPtr = getSettleKernelPtr(updateVelocities, computeVirial);
442 KernelLaunchConfig config;
443 config.blockSize[0] = c_threadsPerBlock;
444 config.blockSize[1] = 1;
445 config.blockSize[2] = 1;
446 config.gridSize[0] = (numSettles_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
447 config.gridSize[1] = 1;
448 config.gridSize[2] = 1;
449 // Shared memory is only used for virial reduction
452 config.sharedMemorySize = c_threadsPerBlock * 6 * sizeof(float);
456 config.sharedMemorySize = 0;
458 config.stream = commandStream_;
460 const auto kernelArgs = prepareGpuKernelArguments(kernelPtr, config, &numSettles_, &d_atomIds_,
461 &settleParameters_, &d_x, &d_xp, &invdt, &d_v,
462 &d_virialScaled_, &pbcAiuc);
464 launchGpuKernel(kernelPtr, config, nullptr, "settle_kernel<updateVelocities, computeVirial>", kernelArgs);
468 copyFromDeviceBuffer(h_virialScaled_.data(), &d_virialScaled_, 0, 6, commandStream_,
469 GpuApiCallBehavior::Sync, nullptr);
471 // Mapping [XX, XY, XZ, YY, YZ, ZZ] internal format to a tensor object
472 virialScaled[XX][XX] += h_virialScaled_[0];
473 virialScaled[XX][YY] += h_virialScaled_[1];
474 virialScaled[XX][ZZ] += h_virialScaled_[2];
476 virialScaled[YY][XX] += h_virialScaled_[1];
477 virialScaled[YY][YY] += h_virialScaled_[3];
478 virialScaled[YY][ZZ] += h_virialScaled_[4];
480 virialScaled[ZZ][XX] += h_virialScaled_[2];
481 virialScaled[ZZ][YY] += h_virialScaled_[4];
482 virialScaled[ZZ][ZZ] += h_virialScaled_[5];
488 SettleGpu::SettleGpu(const gmx_mtop_t& mtop, CommandStream commandStream) :
489 commandStream_(commandStream)
491 static_assert(sizeof(real) == sizeof(float),
492 "Real numbers should be in single precision in GPU code.");
494 c_threadsPerBlock > 0 && ((c_threadsPerBlock & (c_threadsPerBlock - 1)) == 0),
495 "Number of threads per block should be a power of two in order for reduction to work.");
497 // This is to prevent the assertion failure for the systems without water
498 int totalSettles = 0;
499 for (unsigned mt = 0; mt < mtop.moltype.size(); mt++)
501 const int nral1 = 1 + NRAL(F_SETTLE);
502 InteractionList interactionList = mtop.moltype.at(mt).ilist[F_SETTLE];
503 std::vector<int> iatoms = interactionList.iatoms;
504 totalSettles += iatoms.size() / nral1;
506 if (totalSettles == 0)
511 // TODO This should be lifted to a separate subroutine that gets the values of Oxygen and
512 // Hydrogen masses, checks if they are consistent across the topology and if there is no more
513 // than two values for each mass if the free energy perturbation is enabled. In later case,
514 // masses may need to be updated on a regular basis (i.e. in set(...) method).
515 // TODO Do the checks for FEP
519 for (unsigned mt = 0; mt < mtop.moltype.size(); mt++)
521 const int nral1 = 1 + NRAL(F_SETTLE);
522 InteractionList interactionList = mtop.moltype.at(mt).ilist[F_SETTLE];
523 std::vector<int> iatoms = interactionList.iatoms;
524 for (unsigned i = 0; i < iatoms.size() / nral1; i++)
527 settler.x = iatoms[i * nral1 + 1]; // Oxygen index
528 settler.y = iatoms[i * nral1 + 2]; // First hydrogen index
529 settler.z = iatoms[i * nral1 + 3]; // Second hydrogen index
530 t_atom ow1 = mtop.moltype.at(mt).atoms.atom[settler.x];
531 t_atom hw2 = mtop.moltype.at(mt).atoms.atom[settler.y];
532 t_atom hw3 = mtop.moltype.at(mt).atoms.atom[settler.z];
542 GMX_RELEASE_ASSERT(mO == ow1.m,
543 "Topology has different values for oxygen mass. Should be identical "
544 "in order to use SETTLE.");
545 GMX_RELEASE_ASSERT(hw2.m == hw3.m && hw2.m == mH,
546 "Topology has different values for hydrogen mass. Should be "
547 "identical in order to use SETTLE.");
551 GMX_RELEASE_ASSERT(mO > 0, "Could not find oxygen mass in the topology. Needed in SETTLE.");
552 GMX_RELEASE_ASSERT(mH > 0, "Could not find hydrogen mass in the topology. Needed in SETTLE.");
554 // TODO Very similar to SETTLE initialization on CPU. Should be lifted to a separate method
555 // (one that gets dOH and dHH values and checks them for consistency)
556 int settle_type = -1;
557 for (unsigned mt = 0; mt < mtop.moltype.size(); mt++)
559 const int nral1 = 1 + NRAL(F_SETTLE);
560 InteractionList interactionList = mtop.moltype.at(mt).ilist[F_SETTLE];
561 for (int i = 0; i < interactionList.size(); i += nral1)
563 if (settle_type == -1)
565 settle_type = interactionList.iatoms[i];
567 else if (interactionList.iatoms[i] != settle_type)
570 "The [molecules] section of your topology specifies more than one block "
572 "a [moleculetype] with a [settles] block. Only one such is allowed.\n"
573 "If you are trying to partition your solvent into different *groups*\n"
574 "(e.g. for freezing, T-coupling, etc.), you are using the wrong "
576 "files specify groups. Otherwise, you may wish to change the least-used\n"
577 "block of molecules with SETTLE constraints into 3 normal constraints.");
582 GMX_RELEASE_ASSERT(settle_type >= 0, "settle_init called without settles");
584 real dOH = mtop.ffparams.iparams[settle_type].settle.doh;
585 real dHH = mtop.ffparams.iparams[settle_type].settle.dhh;
587 initSettleParameters(&settleParameters_, mO, mH, dOH, dHH);
589 allocateDeviceBuffer(&d_virialScaled_, 6, nullptr);
590 h_virialScaled_.resize(6);
593 SettleGpu::~SettleGpu()
595 // Early exit if there is no settles
596 if (numSettles_ == 0)
600 freeDeviceBuffer(&d_virialScaled_);
601 if (numAtomIdsAlloc_ > 0)
603 freeDeviceBuffer(&d_atomIds_);
607 void SettleGpu::set(const t_idef& idef, const t_mdatoms gmx_unused& md)
609 const int nral1 = 1 + NRAL(F_SETTLE);
610 t_ilist il_settle = idef.il[F_SETTLE];
611 t_iatom* iatoms = il_settle.iatoms;
612 numSettles_ = il_settle.nr / nral1;
614 reallocateDeviceBuffer(&d_atomIds_, numSettles_, &numAtomIds_, &numAtomIdsAlloc_, nullptr);
615 h_atomIds_.resize(numSettles_);
616 for (int i = 0; i < numSettles_; i++)
619 settler.x = iatoms[i * nral1 + 1]; // Oxygen index
620 settler.y = iatoms[i * nral1 + 2]; // First hydrogen index
621 settler.z = iatoms[i * nral1 + 3]; // Second hydrogen index
622 h_atomIds_.at(i) = settler;
624 copyToDeviceBuffer(&d_atomIds_, h_atomIds_.data(), 0, numSettles_, commandStream_,
625 GpuApiCallBehavior::Sync, nullptr);