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