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