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