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