Wrap more device pointers in DeviceBuffer
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_gpu_impl.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 update and constraints class using CUDA.
38  *
39  * The class combines Leap-Frog integrator with LINCS and SETTLE constraints.
40  *
41  * \todo The computational procedures in members should be integrated to improve
42  *       computational performance.
43  *
44  * \author Artem Zhmurov <zhmurov@gmail.com>
45  *
46  * \ingroup module_mdlib
47  */
48 #include "gmxpre.h"
49
50 #include "update_constrain_gpu_impl.h"
51
52 #include <assert.h>
53 #include <stdio.h>
54
55 #include <cmath>
56
57 #include <algorithm>
58
59 #include "gromacs/gpu_utils/cudautils.cuh"
60 #include "gromacs/gpu_utils/device_context.h"
61 #include "gromacs/gpu_utils/device_stream.h"
62 #include "gromacs/gpu_utils/devicebuffer.h"
63 #include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
64 #include "gromacs/gpu_utils/gputraits.cuh"
65 #include "gromacs/gpu_utils/vectype_ops.cuh"
66 #include "gromacs/mdlib/leapfrog_gpu.h"
67 #include "gromacs/mdlib/lincs_gpu.cuh"
68 #include "gromacs/mdlib/settle_gpu.cuh"
69 #include "gromacs/mdlib/update_constrain_gpu.h"
70 #include "gromacs/mdtypes/mdatom.h"
71 #include "gromacs/timing/wallcycle.h"
72
73 namespace gmx
74 {
75 /*!\brief Number of CUDA threads in a block
76  *
77  * \todo Check if using smaller block size will lead to better performance.
78  */
79 constexpr static int c_threadsPerBlock = 256;
80 //! Maximum number of threads in a block (for __launch_bounds__)
81 constexpr static int c_maxThreadsPerBlock = c_threadsPerBlock;
82
83 /*! \brief Scaling matrix struct.
84  *
85  * \todo Should be generalized.
86  */
87 struct ScalingMatrix
88 {
89     float xx, yy, zz, yx, zx, zy;
90 };
91
92 __launch_bounds__(c_maxThreadsPerBlock) __global__
93         static void scaleCoordinates_kernel(const int numAtoms,
94                                             float3* __restrict__ gm_x,
95                                             const ScalingMatrix scalingMatrix)
96 {
97     int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
98     if (threadIndex < numAtoms)
99     {
100         float3 x = gm_x[threadIndex];
101
102         x.x = scalingMatrix.xx * x.x + scalingMatrix.yx * x.y + scalingMatrix.zx * x.z;
103         x.y = scalingMatrix.yy * x.y + scalingMatrix.zy * x.z;
104         x.z = scalingMatrix.zz * x.z;
105
106         gm_x[threadIndex] = x;
107     }
108 }
109
110 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer*             fReadyOnDevice,
111                                          const real                        dt,
112                                          const bool                        updateVelocities,
113                                          const bool                        computeVirial,
114                                          tensor                            virial,
115                                          const bool                        doTemperatureScaling,
116                                          gmx::ArrayRef<const t_grp_tcstat> tcstat,
117                                          const bool                        doParrinelloRahman,
118                                          const float                       dtPressureCouple,
119                                          const matrix                      prVelocityScalingMatrix)
120 {
121     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
122     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
123
124     // Clearing virial matrix
125     // TODO There is no point in having separate virial matrix for constraints
126     clear_mat(virial);
127
128     // Make sure that the forces are ready on device before proceeding with the update.
129     fReadyOnDevice->enqueueWaitEvent(deviceStream_);
130
131     // The integrate should save a copy of the current coordinates in d_xp_ and write updated
132     // once into d_x_. The d_xp_ is only needed by constraints.
133     integrator_->integrate(
134             d_x_, d_xp_, d_v_, d_f_, dt, doTemperatureScaling, tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
135     // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
136     // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
137     // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
138     lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
139     settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
140
141     // scaledVirial -> virial (methods above returns scaled values)
142     float scaleFactor = 0.5f / (dt * dt);
143     for (int i = 0; i < DIM; i++)
144     {
145         for (int j = 0; j < DIM; j++)
146         {
147             virial[i][j] = scaleFactor * virial[i][j];
148         }
149     }
150
151     coordinatesReady_->markEvent(deviceStream_);
152
153     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
154     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
155
156     return;
157 }
158
159 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
160 {
161     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
162     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
163
164     ScalingMatrix mu;
165     mu.xx = scalingMatrix[XX][XX];
166     mu.yy = scalingMatrix[YY][YY];
167     mu.zz = scalingMatrix[ZZ][ZZ];
168     mu.yx = scalingMatrix[YY][XX];
169     mu.zx = scalingMatrix[ZZ][XX];
170     mu.zy = scalingMatrix[ZZ][YY];
171
172     const auto kernelArgs = prepareGpuKernelArguments(
173             scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_x_, &mu);
174     launchGpuKernel(scaleCoordinates_kernel,
175                     coordinateScalingKernelLaunchConfig_,
176                     deviceStream_,
177                     nullptr,
178                     "scaleCoordinates_kernel",
179                     kernelArgs);
180     // TODO: Although this only happens on the pressure coupling steps, this synchronization
181     //       can affect the performance if nstpcouple is small.
182     deviceStream_.synchronize();
183
184     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
185     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
186 }
187
188 void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
189 {
190     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
191     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
192
193     ScalingMatrix mu;
194     mu.xx = scalingMatrix[XX][XX];
195     mu.yy = scalingMatrix[YY][YY];
196     mu.zz = scalingMatrix[ZZ][ZZ];
197     mu.yx = scalingMatrix[YY][XX];
198     mu.zx = scalingMatrix[ZZ][XX];
199     mu.zy = scalingMatrix[ZZ][YY];
200
201     const auto kernelArgs = prepareGpuKernelArguments(
202             scaleCoordinates_kernel, coordinateScalingKernelLaunchConfig_, &numAtoms_, &d_v_, &mu);
203     launchGpuKernel(scaleCoordinates_kernel,
204                     coordinateScalingKernelLaunchConfig_,
205                     deviceStream_,
206                     nullptr,
207                     "scaleCoordinates_kernel",
208                     kernelArgs);
209     // TODO: Although this only happens on the pressure coupling steps, this synchronization
210     //       can affect the performance if nstpcouple is small.
211     deviceStream_.synchronize();
212
213     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
214     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
215 }
216
217 UpdateConstrainGpu::Impl::Impl(const t_inputrec&     ir,
218                                const gmx_mtop_t&     mtop,
219                                const int             numTempScaleValues,
220                                const DeviceContext&  deviceContext,
221                                const DeviceStream&   deviceStream,
222                                GpuEventSynchronizer* xUpdatedOnDevice,
223                                gmx_wallcycle*        wcycle) :
224     deviceContext_(deviceContext),
225     deviceStream_(deviceStream),
226     coordinatesReady_(xUpdatedOnDevice),
227     wcycle_(wcycle)
228 {
229     GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
230
231
232     integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_, numTempScaleValues);
233     lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
234     settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
235
236     coordinateScalingKernelLaunchConfig_.blockSize[0]     = c_threadsPerBlock;
237     coordinateScalingKernelLaunchConfig_.blockSize[1]     = 1;
238     coordinateScalingKernelLaunchConfig_.blockSize[2]     = 1;
239     coordinateScalingKernelLaunchConfig_.sharedMemorySize = 0;
240 }
241
242 UpdateConstrainGpu::Impl::~Impl() {}
243
244 void UpdateConstrainGpu::Impl::set(DeviceBuffer<RVec>            d_x,
245                                    DeviceBuffer<RVec>            d_v,
246                                    const DeviceBuffer<RVec>      d_f,
247                                    const InteractionDefinitions& idef,
248                                    const t_mdatoms&              md)
249 {
250     // TODO wallcycle
251     wallcycle_start_nocount(wcycle_, ewcLAUNCH_GPU);
252     wallcycle_sub_start(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
253
254     GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
255     GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
256     GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
257
258     d_x_ = reinterpret_cast<float3*>(d_x);
259     d_v_ = reinterpret_cast<float3*>(d_v);
260     d_f_ = reinterpret_cast<float3*>(d_f);
261
262     numAtoms_ = md.nr;
263
264     reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
265
266     reallocateDeviceBuffer(
267             &d_inverseMasses_, numAtoms_, &numInverseMasses_, &numInverseMassesAlloc_, deviceContext_);
268
269     // Integrator should also update something, but it does not even have a method yet
270     integrator_->set(numAtoms_, md.invmass, md.cTC);
271     lincsGpu_->set(idef, numAtoms_, md.invmass);
272     settleGpu_->set(idef);
273
274     coordinateScalingKernelLaunchConfig_.gridSize[0] =
275             (numAtoms_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
276
277     wallcycle_sub_stop(wcycle_, ewcsLAUNCH_GPU_UPDATE_CONSTRAIN);
278     wallcycle_stop(wcycle_, ewcLAUNCH_GPU);
279 }
280
281 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
282 {
283     // TODO wallcycle
284     setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
285 }
286
287 GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
288 {
289     return coordinatesReady_;
290 }
291
292 UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec&     ir,
293                                        const gmx_mtop_t&     mtop,
294                                        const int             numTempScaleValues,
295                                        const DeviceContext&  deviceContext,
296                                        const DeviceStream&   deviceStream,
297                                        GpuEventSynchronizer* xUpdatedOnDevice,
298                                        gmx_wallcycle*        wcycle) :
299     impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
300 {
301 }
302
303 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
304
305 void UpdateConstrainGpu::integrate(GpuEventSynchronizer*             fReadyOnDevice,
306                                    const real                        dt,
307                                    const bool                        updateVelocities,
308                                    const bool                        computeVirial,
309                                    tensor                            virialScaled,
310                                    const bool                        doTemperatureScaling,
311                                    gmx::ArrayRef<const t_grp_tcstat> tcstat,
312                                    const bool                        doParrinelloRahman,
313                                    const float                       dtPressureCouple,
314                                    const matrix                      prVelocityScalingMatrix)
315 {
316     impl_->integrate(fReadyOnDevice,
317                      dt,
318                      updateVelocities,
319                      computeVirial,
320                      virialScaled,
321                      doTemperatureScaling,
322                      tcstat,
323                      doParrinelloRahman,
324                      dtPressureCouple,
325                      prVelocityScalingMatrix);
326 }
327
328 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
329 {
330     impl_->scaleCoordinates(scalingMatrix);
331 }
332
333 void UpdateConstrainGpu::scaleVelocities(const matrix scalingMatrix)
334 {
335     impl_->scaleVelocities(scalingMatrix);
336 }
337
338 void UpdateConstrainGpu::set(DeviceBuffer<RVec>            d_x,
339                              DeviceBuffer<RVec>            d_v,
340                              const DeviceBuffer<RVec>      d_f,
341                              const InteractionDefinitions& idef,
342                              const t_mdatoms&              md)
343 {
344     impl_->set(d_x, d_v, d_f, idef, md);
345 }
346
347 void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
348 {
349     impl_->setPbc(pbcType, box);
350 }
351
352 GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
353 {
354     return impl_->getCoordinatesReadySync();
355 }
356
357 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
358 {
359     return LincsGpu::isNumCoupledConstraintsSupported(mtop);
360 }
361
362 } // namespace gmx