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