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