c5e92b7c296c812e114f95072ecbf06957c8088a
[alexxy/gromacs.git] / src / gromacs / mdlib / lincs_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 LINCS using CUDA
38  *
39  * This file contains implementation of LINCS constraints algorithm
40  * using CUDA, including class initialization, data-structures management
41  * and GPU kernel.
42  *
43  * \author Artem Zhmurov <zhmurov@gmail.com>
44  * \author Alan Gray <alang@nvidia.com>
45  *
46  * \ingroup module_mdlib
47  */
48 #include "gmxpre.h"
49
50 #include "lincs_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.cuh"
62 #include "gromacs/gpu_utils/gputraits.cuh"
63 #include "gromacs/gpu_utils/vectype_ops.cuh"
64 #include "gromacs/math/vec.h"
65 #include "gromacs/mdlib/constr.h"
66 #include "gromacs/pbcutil/pbc.h"
67 #include "gromacs/pbcutil/pbc_aiuc_cuda.cuh"
68 #include "gromacs/topology/ifunc.h"
69 #include "gromacs/topology/topology.h"
70
71 namespace gmx
72 {
73
74 //! Number of CUDA threads in a block
75 constexpr static int c_threadsPerBlock = 256;
76 //! Maximum number of threads in a block (for __launch_bounds__)
77 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
78
79 /*! \brief Main kernel for LINCS constraints.
80  *
81  * See Hess et al., J. Comput. Chem. 18: 1463-1472 (1997) for the description of the algorithm.
82  *
83  * In CUDA version, one thread is responsible for all computations for one constraint. The blocks are
84  * filled in a way that no constraint is coupled to the constraint from the next block. This is achieved
85  * by moving active threads to the next block, if the correspondent group of coupled constraints is to big
86  * to fit the current thread block. This may leave some 'dummy' threads in the end of the thread block, i.e.
87  * threads that are not required to do actual work. Since constraints from different blocks are not coupled,
88  * there is no need to synchronize across the device. However, extensive communication in a thread block
89  * are still needed.
90  *
91  * \todo Reduce synchronization overhead. Some ideas are:
92  *        1. Consider going to warp-level synchronization for the coupled constraints.
93  *        2. Move more data to local/shared memory and try to get rid of atomic operations (at least on
94  *           the device level).
95  *        3. Use analytical solution for matrix A inversion.
96  *        4. Introduce mapping of thread id to both single constraint and single atom, thus designating
97  *           Nth threads to deal with Nat <= Nth coupled atoms and Nc <= Nth coupled constraints.
98  *       See Redmine issue #2885 for details (https://redmine.gromacs.org/issues/2885)
99  * \todo The use of __restrict__  for gm_xp and gm_v causes failure, probably because of the atomic
100          operations. Investigate this issue further.
101  *
102  * \param[in,out] kernelParams  All parameters and pointers for the kernel condensed in single struct.
103  * \param[in]     invdt         Inverse timestep (needed to update velocities).
104  */
105 template<bool updateVelocities, bool computeVirial>
106 __launch_bounds__(c_maxThreadsPerBlock) __global__
107         void lincs_kernel(LincsGpuKernelParameters kernelParams,
108                           const float3* __restrict__ gm_x,
109                           float3*     gm_xp,
110                           float3*     gm_v,
111                           const float invdt)
112 {
113     const PbcAiuc pbcAiuc                                 = kernelParams.pbcAiuc;
114     const int     numConstraintsThreads                   = kernelParams.numConstraintsThreads;
115     const int     numIterations                           = kernelParams.numIterations;
116     const int     expansionOrder                          = kernelParams.expansionOrder;
117     const int2* __restrict__ gm_constraints               = kernelParams.d_constraints;
118     const float* __restrict__ gm_constraintsTargetLengths = kernelParams.d_constraintsTargetLengths;
119     const int* __restrict__ gm_coupledConstraintsCounts   = kernelParams.d_coupledConstraintsCounts;
120     const int* __restrict__ gm_coupledConstraintsIdxes = kernelParams.d_coupledConstraintsIndices;
121     const float* __restrict__ gm_massFactors           = kernelParams.d_massFactors;
122     float* __restrict__ gm_matrixA                     = kernelParams.d_matrixA;
123     const float* __restrict__ gm_inverseMasses         = kernelParams.d_inverseMasses;
124     float* __restrict__ gm_virialScaled                = kernelParams.d_virialScaled;
125
126     int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
127
128     // numConstraintsThreads should be a integer multiple of blockSize (numConstraintsThreads = numBlocks*blockSize).
129     // This is to ensure proper synchronizations and reduction. All array are padded to the required size.
130     assert(threadIndex < numConstraintsThreads);
131
132     // Vectors connecting constrained atoms before algorithm was applied.
133     // Needed to construct constrain matrix A
134     extern __shared__ float3 sm_r[];
135
136     int2 pair = gm_constraints[threadIndex];
137     int  i    = pair.x;
138     int  j    = pair.y;
139
140     // Mass-scaled Lagrange multiplier
141     float lagrangeScaled = 0.0f;
142
143     float targetLength;
144     float inverseMassi;
145     float inverseMassj;
146     float sqrtReducedMass;
147
148     float3 xi;
149     float3 xj;
150     float3 rc;
151
152     // i == -1 indicates dummy constraint at the end of the thread block.
153     bool isDummyThread = (i == -1);
154
155     // Everything computed for these dummies will be equal to zero
156     if (isDummyThread)
157     {
158         targetLength    = 0.0f;
159         inverseMassi    = 0.0f;
160         inverseMassj    = 0.0f;
161         sqrtReducedMass = 0.0f;
162
163         xi = make_float3(0.0f, 0.0f, 0.0f);
164         xj = make_float3(0.0f, 0.0f, 0.0f);
165         rc = make_float3(0.0f, 0.0f, 0.0f);
166     }
167     else
168     {
169         // Collecting data
170         targetLength    = gm_constraintsTargetLengths[threadIndex];
171         inverseMassi    = gm_inverseMasses[i];
172         inverseMassj    = gm_inverseMasses[j];
173         sqrtReducedMass = rsqrt(inverseMassi + inverseMassj);
174
175         xi = gm_x[i];
176         xj = gm_x[j];
177
178         float3 dx = pbcDxAiuc(pbcAiuc, xi, xj);
179
180         float rlen = rsqrtf(dx.x * dx.x + dx.y * dx.y + dx.z * dx.z);
181         rc         = rlen * dx;
182     }
183
184     sm_r[threadIdx.x] = rc;
185     // Make sure that all r's are saved into shared memory
186     // before they are accessed in the loop below
187     __syncthreads();
188
189     /*
190      * Constructing LINCS matrix (A)
191      */
192
193     // Only non-zero values are saved (for coupled constraints)
194     int coupledConstraintsCount = gm_coupledConstraintsCounts[threadIndex];
195     for (int n = 0; n < coupledConstraintsCount; n++)
196     {
197         int index = n * numConstraintsThreads + threadIndex;
198         int c1    = gm_coupledConstraintsIdxes[index];
199
200         float3 rc1        = sm_r[c1 - blockIdx.x * blockDim.x];
201         gm_matrixA[index] = gm_massFactors[index] * (rc.x * rc1.x + rc.y * rc1.y + rc.z * rc1.z);
202     }
203
204     // Skipping in dummy threads
205     if (!isDummyThread)
206     {
207         xi = gm_xp[i];
208         xj = gm_xp[j];
209     }
210
211     float3 dx = pbcDxAiuc(pbcAiuc, xi, xj);
212
213     float sol = sqrtReducedMass * ((rc.x * dx.x + rc.y * dx.y + rc.z * dx.z) - targetLength);
214
215     /*
216      *  Inverse matrix using a set of expansionOrder matrix multiplications
217      */
218
219     // This will use the same memory space as sm_r, which is no longer needed.
220     extern __shared__ float sm_rhs[];
221     // Save current right-hand-side vector in the shared memory
222     sm_rhs[threadIdx.x] = sol;
223
224     for (int rec = 0; rec < expansionOrder; rec++)
225     {
226         // Making sure that all sm_rhs are saved before they are accessed in a loop below
227         __syncthreads();
228         float mvb = 0.0f;
229
230         for (int n = 0; n < coupledConstraintsCount; n++)
231         {
232             int index = n * numConstraintsThreads + threadIndex;
233             int c1    = gm_coupledConstraintsIdxes[index];
234             // Convolute current right-hand-side with A
235             // Different, non overlapping parts of sm_rhs[..] are read during odd and even iterations
236             mvb = mvb + gm_matrixA[index] * sm_rhs[c1 - blockIdx.x * blockDim.x + blockDim.x * (rec % 2)];
237         }
238         // 'Switch' rhs vectors, save current result
239         // These values will be accessed in the loop above during the next iteration.
240         sm_rhs[threadIdx.x + blockDim.x * ((rec + 1) % 2)] = mvb;
241         sol                                                = sol + mvb;
242     }
243
244     // Current mass-scaled Lagrange multipliers
245     lagrangeScaled = sqrtReducedMass * sol;
246
247     // Save updated coordinates before correction for the rotational lengthening
248     float3 tmp = rc * lagrangeScaled;
249
250     // Writing for all but dummy constraints
251     if (!isDummyThread)
252     {
253         atomicAdd(&gm_xp[i], -tmp * inverseMassi);
254         atomicAdd(&gm_xp[j], tmp * inverseMassj);
255     }
256
257     /*
258      *  Correction for centripetal effects
259      */
260     for (int iter = 0; iter < numIterations; iter++)
261     {
262         // Make sure that all xp's are saved: atomic operation calls before are
263         // communicating current xp[..] values across thread block.
264         __syncthreads();
265
266         if (!isDummyThread)
267         {
268             xi = gm_xp[i];
269             xj = gm_xp[j];
270         }
271
272         float3 dx = pbcDxAiuc(pbcAiuc, xi, xj);
273
274         float len2  = targetLength * targetLength;
275         float dlen2 = 2.0f * len2 - norm2(dx);
276
277         // TODO A little bit more effective but slightly less readable version of the below would be:
278         //      float proj = sqrtReducedMass*(targetLength - (dlen2 > 0.0f ? 1.0f : 0.0f)*dlen2*rsqrt(dlen2));
279         float proj;
280         if (dlen2 > 0.0f)
281         {
282             proj = sqrtReducedMass * (targetLength - dlen2 * rsqrt(dlen2));
283         }
284         else
285         {
286             proj = sqrtReducedMass * targetLength;
287         }
288
289         sm_rhs[threadIdx.x] = proj;
290         float sol           = proj;
291
292         /*
293          * Same matrix inversion as above is used for updated data
294          */
295         for (int rec = 0; rec < expansionOrder; rec++)
296         {
297             // Make sure that all elements of rhs are saved into shared memory
298             __syncthreads();
299             float mvb = 0;
300
301             for (int n = 0; n < coupledConstraintsCount; n++)
302             {
303                 int index = n * numConstraintsThreads + threadIndex;
304                 int c1    = gm_coupledConstraintsIdxes[index];
305
306                 mvb = mvb + gm_matrixA[index] * sm_rhs[c1 - blockIdx.x * blockDim.x + blockDim.x * (rec % 2)];
307             }
308             sm_rhs[threadIdx.x + blockDim.x * ((rec + 1) % 2)] = mvb;
309             sol                                                = sol + mvb;
310         }
311
312         // Add corrections to Lagrange multipliers
313         float sqrtmu_sol = sqrtReducedMass * sol;
314         lagrangeScaled += sqrtmu_sol;
315
316         // Save updated coordinates for the next iteration
317         // Dummy constraints are skipped
318         if (!isDummyThread)
319         {
320             float3 tmp = rc * sqrtmu_sol;
321             atomicAdd(&gm_xp[i], -tmp * inverseMassi);
322             atomicAdd(&gm_xp[j], tmp * inverseMassj);
323         }
324     }
325
326     // Updating particle velocities for all but dummy threads
327     if (updateVelocities && !isDummyThread)
328     {
329         float3 tmp = rc * invdt * lagrangeScaled;
330         atomicAdd(&gm_v[i], -tmp * inverseMassi);
331         atomicAdd(&gm_v[j], tmp * inverseMassj);
332     }
333
334
335     if (computeVirial)
336     {
337         // Virial is computed from Lagrange multiplier (lagrangeScaled), target constrain length
338         // (targetLength) and the normalized vector connecting constrained atoms before
339         // the algorithm was applied (rc). The evaluation of virial in each thread is
340         // followed by basic reduction for the values inside single thread block.
341         // Then, the values are reduced across grid by atomicAdd(...).
342         //
343         // TODO Shuffle reduction.
344         // TODO Should be unified and/or done once when virial is actually needed.
345         // TODO Recursive version that removes atomicAdd(...)'s entirely is needed. Ideally,
346         //      one that works for any datatype.
347
348         // Save virial for each thread into the shared memory. Tensor is symmetrical, hence only
349         // 6 values are saved. Dummy threads will have zeroes in their virial: targetLength,
350         // lagrangeScaled and rc are all set to zero for them in the beginning of the kernel.
351         // The sm_threadVirial[..] will overlap with the sm_r[..] and sm_rhs[..], but the latter
352         // two are no longer in use.
353         extern __shared__ float sm_threadVirial[];
354         float                   mult                  = targetLength * lagrangeScaled;
355         sm_threadVirial[0 * blockDim.x + threadIdx.x] = mult * rc.x * rc.x;
356         sm_threadVirial[1 * blockDim.x + threadIdx.x] = mult * rc.x * rc.y;
357         sm_threadVirial[2 * blockDim.x + threadIdx.x] = mult * rc.x * rc.z;
358         sm_threadVirial[3 * blockDim.x + threadIdx.x] = mult * rc.y * rc.y;
359         sm_threadVirial[4 * blockDim.x + threadIdx.x] = mult * rc.y * rc.z;
360         sm_threadVirial[5 * blockDim.x + threadIdx.x] = mult * rc.z * rc.z;
361
362         __syncthreads();
363
364         // Reduce up to one virial per thread block. All blocks are divided by half, the first
365         // half of threads sums two virials. Then the first half is divided by two and the first
366         // half of it sums two values. This procedure is repeated until only one thread is left.
367         // Only works if the threads per blocks is a power of two (hence static_assert
368         // in the beginning of the kernel).
369         for (int divideBy = 2; divideBy <= static_cast<int>(blockDim.x); divideBy *= 2)
370         {
371             int dividedAt = blockDim.x / divideBy;
372             if (static_cast<int>(threadIdx.x) < dividedAt)
373             {
374                 for (int d = 0; d < 6; d++)
375                 {
376                     sm_threadVirial[d * blockDim.x + threadIdx.x] +=
377                             sm_threadVirial[d * blockDim.x + (threadIdx.x + dividedAt)];
378                 }
379             }
380             // Syncronize if not within one warp
381             if (dividedAt > warpSize / 2)
382             {
383                 __syncthreads();
384             }
385         }
386         // First 6 threads in the block add the results of 6 tensor components to the global memory address.
387         if (threadIdx.x < 6)
388         {
389             atomicAdd(&(gm_virialScaled[threadIdx.x]), sm_threadVirial[threadIdx.x * blockDim.x]);
390         }
391     }
392
393     return;
394 }
395
396 /*! \brief Select templated kernel.
397  *
398  * Returns pointer to a CUDA kernel based on provided booleans.
399  *
400  * \param[in] updateVelocities  If the velocities should be constrained.
401  * \param[in] computeVirial     If virial should be updated.
402  *
403  * \return                      Pointer to CUDA kernel
404  */
405 inline auto getLincsKernelPtr(const bool updateVelocities, const bool computeVirial)
406 {
407
408     auto kernelPtr = lincs_kernel<true, true>;
409     if (updateVelocities && computeVirial)
410     {
411         kernelPtr = lincs_kernel<true, true>;
412     }
413     else if (updateVelocities && !computeVirial)
414     {
415         kernelPtr = lincs_kernel<true, false>;
416     }
417     else if (!updateVelocities && computeVirial)
418     {
419         kernelPtr = lincs_kernel<false, true>;
420     }
421     else if (!updateVelocities && !computeVirial)
422     {
423         kernelPtr = lincs_kernel<false, false>;
424     }
425     return kernelPtr;
426 }
427
428 void LincsGpu::apply(const float3* d_x,
429                      float3*       d_xp,
430                      const bool    updateVelocities,
431                      float3*       d_v,
432                      const real    invdt,
433                      const bool    computeVirial,
434                      tensor        virialScaled,
435                      const PbcAiuc pbcAiuc)
436 {
437     ensureNoPendingCudaError("In CUDA version of LINCS");
438
439     // Early exit if no constraints
440     if (kernelParams_.numConstraintsThreads == 0)
441     {
442         return;
443     }
444
445     if (computeVirial)
446     {
447         // Fill with zeros so the values can be reduced to it
448         // Only 6 values are needed because virial is symmetrical
449         clearDeviceBufferAsync(&kernelParams_.d_virialScaled, 0, 6, commandStream_);
450     }
451
452     auto kernelPtr = getLincsKernelPtr(updateVelocities, computeVirial);
453
454     KernelLaunchConfig config;
455     config.blockSize[0] = c_threadsPerBlock;
456     config.blockSize[1] = 1;
457     config.blockSize[2] = 1;
458     config.gridSize[0] = (kernelParams_.numConstraintsThreads + c_threadsPerBlock - 1) / c_threadsPerBlock;
459     config.gridSize[1] = 1;
460     config.gridSize[2] = 1;
461
462     // Shared memory is used to store:
463     // -- Current coordinates (3 floats per thread)
464     // -- Right-hand-sides for matrix inversion (2 floats per thread)
465     // -- Virial tensor components (6 floats per thread)
466     // Since none of these three are needed simultaneously, they can be saved at the same shared memory address
467     // (i.e. correspondent arrays are intentionally overlapped in address space). Consequently, only
468     // max{3, 2, 6} = 6 floats per thread are needed in case virial is computed, or max{3, 2} = 3 if not.
469     if (computeVirial)
470     {
471         config.sharedMemorySize = c_threadsPerBlock * 6 * sizeof(float);
472     }
473     else
474     {
475         config.sharedMemorySize = c_threadsPerBlock * 3 * sizeof(float);
476     }
477     config.stream = commandStream_;
478
479     kernelParams_.pbcAiuc = pbcAiuc;
480
481     const auto kernelArgs =
482             prepareGpuKernelArguments(kernelPtr, config, &kernelParams_, &d_x, &d_xp, &d_v, &invdt);
483
484     launchGpuKernel(kernelPtr, config, nullptr, "lincs_kernel<updateVelocities, computeVirial>", kernelArgs);
485
486     if (computeVirial)
487     {
488         // Copy LINCS virial data and add it to the common virial
489         copyFromDeviceBuffer(h_virialScaled_.data(), &kernelParams_.d_virialScaled, 0, 6,
490                              commandStream_, GpuApiCallBehavior::Sync, nullptr);
491
492         // Mapping [XX, XY, XZ, YY, YZ, ZZ] internal format to a tensor object
493         virialScaled[XX][XX] += h_virialScaled_[0];
494         virialScaled[XX][YY] += h_virialScaled_[1];
495         virialScaled[XX][ZZ] += h_virialScaled_[2];
496
497         virialScaled[YY][XX] += h_virialScaled_[1];
498         virialScaled[YY][YY] += h_virialScaled_[3];
499         virialScaled[YY][ZZ] += h_virialScaled_[4];
500
501         virialScaled[ZZ][XX] += h_virialScaled_[2];
502         virialScaled[ZZ][YY] += h_virialScaled_[4];
503         virialScaled[ZZ][ZZ] += h_virialScaled_[5];
504     }
505
506     return;
507 }
508
509 LincsGpu::LincsGpu(int numIterations, int expansionOrder, CommandStream commandStream) :
510     commandStream_(commandStream)
511 {
512     kernelParams_.numIterations  = numIterations;
513     kernelParams_.expansionOrder = expansionOrder;
514
515     static_assert(sizeof(real) == sizeof(float),
516                   "Real numbers should be in single precision in GPU code.");
517     static_assert(
518             c_threadsPerBlock > 0 && ((c_threadsPerBlock & (c_threadsPerBlock - 1)) == 0),
519             "Number of threads per block should be a power of two in order for reduction to work.");
520
521     allocateDeviceBuffer(&kernelParams_.d_virialScaled, 6, nullptr);
522     h_virialScaled_.resize(6);
523
524     // The data arrays should be expanded/reallocated on first call of set() function.
525     numConstraintsThreadsAlloc_ = 0;
526     numAtomsAlloc_              = 0;
527 }
528
529 LincsGpu::~LincsGpu()
530 {
531     freeDeviceBuffer(&kernelParams_.d_virialScaled);
532
533     if (numConstraintsThreadsAlloc_ > 0)
534     {
535         freeDeviceBuffer(&kernelParams_.d_constraints);
536         freeDeviceBuffer(&kernelParams_.d_constraintsTargetLengths);
537
538         freeDeviceBuffer(&kernelParams_.d_coupledConstraintsCounts);
539         freeDeviceBuffer(&kernelParams_.d_coupledConstraintsIndices);
540         freeDeviceBuffer(&kernelParams_.d_massFactors);
541         freeDeviceBuffer(&kernelParams_.d_matrixA);
542     }
543     if (numAtomsAlloc_ > 0)
544     {
545         freeDeviceBuffer(&kernelParams_.d_inverseMasses);
546     }
547 }
548
549 //! Helper type for discovering coupled constraints
550 struct AtomsAdjacencyListElement
551 {
552     AtomsAdjacencyListElement(const int indexOfSecondConstrainedAtom,
553                               const int indexOfConstraint,
554                               const int signFactor) :
555         indexOfSecondConstrainedAtom_(indexOfSecondConstrainedAtom),
556         indexOfConstraint_(indexOfConstraint),
557         signFactor_(signFactor)
558     {
559     }
560     //! The index of the other atom constrained to this atom.
561     int indexOfSecondConstrainedAtom_;
562     //! The index of this constraint in the container of constraints.
563     int indexOfConstraint_;
564     /*! \brief A multiplicative factor that indicates the relative
565      * order of the atoms in the atom list.
566      *
567      * Used for computing the mass factor of this constraint
568      * relative to any coupled constraints. */
569     int signFactor_;
570 };
571 //! Constructs and returns an atom constraint adjacency list
572 static std::vector<std::vector<AtomsAdjacencyListElement>>
573 constructAtomsAdjacencyList(const int numAtoms, ArrayRef<const int> iatoms)
574 {
575     const int                                           stride         = 1 + NRAL(F_CONSTR);
576     const int                                           numConstraints = iatoms.ssize() / stride;
577     std::vector<std::vector<AtomsAdjacencyListElement>> atomsAdjacencyList(numAtoms);
578     for (int c = 0; c < numConstraints; c++)
579     {
580         int a1 = iatoms[stride * c + 1];
581         int a2 = iatoms[stride * c + 2];
582
583         // Each constraint will be represented as a tuple, containing index of the second
584         // constrained atom, index of the constraint and a sign that indicates the order of atoms in
585         // which they are listed. Sign is needed to compute the mass factors.
586         atomsAdjacencyList[a1].emplace_back(a2, c, +1);
587         atomsAdjacencyList[a2].emplace_back(a1, c, -1);
588     }
589
590     return atomsAdjacencyList;
591 }
592
593 /*! \brief Helper function to go through constraints recursively.
594  *
595  *  For each constraint, counts the number of coupled constraints and stores the value in \p numCoupledConstraints array.
596  *  This information is used to split the array of constraints between thread blocks on a GPU so there is no
597  *  coupling between constraints from different thread blocks. After the \p numCoupledConstraints array is filled, the
598  *  value \p numCoupledConstraints[c] should be equal to the number of constraints that are coupled to \p c and located
599  *  after it in the constraints array.
600  *
601  * \param[in]     a                   Atom index.
602  * \param[in,out] numCoupledConstraints  Indicates if the constraint was already counted and stores
603  *                                    the number of constraints (i) connected to it and (ii) located
604  *                                    after it in memory. This array is filled by this recursive function.
605  *                                    For a set of coupled constraints, only for the first one in this list
606  *                                    the number of consecutive coupled constraints is needed: if there is
607  *                                    not enough space for this set of constraints in the thread block,
608  *                                    the group has to be moved to the next one.
609  * \param[in]     atomsAdjacencyList  Stores information about connections between atoms.
610  */
611 inline int countCoupled(int           a,
612                         ArrayRef<int> numCoupledConstraints,
613                         ArrayRef<const std::vector<AtomsAdjacencyListElement>> atomsAdjacencyList)
614
615 {
616     int counted = 0;
617     for (const auto& adjacentAtom : atomsAdjacencyList[a])
618     {
619         const int c2 = adjacentAtom.indexOfConstraint_;
620         if (numCoupledConstraints[c2] == -1)
621         {
622             numCoupledConstraints[c2] = 0; // To indicate we've been here
623             counted += 1
624                        + countCoupled(adjacentAtom.indexOfSecondConstrainedAtom_,
625                                       numCoupledConstraints, atomsAdjacencyList);
626         }
627     }
628     return counted;
629 }
630
631 /*! \brief Add constraint to \p splitMap with all constraints coupled to it.
632  *
633  *  Adds the constraint \p c from the constrain list \p iatoms to the map \p splitMap
634  *  if it was not yet added. Then goes through all the constraints coupled to \p c
635  *  and calls itself recursively. This ensures that all the coupled constraints will
636  *  be added to neighboring locations in the final data structures on the device,
637  *  hence mapping all coupled constraints to the same thread block. A value of -1 in
638  *  the \p splitMap is used to flag that constraint was not yet added to the \p splitMap.
639  *
640  * \param[in]     iatoms              The list of constraints.
641  * \param[in]     stride              Number of elements per constraint in \p iatoms.
642  * \param[in]     atomsAdjacencyList  Information about connections between atoms.
643  * \param[our]    splitMap            Map of sequential constraint indexes to indexes to be on the device
644  * \param[in]     c                   Sequential index for constraint to consider adding.
645  * \param[in,out] currentMapIndex     The rolling index for the constraints mapping.
646  */
647 inline void addWithCoupled(ArrayRef<const int>                                    iatoms,
648                            const int                                              stride,
649                            ArrayRef<const std::vector<AtomsAdjacencyListElement>> atomsAdjacencyList,
650                            ArrayRef<int>                                          splitMap,
651                            const int                                              c,
652                            int*                                                   currentMapIndex)
653 {
654     if (splitMap[c] == -1)
655     {
656         splitMap[c] = *currentMapIndex;
657         (*currentMapIndex)++;
658
659         // Constraints, coupled through both atoms.
660         for (int atomIndexInConstraint = 0; atomIndexInConstraint < 2; atomIndexInConstraint++)
661         {
662             const int a1 = iatoms[stride * c + 1 + atomIndexInConstraint];
663             for (const auto& adjacentAtom : atomsAdjacencyList[a1])
664             {
665                 const int c2 = adjacentAtom.indexOfConstraint_;
666                 if (c2 != c)
667                 {
668                     addWithCoupled(iatoms, stride, atomsAdjacencyList, splitMap, c2, currentMapIndex);
669                 }
670             }
671         }
672     }
673 }
674
675 /*! \brief Computes and returns how many constraints are coupled to each constraint
676  *
677  * Needed to introduce splits in data so that all coupled constraints will be computed in a
678  * single GPU block. The position \p c of the vector \p numCoupledConstraints should have the number
679  * of constraints that are coupled to a constraint \p c and are after \p c in the vector. Only
680  * first index of the connected group of the constraints is needed later in the code, hence the
681  * numCoupledConstraints vector is also used to keep track if the constrain was already counted.
682  */
683 static std::vector<int> countNumCoupledConstraints(ArrayRef<const int> iatoms,
684                                                    ArrayRef<const std::vector<AtomsAdjacencyListElement>> atomsAdjacencyList)
685 {
686     const int        stride         = 1 + NRAL(F_CONSTR);
687     const int        numConstraints = iatoms.ssize() / stride;
688     std::vector<int> numCoupledConstraints(numConstraints, -1);
689     for (int c = 0; c < numConstraints; c++)
690     {
691         const int a1 = iatoms[stride * c + 1];
692         const int a2 = iatoms[stride * c + 2];
693         if (numCoupledConstraints[c] == -1)
694         {
695             numCoupledConstraints[c] = countCoupled(a1, numCoupledConstraints, atomsAdjacencyList)
696                                        + countCoupled(a2, numCoupledConstraints, atomsAdjacencyList);
697         }
698     }
699
700     return numCoupledConstraints;
701 }
702
703 bool LincsGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
704 {
705     for (const gmx_moltype_t& molType : mtop.moltype)
706     {
707         ArrayRef<const int> iatoms    = molType.ilist[F_CONSTR].iatoms;
708         const auto atomsAdjacencyList = constructAtomsAdjacencyList(molType.atoms.nr, iatoms);
709         // Compute, how many constraints are coupled to each constraint
710         const auto numCoupledConstraints = countNumCoupledConstraints(iatoms, atomsAdjacencyList);
711         for (const int numCoupled : numCoupledConstraints)
712         {
713             if (numCoupled > c_threadsPerBlock)
714             {
715                 return false;
716             }
717         }
718     }
719
720     return true;
721 }
722
723 void LincsGpu::set(const t_idef& idef, const t_mdatoms& md)
724 {
725     int numAtoms = md.nr;
726     // List of constrained atoms (CPU memory)
727     std::vector<int2> constraintsHost;
728     // Equilibrium distances for the constraints (CPU)
729     std::vector<float> constraintsTargetLengthsHost;
730     // Number of constraints, coupled with the current one (CPU)
731     std::vector<int> coupledConstraintsCountsHost;
732     // List of coupled with the current one (CPU)
733     std::vector<int> coupledConstraintsIndicesHost;
734     // Mass factors (CPU)
735     std::vector<float> massFactorsHost;
736
737     // List of constrained atoms in local topology
738     gmx::ArrayRef<const int> iatoms =
739             constArrayRefFromArray(idef.il[F_CONSTR].iatoms, idef.il[F_CONSTR].nr);
740     const int stride         = NRAL(F_CONSTR) + 1;
741     const int numConstraints = idef.il[F_CONSTR].nr / stride;
742
743     // Early exit if no constraints
744     if (numConstraints == 0)
745     {
746         kernelParams_.numConstraintsThreads = 0;
747         return;
748     }
749
750     // Construct the adjacency list, a useful intermediate structure
751     const auto atomsAdjacencyList = constructAtomsAdjacencyList(numAtoms, iatoms);
752
753     // Compute, how many constraints are coupled to each constraint
754     const auto numCoupledConstraints = countNumCoupledConstraints(iatoms, atomsAdjacencyList);
755
756     // Map of splits in the constraints data. For each 'old' constraint index gives 'new' which
757     // takes into account the empty spaces which might be needed in the end of each thread block.
758     std::vector<int> splitMap(numConstraints, -1);
759     int              currentMapIndex = 0;
760     for (int c = 0; c < numConstraints; c++)
761     {
762         // Check if coupled constraints all fit in one block
763         if (numCoupledConstraints.at(c) > c_threadsPerBlock)
764         {
765             gmx_fatal(FARGS,
766                       "Maximum number of coupled constraints (%d) exceeds the size of the CUDA "
767                       "thread block (%d). Most likely, you are trying to use the GPU version of "
768                       "LINCS with constraints on all-bonds, which is not supported for large "
769                       "molecules. When compatible with the force field and integration settings, "
770                       "using constraints on H-bonds only.",
771                       numCoupledConstraints.at(c), c_threadsPerBlock);
772         }
773         if (currentMapIndex / c_threadsPerBlock
774             != (currentMapIndex + numCoupledConstraints.at(c)) / c_threadsPerBlock)
775         {
776             currentMapIndex = ((currentMapIndex / c_threadsPerBlock) + 1) * c_threadsPerBlock;
777         }
778         addWithCoupled(iatoms, stride, atomsAdjacencyList, splitMap, c, &currentMapIndex);
779     }
780
781     kernelParams_.numConstraintsThreads =
782             currentMapIndex + c_threadsPerBlock - currentMapIndex % c_threadsPerBlock;
783     GMX_RELEASE_ASSERT(kernelParams_.numConstraintsThreads % c_threadsPerBlock == 0,
784                        "Number of threads should be a multiple of the block size");
785
786     // Initialize constraints and their target indexes taking into account the splits in the data arrays.
787     int2 pair;
788     pair.x = -1;
789     pair.y = -1;
790     constraintsHost.resize(kernelParams_.numConstraintsThreads, pair);
791     std::fill(constraintsHost.begin(), constraintsHost.end(), pair);
792     constraintsTargetLengthsHost.resize(kernelParams_.numConstraintsThreads, 0.0);
793     std::fill(constraintsTargetLengthsHost.begin(), constraintsTargetLengthsHost.end(), 0.0);
794     for (int c = 0; c < numConstraints; c++)
795     {
796         int a1   = iatoms[stride * c + 1];
797         int a2   = iatoms[stride * c + 2];
798         int type = iatoms[stride * c];
799
800         int2 pair;
801         pair.x                                          = a1;
802         pair.y                                          = a2;
803         constraintsHost.at(splitMap.at(c))              = pair;
804         constraintsTargetLengthsHost.at(splitMap.at(c)) = idef.iparams[type].constr.dA;
805     }
806
807     // The adjacency list of constraints (i.e. the list of coupled constraints for each constraint).
808     // We map a single thread to a single constraint, hence each thread 'c' will be using one
809     // element from coupledConstraintsCountsHost array, which is the number of constraints coupled
810     // to the constraint 'c'. The coupled constraints indexes are placed into the
811     // coupledConstraintsIndicesHost array. Latter is organized as a one-dimensional array to ensure
812     // good memory alignment. It is addressed as [c + i*numConstraintsThreads], where 'i' goes from
813     // zero to the number of constraints coupled to 'c'. 'numConstraintsThreads' is the width of the
814     // array --- a number, greater then total number of constraints, taking into account the splits
815     // in the constraints array due to the GPU block borders. This number can be adjusted to improve
816     // memory access pattern. Mass factors are saved in a similar data structure.
817     int maxCoupledConstraints = 0;
818     for (int c = 0; c < numConstraints; c++)
819     {
820         int a1 = iatoms[stride * c + 1];
821         int a2 = iatoms[stride * c + 2];
822
823         // Constraint 'c' is counted twice, but it should be excluded altogether. Hence '-2'.
824         int nCoupedConstraints = atomsAdjacencyList.at(a1).size() + atomsAdjacencyList.at(a2).size() - 2;
825
826         if (nCoupedConstraints > maxCoupledConstraints)
827         {
828             maxCoupledConstraints = nCoupedConstraints;
829         }
830     }
831
832     coupledConstraintsCountsHost.resize(kernelParams_.numConstraintsThreads, 0);
833     coupledConstraintsIndicesHost.resize(maxCoupledConstraints * kernelParams_.numConstraintsThreads, -1);
834     massFactorsHost.resize(maxCoupledConstraints * kernelParams_.numConstraintsThreads, -1);
835
836     for (int c1 = 0; c1 < numConstraints; c1++)
837     {
838         coupledConstraintsCountsHost.at(splitMap.at(c1)) = 0;
839         int c1a1                                         = iatoms[stride * c1 + 1];
840         int c1a2                                         = iatoms[stride * c1 + 2];
841
842         // Constraints, coupled through the first atom.
843         int c2a1 = c1a1;
844         for (const auto& atomAdjacencyList : atomsAdjacencyList[c1a1])
845         {
846             int c2 = atomAdjacencyList.indexOfConstraint_;
847
848             if (c1 != c2)
849             {
850                 int c2a2  = atomAdjacencyList.indexOfSecondConstrainedAtom_;
851                 int sign  = atomAdjacencyList.signFactor_;
852                 int index = kernelParams_.numConstraintsThreads
853                                     * coupledConstraintsCountsHost.at(splitMap.at(c1))
854                             + splitMap.at(c1);
855
856                 coupledConstraintsIndicesHost.at(index) = splitMap.at(c2);
857
858                 int center = c1a1;
859
860                 float sqrtmu1 = 1.0 / sqrt(md.invmass[c1a1] + md.invmass[c1a2]);
861                 float sqrtmu2 = 1.0 / sqrt(md.invmass[c2a1] + md.invmass[c2a2]);
862
863                 massFactorsHost.at(index) = -sign * md.invmass[center] * sqrtmu1 * sqrtmu2;
864
865                 coupledConstraintsCountsHost.at(splitMap.at(c1))++;
866             }
867         }
868
869         // Constraints, coupled through the second atom.
870         c2a1 = c1a2;
871         for (const auto& atomAdjacencyList : atomsAdjacencyList[c1a2])
872         {
873             int c2 = atomAdjacencyList.indexOfConstraint_;
874
875             if (c1 != c2)
876             {
877                 int c2a2  = atomAdjacencyList.indexOfSecondConstrainedAtom_;
878                 int sign  = atomAdjacencyList.signFactor_;
879                 int index = kernelParams_.numConstraintsThreads
880                                     * coupledConstraintsCountsHost.at(splitMap.at(c1))
881                             + splitMap.at(c1);
882
883                 coupledConstraintsIndicesHost.at(index) = splitMap.at(c2);
884
885                 int center = c1a2;
886
887                 float sqrtmu1 = 1.0 / sqrt(md.invmass[c1a1] + md.invmass[c1a2]);
888                 float sqrtmu2 = 1.0 / sqrt(md.invmass[c2a1] + md.invmass[c2a2]);
889
890                 massFactorsHost.at(index) = sign * md.invmass[center] * sqrtmu1 * sqrtmu2;
891
892                 coupledConstraintsCountsHost.at(splitMap.at(c1))++;
893             }
894         }
895     }
896
897     // (Re)allocate the memory, if the number of constraints has increased.
898     if (kernelParams_.numConstraintsThreads > numConstraintsThreadsAlloc_)
899     {
900         // Free memory if it was allocated before (i.e. if not the first time here).
901         if (numConstraintsThreadsAlloc_ > 0)
902         {
903             freeDeviceBuffer(&kernelParams_.d_constraints);
904             freeDeviceBuffer(&kernelParams_.d_constraintsTargetLengths);
905
906             freeDeviceBuffer(&kernelParams_.d_coupledConstraintsCounts);
907             freeDeviceBuffer(&kernelParams_.d_coupledConstraintsIndices);
908             freeDeviceBuffer(&kernelParams_.d_massFactors);
909             freeDeviceBuffer(&kernelParams_.d_matrixA);
910         }
911
912         numConstraintsThreadsAlloc_ = kernelParams_.numConstraintsThreads;
913
914         allocateDeviceBuffer(&kernelParams_.d_constraints, kernelParams_.numConstraintsThreads, nullptr);
915         allocateDeviceBuffer(&kernelParams_.d_constraintsTargetLengths,
916                              kernelParams_.numConstraintsThreads, nullptr);
917
918         allocateDeviceBuffer(&kernelParams_.d_coupledConstraintsCounts,
919                              kernelParams_.numConstraintsThreads, nullptr);
920         allocateDeviceBuffer(&kernelParams_.d_coupledConstraintsIndices,
921                              maxCoupledConstraints * kernelParams_.numConstraintsThreads, nullptr);
922         allocateDeviceBuffer(&kernelParams_.d_massFactors,
923                              maxCoupledConstraints * kernelParams_.numConstraintsThreads, nullptr);
924         allocateDeviceBuffer(&kernelParams_.d_matrixA,
925                              maxCoupledConstraints * kernelParams_.numConstraintsThreads, nullptr);
926     }
927
928     // (Re)allocate the memory, if the number of atoms has increased.
929     if (numAtoms > numAtomsAlloc_)
930     {
931         if (numAtomsAlloc_ > 0)
932         {
933             freeDeviceBuffer(&kernelParams_.d_inverseMasses);
934         }
935         numAtomsAlloc_ = numAtoms;
936         allocateDeviceBuffer(&kernelParams_.d_inverseMasses, numAtoms, nullptr);
937     }
938
939     // Copy data to GPU.
940     copyToDeviceBuffer(&kernelParams_.d_constraints, constraintsHost.data(), 0,
941                        kernelParams_.numConstraintsThreads, commandStream_,
942                        GpuApiCallBehavior::Sync, nullptr);
943     copyToDeviceBuffer(&kernelParams_.d_constraintsTargetLengths,
944                        constraintsTargetLengthsHost.data(), 0, kernelParams_.numConstraintsThreads,
945                        commandStream_, GpuApiCallBehavior::Sync, nullptr);
946     copyToDeviceBuffer(&kernelParams_.d_coupledConstraintsCounts,
947                        coupledConstraintsCountsHost.data(), 0, kernelParams_.numConstraintsThreads,
948                        commandStream_, GpuApiCallBehavior::Sync, nullptr);
949     copyToDeviceBuffer(&kernelParams_.d_coupledConstraintsIndices, coupledConstraintsIndicesHost.data(),
950                        0, maxCoupledConstraints * kernelParams_.numConstraintsThreads,
951                        commandStream_, GpuApiCallBehavior::Sync, nullptr);
952     copyToDeviceBuffer(&kernelParams_.d_massFactors, massFactorsHost.data(), 0,
953                        maxCoupledConstraints * kernelParams_.numConstraintsThreads, commandStream_,
954                        GpuApiCallBehavior::Sync, nullptr);
955
956     GMX_RELEASE_ASSERT(md.invmass != nullptr, "Masses of atoms should be specified.\n");
957     copyToDeviceBuffer(&kernelParams_.d_inverseMasses, md.invmass, 0, numAtoms, commandStream_,
958                        GpuApiCallBehavior::Sync, nullptr);
959 }
960
961 } // namespace gmx