Add gmx::isPowerOfTwo function
[alexxy/gromacs.git] / src / gromacs / mdlib / settle_gpu.cu
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Implements SETTLE using CUDA
38  *
39  * This file contains implementation of SETTLE constraints algorithm
40  * using CUDA, including class initialization, data-structures management
41  * and GPU kernel.
42  *
43  *
44  * \author Artem Zhmurov <zhmurov@gmail.com>
45  *
46  * \ingroup module_mdlib
47  */
48 #include "gmxpre.h"
49
50 #include "settle_gpu.cuh"
51
52 #include <assert.h>
53 #include <stdio.h>
54
55 #include <cmath>
56
57 #include <algorithm>
58
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/functions.h"
65 #include "gromacs/math/vec.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/pbcutil/pbc_aiuc_cuda.cuh"
68
69 namespace gmx
70 {
71
72 //! Number of CUDA threads in a block
73 constexpr static int c_threadsPerBlock = 256;
74 //! Maximum number of threads in a block (for __launch_bounds__)
75 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
76
77 /*! \brief SETTLE constraints kernel
78  *
79  * Each thread corresponds to a single constraints triangle (i.e. single water molecule).
80  *
81  * See original CPU version in settle.cpp
82  *
83  * \param [in]      numSettles       Number of constraints triangles (water molecules).
84  * \param [in]      gm_settles       Indexes of three atoms in the constraints triangle. The field .x of int3
85  *                                   data type corresponds to Oxygen, fields .y and .z are two hydrogen atoms.
86  * \param [in]      pars             Parameters for the algorithm (i.e. masses, target distances, etc.).
87  * \param [in]      gm_x             Coordinates of atoms before the timestep.
88  * \param [in,out]  gm_x             Coordinates of atoms after the timestep (constrained coordinates will be
89  *                                   saved here).
90  * \param [in]      invdt            Reciprocal timestep.
91  * \param [in]      gm_v             Velocities of the particles.
92  * \param [in]      gm_virialScaled  Virial tensor.
93  * \param [in]      pbcAiuc          Periodic boundary conditions data.
94  */
95 template<bool updateVelocities, bool computeVirial>
96 __launch_bounds__(c_maxThreadsPerBlock) __global__
97         void settle_kernel(const int numSettles,
98                            const int3* __restrict__ gm_settles,
99                            const SettleParameters pars,
100                            const float3* __restrict__ gm_x,
101                            float3* __restrict__ gm_xprime,
102                            float invdt,
103                            float3* __restrict__ gm_v,
104                            float* __restrict__ gm_virialScaled,
105                            const PbcAiuc pbcAiuc)
106 {
107     /* ******************************************************************* */
108     /*                                                                  ** */
109     /*    Original code by Shuichi Miyamoto, last update Oct. 1, 1992   ** */
110     /*                                                                  ** */
111     /*    Algorithm changes by Berk Hess:                               ** */
112     /*    2004-07-15 Convert COM to double precision to avoid drift     ** */
113     /*    2006-10-16 Changed velocity update to use differences         ** */
114     /*    2012-09-24 Use oxygen as reference instead of COM             ** */
115     /*    2016-02    Complete rewrite of the code for SIMD              ** */
116     /*    2020-06    Completely remove use of COM to minimize drift     ** */
117     /*                                                                  ** */
118     /*    Reference for the SETTLE algorithm                            ** */
119     /*           S. Miyamoto et al., J. Comp. Chem., 13, 952 (1992).    ** */
120     /*                                                                  ** */
121     /* ******************************************************************* */
122
123     constexpr float almost_zero = real(1e-12);
124
125     extern __shared__ float sm_threadVirial[];
126
127     int tid = static_cast<int>(blockIdx.x * blockDim.x + threadIdx.x);
128
129     if (tid < numSettles)
130     {
131         // These are the indexes of three atoms in a single 'water' molecule.
132         // TODO Can be reduced to one integer if atoms are consecutive in memory.
133         int3 indices = gm_settles[tid];
134
135         float3 x_ow1 = gm_x[indices.x];
136         float3 x_hw2 = gm_x[indices.y];
137         float3 x_hw3 = gm_x[indices.z];
138
139         float3 xprime_ow1 = gm_xprime[indices.x];
140         float3 xprime_hw2 = gm_xprime[indices.y];
141         float3 xprime_hw3 = gm_xprime[indices.z];
142
143         float3 dist21 = pbcDxAiuc(pbcAiuc, x_hw2, x_ow1);
144         float3 dist31 = pbcDxAiuc(pbcAiuc, x_hw3, x_ow1);
145         float3 doh2   = pbcDxAiuc(pbcAiuc, xprime_hw2, xprime_ow1);
146
147         float3 doh3 = pbcDxAiuc(pbcAiuc, xprime_hw3, xprime_ow1);
148
149         float3 a1 = (-doh2 - doh3) * pars.wh;
150
151         float3 b1 = doh2 + a1;
152
153         float3 c1 = doh3 + a1;
154
155         float xakszd = dist21.y * dist31.z - dist21.z * dist31.y;
156         float yakszd = dist21.z * dist31.x - dist21.x * dist31.z;
157         float zakszd = dist21.x * dist31.y - dist21.y * dist31.x;
158
159         float xaksxd = a1.y * zakszd - a1.z * yakszd;
160         float yaksxd = a1.z * xakszd - a1.x * zakszd;
161         float zaksxd = a1.x * yakszd - a1.y * xakszd;
162
163         float xaksyd = yakszd * zaksxd - zakszd * yaksxd;
164         float yaksyd = zakszd * xaksxd - xakszd * zaksxd;
165         float zaksyd = xakszd * yaksxd - yakszd * xaksxd;
166
167         float axlng = rsqrt(xaksxd * xaksxd + yaksxd * yaksxd + zaksxd * zaksxd);
168         float aylng = rsqrt(xaksyd * xaksyd + yaksyd * yaksyd + zaksyd * zaksyd);
169         float azlng = rsqrt(xakszd * xakszd + yakszd * yakszd + zakszd * zakszd);
170
171         // TODO {1,2,3} indexes should be swapped with {.x, .y, .z} components.
172         //      This way, we will be able to use vector ops more.
173         float3 trns1, trns2, trns3;
174
175         trns1.x = xaksxd * axlng;
176         trns2.x = yaksxd * axlng;
177         trns3.x = zaksxd * axlng;
178
179         trns1.y = xaksyd * aylng;
180         trns2.y = yaksyd * aylng;
181         trns3.y = zaksyd * aylng;
182
183         trns1.z = xakszd * azlng;
184         trns2.z = yakszd * azlng;
185         trns3.z = zakszd * azlng;
186
187
188         float2 b0d, c0d;
189
190         b0d.x = trns1.x * dist21.x + trns2.x * dist21.y + trns3.x * dist21.z;
191         b0d.y = trns1.y * dist21.x + trns2.y * dist21.y + trns3.y * dist21.z;
192
193         c0d.x = trns1.x * dist31.x + trns2.x * dist31.y + trns3.x * dist31.z;
194         c0d.y = trns1.y * dist31.x + trns2.y * dist31.y + trns3.y * dist31.z;
195
196         float3 b1d, c1d;
197
198         float a1d_z = trns1.z * a1.x + trns2.z * a1.y + trns3.z * a1.z;
199
200         b1d.x = trns1.x * b1.x + trns2.x * b1.y + trns3.x * b1.z;
201         b1d.y = trns1.y * b1.x + trns2.y * b1.y + trns3.y * b1.z;
202         b1d.z = trns1.z * b1.x + trns2.z * b1.y + trns3.z * b1.z;
203
204         c1d.x = trns1.x * c1.x + trns2.x * c1.y + trns3.x * c1.z;
205         c1d.y = trns1.y * c1.x + trns2.y * c1.y + trns3.y * c1.z;
206         c1d.z = trns1.z * c1.x + trns2.z * c1.y + trns3.z * c1.z;
207
208
209         float sinphi = a1d_z * rsqrt(pars.ra * pars.ra);
210         float tmp2   = 1.0f - sinphi * sinphi;
211
212         if (almost_zero > tmp2)
213         {
214             tmp2 = almost_zero;
215         }
216
217         float tmp    = rsqrt(tmp2);
218         float cosphi = tmp2 * tmp;
219         float sinpsi = (b1d.z - c1d.z) * pars.irc2 * tmp;
220         tmp2         = 1.0f - sinpsi * sinpsi;
221
222         float cospsi = tmp2 * rsqrt(tmp2);
223
224         float a2d_y = pars.ra * cosphi;
225         float b2d_x = -pars.rc * cospsi;
226         float t1    = -pars.rb * cosphi;
227         float t2    = pars.rc * sinpsi * sinphi;
228         float b2d_y = t1 - t2;
229         float c2d_y = t1 + t2;
230
231         /*     --- Step3  al,be,ga            --- */
232         float alpha  = b2d_x * (b0d.x - c0d.x) + b0d.y * b2d_y + c0d.y * c2d_y;
233         float beta   = b2d_x * (c0d.y - b0d.y) + b0d.x * b2d_y + c0d.x * c2d_y;
234         float gamma  = b0d.x * b1d.y - b1d.x * b0d.y + c0d.x * c1d.y - c1d.x * c0d.y;
235         float al2be2 = alpha * alpha + beta * beta;
236         tmp2         = (al2be2 - gamma * gamma);
237         float sinthe = (alpha * gamma - beta * tmp2 * rsqrt(tmp2)) * rsqrt(al2be2 * al2be2);
238
239         /*  --- Step4  A3' --- */
240         tmp2         = 1.0f - sinthe * sinthe;
241         float costhe = tmp2 * rsqrt(tmp2);
242
243         float3 a3d, b3d, c3d;
244
245         a3d.x = -a2d_y * sinthe;
246         a3d.y = a2d_y * costhe;
247         a3d.z = a1d_z;
248         b3d.x = b2d_x * costhe - b2d_y * sinthe;
249         b3d.y = b2d_x * sinthe + b2d_y * costhe;
250         b3d.z = b1d.z;
251         c3d.x = -b2d_x * costhe - c2d_y * sinthe;
252         c3d.y = -b2d_x * sinthe + c2d_y * costhe;
253         c3d.z = c1d.z;
254
255         /*    --- Step5  A3 --- */
256         float3 a3, b3, c3;
257
258         a3.x = trns1.x * a3d.x + trns1.y * a3d.y + trns1.z * a3d.z;
259         a3.y = trns2.x * a3d.x + trns2.y * a3d.y + trns2.z * a3d.z;
260         a3.z = trns3.x * a3d.x + trns3.y * a3d.y + trns3.z * a3d.z;
261
262         b3.x = trns1.x * b3d.x + trns1.y * b3d.y + trns1.z * b3d.z;
263         b3.y = trns2.x * b3d.x + trns2.y * b3d.y + trns2.z * b3d.z;
264         b3.z = trns3.x * b3d.x + trns3.y * b3d.y + trns3.z * b3d.z;
265
266         c3.x = trns1.x * c3d.x + trns1.y * c3d.y + trns1.z * c3d.z;
267         c3.y = trns2.x * c3d.x + trns2.y * c3d.y + trns2.z * c3d.z;
268         c3.z = trns3.x * c3d.x + trns3.y * c3d.y + trns3.z * c3d.z;
269
270
271         /* Compute and store the corrected new coordinate */
272         const float3 dxOw1 = a3 - a1;
273         const float3 dxHw2 = b3 - b1;
274         const float3 dxHw3 = c3 - c1;
275
276         gm_xprime[indices.x] = xprime_ow1 + dxOw1;
277         gm_xprime[indices.y] = xprime_hw2 + dxHw2;
278         gm_xprime[indices.z] = xprime_hw3 + dxHw3;
279
280         if (updateVelocities)
281         {
282             float3 v_ow1 = gm_v[indices.x];
283             float3 v_hw2 = gm_v[indices.y];
284             float3 v_hw3 = gm_v[indices.z];
285
286             /* Add the position correction divided by dt to the velocity */
287             v_ow1 = dxOw1 * invdt + v_ow1;
288             v_hw2 = dxHw2 * invdt + v_hw2;
289             v_hw3 = dxHw3 * invdt + v_hw3;
290
291             gm_v[indices.x] = v_ow1;
292             gm_v[indices.y] = v_hw2;
293             gm_v[indices.z] = v_hw3;
294         }
295
296         if (computeVirial)
297         {
298             float3 mdb = pars.mH * dxHw2;
299             float3 mdc = pars.mH * dxHw3;
300             float3 mdo = pars.mO * dxOw1 + mdb + mdc;
301
302             sm_threadVirial[0 * blockDim.x + threadIdx.x] =
303                     -(x_ow1.x * mdo.x + dist21.x * mdb.x + dist31.x * mdc.x);
304             sm_threadVirial[1 * blockDim.x + threadIdx.x] =
305                     -(x_ow1.x * mdo.y + dist21.x * mdb.y + dist31.x * mdc.y);
306             sm_threadVirial[2 * blockDim.x + threadIdx.x] =
307                     -(x_ow1.x * mdo.z + dist21.x * mdb.z + dist31.x * mdc.z);
308             sm_threadVirial[3 * blockDim.x + threadIdx.x] =
309                     -(x_ow1.y * mdo.y + dist21.y * mdb.y + dist31.y * mdc.y);
310             sm_threadVirial[4 * blockDim.x + threadIdx.x] =
311                     -(x_ow1.y * mdo.z + dist21.y * mdb.z + dist31.y * mdc.z);
312             sm_threadVirial[5 * blockDim.x + threadIdx.x] =
313                     -(x_ow1.z * mdo.z + dist21.z * mdb.z + dist31.z * mdc.z);
314         }
315     }
316     else
317     {
318         // Filling data for dummy threads with zeroes
319         if (computeVirial)
320         {
321             for (int d = 0; d < 6; d++)
322             {
323                 sm_threadVirial[d * blockDim.x + threadIdx.x] = 0.0f;
324             }
325         }
326     }
327     // Basic reduction for the values inside single thread block
328     // TODO what follows should be separated out as a standard virial reduction subroutine
329     if (computeVirial)
330     {
331         // This is to ensure that all threads saved the data before reduction starts
332         __syncthreads();
333         // This casts unsigned into signed integers to avoid clang warnings
334         int tib       = static_cast<int>(threadIdx.x);
335         int blockSize = static_cast<int>(blockDim.x);
336         // Reduce up to one virial per thread block
337         // All blocks are divided by half, the first half of threads sums
338         // two virials. Then the first half is divided by two and the first half
339         // of it sums two values... The procedure continues until only one thread left.
340         // Only works if the threads per blocks is a power of two.
341         for (int divideBy = 2; divideBy <= blockSize; divideBy *= 2)
342         {
343             int dividedAt = blockSize / divideBy;
344             if (tib < dividedAt)
345             {
346                 for (int d = 0; d < 6; d++)
347                 {
348                     sm_threadVirial[d * blockSize + tib] +=
349                             sm_threadVirial[d * blockSize + (tib + dividedAt)];
350                 }
351             }
352             if (dividedAt > warpSize / 2)
353             {
354                 __syncthreads();
355             }
356         }
357         // First 6 threads in the block add the 6 components of virial to the global memory address
358         if (tib < 6)
359         {
360             atomicAdd(&(gm_virialScaled[tib]), sm_threadVirial[tib * blockSize]);
361         }
362     }
363
364     return;
365 }
366
367 /*! \brief Select templated kernel.
368  *
369  * Returns pointer to a CUDA kernel based on provided booleans.
370  *
371  * \param[in] updateVelocities  If the velocities should be constrained.
372  * \param[in] bCalcVir          If virial should be updated.
373  *
374  * \retrun                      Pointer to CUDA kernel
375  */
376 inline auto getSettleKernelPtr(const bool updateVelocities, const bool computeVirial)
377 {
378
379     auto kernelPtr = settle_kernel<true, true>;
380     if (updateVelocities && computeVirial)
381     {
382         kernelPtr = settle_kernel<true, true>;
383     }
384     else if (updateVelocities && !computeVirial)
385     {
386         kernelPtr = settle_kernel<true, false>;
387     }
388     else if (!updateVelocities && computeVirial)
389     {
390         kernelPtr = settle_kernel<false, true>;
391     }
392     else if (!updateVelocities && !computeVirial)
393     {
394         kernelPtr = settle_kernel<false, false>;
395     }
396     return kernelPtr;
397 }
398
399 void SettleGpu::apply(const float3* d_x,
400                       float3*       d_xp,
401                       const bool    updateVelocities,
402                       float3*       d_v,
403                       const real    invdt,
404                       const bool    computeVirial,
405                       tensor        virialScaled,
406                       const PbcAiuc pbcAiuc)
407 {
408
409     ensureNoPendingDeviceError("In CUDA version SETTLE");
410
411     // Early exit if no settles
412     if (numSettles_ == 0)
413     {
414         return;
415     }
416
417     if (computeVirial)
418     {
419         // Fill with zeros so the values can be reduced to it
420         // Only 6 values are needed because virial is symmetrical
421         clearDeviceBufferAsync(&d_virialScaled_, 0, 6, deviceStream_);
422     }
423
424     auto kernelPtr = getSettleKernelPtr(updateVelocities, computeVirial);
425
426     KernelLaunchConfig config;
427     config.blockSize[0] = c_threadsPerBlock;
428     config.blockSize[1] = 1;
429     config.blockSize[2] = 1;
430     config.gridSize[0]  = (numSettles_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
431     config.gridSize[1]  = 1;
432     config.gridSize[2]  = 1;
433     // Shared memory is only used for virial reduction
434     if (computeVirial)
435     {
436         config.sharedMemorySize = c_threadsPerBlock * 6 * sizeof(float);
437     }
438     else
439     {
440         config.sharedMemorySize = 0;
441     }
442
443     const auto kernelArgs = prepareGpuKernelArguments(
444             kernelPtr, config, &numSettles_, &d_atomIds_, &settleParameters_, &d_x, &d_xp, &invdt, &d_v, &d_virialScaled_, &pbcAiuc);
445
446     launchGpuKernel(kernelPtr,
447                     config,
448                     deviceStream_,
449                     nullptr,
450                     "settle_kernel<updateVelocities, computeVirial>",
451                     kernelArgs);
452
453     if (computeVirial)
454     {
455         copyFromDeviceBuffer(
456                 h_virialScaled_.data(), &d_virialScaled_, 0, 6, deviceStream_, GpuApiCallBehavior::Sync, nullptr);
457
458         // Mapping [XX, XY, XZ, YY, YZ, ZZ] internal format to a tensor object
459         virialScaled[XX][XX] += h_virialScaled_[0];
460         virialScaled[XX][YY] += h_virialScaled_[1];
461         virialScaled[XX][ZZ] += h_virialScaled_[2];
462
463         virialScaled[YY][XX] += h_virialScaled_[1];
464         virialScaled[YY][YY] += h_virialScaled_[3];
465         virialScaled[YY][ZZ] += h_virialScaled_[4];
466
467         virialScaled[ZZ][XX] += h_virialScaled_[2];
468         virialScaled[ZZ][YY] += h_virialScaled_[4];
469         virialScaled[ZZ][ZZ] += h_virialScaled_[5];
470     }
471
472     return;
473 }
474
475 SettleGpu::SettleGpu(const gmx_mtop_t& mtop, const DeviceContext& deviceContext, const DeviceStream& deviceStream) :
476     deviceContext_(deviceContext),
477     deviceStream_(deviceStream)
478 {
479     static_assert(sizeof(real) == sizeof(float),
480                   "Real numbers should be in single precision in GPU code.");
481     static_assert(
482             gmx::isPowerOfTwo(c_threadsPerBlock),
483             "Number of threads per block should be a power of two in order for reduction to work.");
484
485     // This is to prevent the assertion failure for the systems without water
486     int totalSettles = 0;
487     for (unsigned mt = 0; mt < mtop.moltype.size(); mt++)
488     {
489         const int        nral1           = 1 + NRAL(F_SETTLE);
490         InteractionList  interactionList = mtop.moltype.at(mt).ilist[F_SETTLE];
491         std::vector<int> iatoms          = interactionList.iatoms;
492         totalSettles += iatoms.size() / nral1;
493     }
494     if (totalSettles == 0)
495     {
496         return;
497     }
498
499     // TODO This should be lifted to a separate subroutine that gets the values of Oxygen and
500     // Hydrogen masses, checks if they are consistent across the topology and if there is no more
501     // than two values for each mass if the free energy perturbation is enabled. In later case,
502     // masses may need to be updated on a regular basis (i.e. in set(...) method).
503     // TODO Do the checks for FEP
504     real mO = -1.0;
505     real mH = -1.0;
506
507     for (unsigned mt = 0; mt < mtop.moltype.size(); mt++)
508     {
509         const int        nral1           = 1 + NRAL(F_SETTLE);
510         InteractionList  interactionList = mtop.moltype.at(mt).ilist[F_SETTLE];
511         std::vector<int> iatoms          = interactionList.iatoms;
512         for (unsigned i = 0; i < iatoms.size() / nral1; i++)
513         {
514             int3 settler;
515             settler.x  = iatoms[i * nral1 + 1]; // Oxygen index
516             settler.y  = iatoms[i * nral1 + 2]; // First hydrogen index
517             settler.z  = iatoms[i * nral1 + 3]; // Second hydrogen index
518             t_atom ow1 = mtop.moltype.at(mt).atoms.atom[settler.x];
519             t_atom hw2 = mtop.moltype.at(mt).atoms.atom[settler.y];
520             t_atom hw3 = mtop.moltype.at(mt).atoms.atom[settler.z];
521
522             if (mO < 0)
523             {
524                 mO = ow1.m;
525             }
526             if (mH < 0)
527             {
528                 mH = hw2.m;
529             }
530             GMX_RELEASE_ASSERT(mO == ow1.m,
531                                "Topology has different values for oxygen mass. Should be identical "
532                                "in order to use SETTLE.");
533             GMX_RELEASE_ASSERT(hw2.m == hw3.m && hw2.m == mH,
534                                "Topology has different values for hydrogen mass. Should be "
535                                "identical in order to use SETTLE.");
536         }
537     }
538
539     GMX_RELEASE_ASSERT(mO > 0, "Could not find oxygen mass in the topology. Needed in SETTLE.");
540     GMX_RELEASE_ASSERT(mH > 0, "Could not find hydrogen mass in the topology. Needed in SETTLE.");
541
542     // TODO Very similar to SETTLE initialization on CPU. Should be lifted to a separate method
543     // (one that gets dOH and dHH values and checks them for consistency)
544     int settle_type = -1;
545     for (unsigned mt = 0; mt < mtop.moltype.size(); mt++)
546     {
547         const int       nral1           = 1 + NRAL(F_SETTLE);
548         InteractionList interactionList = mtop.moltype.at(mt).ilist[F_SETTLE];
549         for (int i = 0; i < interactionList.size(); i += nral1)
550         {
551             if (settle_type == -1)
552             {
553                 settle_type = interactionList.iatoms[i];
554             }
555             else if (interactionList.iatoms[i] != settle_type)
556             {
557                 gmx_fatal(FARGS,
558                           "The [molecules] section of your topology specifies more than one block "
559                           "of\n"
560                           "a [moleculetype] with a [settles] block. Only one such is allowed.\n"
561                           "If you are trying to partition your solvent into different *groups*\n"
562                           "(e.g. for freezing, T-coupling, etc.), you are using the wrong "
563                           "approach. Index\n"
564                           "files specify groups. Otherwise, you may wish to change the least-used\n"
565                           "block of molecules with SETTLE constraints into 3 normal constraints.");
566             }
567         }
568     }
569
570     GMX_RELEASE_ASSERT(settle_type >= 0, "settle_init called without settles");
571
572     real dOH = mtop.ffparams.iparams[settle_type].settle.doh;
573     real dHH = mtop.ffparams.iparams[settle_type].settle.dhh;
574
575     settleParameters_ = settleParameters(mO, mH, 1.0 / mO, 1.0 / mH, dOH, dHH);
576
577     allocateDeviceBuffer(&d_virialScaled_, 6, deviceContext_);
578     h_virialScaled_.resize(6);
579 }
580
581 SettleGpu::~SettleGpu()
582 {
583     // Early exit if there is no settles
584     if (numSettles_ == 0)
585     {
586         return;
587     }
588     freeDeviceBuffer(&d_virialScaled_);
589     if (numAtomIdsAlloc_ > 0)
590     {
591         freeDeviceBuffer(&d_atomIds_);
592     }
593 }
594
595 void SettleGpu::set(const InteractionDefinitions& idef)
596 {
597     const int              nral1     = 1 + NRAL(F_SETTLE);
598     const InteractionList& il_settle = idef.il[F_SETTLE];
599     ArrayRef<const int>    iatoms    = il_settle.iatoms;
600     numSettles_                      = il_settle.size() / nral1;
601
602     reallocateDeviceBuffer(&d_atomIds_, numSettles_, &numAtomIds_, &numAtomIdsAlloc_, deviceContext_);
603     h_atomIds_.resize(numSettles_);
604     for (int i = 0; i < numSettles_; i++)
605     {
606         int3 settler;
607         settler.x        = iatoms[i * nral1 + 1]; // Oxygen index
608         settler.y        = iatoms[i * nral1 + 2]; // First hydrogen index
609         settler.z        = iatoms[i * nral1 + 3]; // Second hydrogen index
610         h_atomIds_.at(i) = settler;
611     }
612     copyToDeviceBuffer(
613             &d_atomIds_, h_atomIds_.data(), 0, numSettles_, deviceStream_, GpuApiCallBehavior::Sync, nullptr);
614 }
615
616 } // namespace gmx