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