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