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