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