96d447d36c22d4851b6fa87a87600363db245af7
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_gpu_impl.cpp
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.
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/device_context.h"
60 #include "gromacs/gpu_utils/device_stream.h"
61 #include "gromacs/gpu_utils/devicebuffer.h"
62 #if GMX_GPU_CUDA
63 #    include "gromacs/gpu_utils/gpueventsynchronizer.cuh"
64 #elif GMX_GPU_SYCL
65 #    include "gromacs/gpu_utils/gpueventsynchronizer_sycl.h"
66 #endif
67 #include "gromacs/gpu_utils/gputraits.h"
68 #include "gromacs/mdlib/leapfrog_gpu.h"
69 #include "gromacs/mdlib/update_constrain_gpu.h"
70 #include "gromacs/mdlib/update_constrain_gpu_internal.h"
71 #include "gromacs/mdtypes/mdatom.h"
72 #include "gromacs/timing/wallcycle.h"
73 #include "gromacs/topology/mtop_util.h"
74
75 static constexpr bool sc_haveGpuConstraintSupport = GMX_GPU_CUDA;
76
77 namespace gmx
78 {
79
80 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer*             fReadyOnDevice,
81                                          const real                        dt,
82                                          const bool                        updateVelocities,
83                                          const bool                        computeVirial,
84                                          tensor                            virial,
85                                          const bool                        doTemperatureScaling,
86                                          gmx::ArrayRef<const t_grp_tcstat> tcstat,
87                                          const bool                        doParrinelloRahman,
88                                          const float                       dtPressureCouple,
89                                          const matrix                      prVelocityScalingMatrix)
90 {
91     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
92     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
93
94     // Clearing virial matrix
95     // TODO There is no point in having separate virial matrix for constraints
96     clear_mat(virial);
97
98     // Make sure that the forces are ready on device before proceeding with the update.
99     fReadyOnDevice->enqueueWaitEvent(deviceStream_);
100
101     // The integrate should save a copy of the current coordinates in d_xp_ and write updated
102     // once into d_x_. The d_xp_ is only needed by constraints.
103     integrator_->integrate(
104             d_x_, d_xp_, d_v_, d_f_, dt, doTemperatureScaling, tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
105     // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
106     // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
107     // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
108     if (sc_haveGpuConstraintSupport)
109     {
110         lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
111         settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
112     }
113
114     // scaledVirial -> virial (methods above returns scaled values)
115     float scaleFactor = 0.5F / (dt * dt);
116     for (int i = 0; i < DIM; i++)
117     {
118         for (int j = 0; j < DIM; j++)
119         {
120             virial[i][j] = scaleFactor * virial[i][j];
121         }
122     }
123
124     coordinatesReady_->markEvent(deviceStream_);
125
126     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
127     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
128 }
129
130 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
131 {
132     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
133     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
134
135     ScalingMatrix mu(scalingMatrix);
136
137     launchScaleCoordinatesKernel(numAtoms_, d_x_, mu, deviceStream_);
138
139     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
140     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
141 }
142
143 void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
144 {
145     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
146     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
147
148     ScalingMatrix mu(scalingMatrix);
149
150     launchScaleCoordinatesKernel(numAtoms_, d_v_, mu, deviceStream_);
151
152     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
153     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
154 }
155
156 UpdateConstrainGpu::Impl::Impl(const t_inputrec&     ir,
157                                const gmx_mtop_t&     mtop,
158                                const int             numTempScaleValues,
159                                const DeviceContext&  deviceContext,
160                                const DeviceStream&   deviceStream,
161                                GpuEventSynchronizer* xUpdatedOnDevice,
162                                gmx_wallcycle*        wcycle) :
163     deviceContext_(deviceContext),
164     deviceStream_(deviceStream),
165     coordinatesReady_(xUpdatedOnDevice),
166     wcycle_(wcycle)
167 {
168     GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
169
170     integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_, numTempScaleValues);
171     if (sc_haveGpuConstraintSupport)
172     {
173         lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
174         settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
175     }
176 }
177
178 UpdateConstrainGpu::Impl::~Impl() {}
179
180 void UpdateConstrainGpu::Impl::set(DeviceBuffer<Float3>          d_x,
181                                    DeviceBuffer<Float3>          d_v,
182                                    const DeviceBuffer<Float3>    d_f,
183                                    const InteractionDefinitions& idef,
184                                    const t_mdatoms&              md)
185 {
186     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
187     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
188
189     GMX_ASSERT(d_x, "Coordinates device buffer should not be null.");
190     GMX_ASSERT(d_v, "Velocities device buffer should not be null.");
191     GMX_ASSERT(d_f, "Forces device buffer should not be null.");
192
193     d_x_ = d_x;
194     d_v_ = d_v;
195     d_f_ = d_f;
196
197     numAtoms_ = md.nr;
198
199     reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
200
201     reallocateDeviceBuffer(
202             &d_inverseMasses_, numAtoms_, &numInverseMasses_, &numInverseMassesAlloc_, deviceContext_);
203
204     // Integrator should also update something, but it does not even have a method yet
205     integrator_->set(numAtoms_, md.invmass, md.cTC);
206     if (sc_haveGpuConstraintSupport)
207     {
208         lincsGpu_->set(idef, numAtoms_, md.invmass);
209         settleGpu_->set(idef);
210     }
211     else
212     {
213         GMX_ASSERT(idef.il[F_SETTLE].empty(), "SETTLE not supported");
214         GMX_ASSERT(idef.il[F_CONSTR].empty(), "LINCS not supported");
215     }
216
217     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
218     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
219 }
220
221 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
222 {
223     // TODO wallcycle
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 int             numTempScaleValues,
235                                        const DeviceContext&  deviceContext,
236                                        const DeviceStream&   deviceStream,
237                                        GpuEventSynchronizer* xUpdatedOnDevice,
238                                        gmx_wallcycle*        wcycle) :
239     impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
240 {
241 }
242
243 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
244
245 void UpdateConstrainGpu::integrate(GpuEventSynchronizer*             fReadyOnDevice,
246                                    const real                        dt,
247                                    const bool                        updateVelocities,
248                                    const bool                        computeVirial,
249                                    tensor                            virialScaled,
250                                    const bool                        doTemperatureScaling,
251                                    gmx::ArrayRef<const t_grp_tcstat> tcstat,
252                                    const bool                        doParrinelloRahman,
253                                    const float                       dtPressureCouple,
254                                    const matrix                      prVelocityScalingMatrix)
255 {
256     impl_->integrate(fReadyOnDevice,
257                      dt,
258                      updateVelocities,
259                      computeVirial,
260                      virialScaled,
261                      doTemperatureScaling,
262                      tcstat,
263                      doParrinelloRahman,
264                      dtPressureCouple,
265                      prVelocityScalingMatrix);
266 }
267
268 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
269 {
270     impl_->scaleCoordinates(scalingMatrix);
271 }
272
273 void UpdateConstrainGpu::scaleVelocities(const matrix scalingMatrix)
274 {
275     impl_->scaleVelocities(scalingMatrix);
276 }
277
278 void UpdateConstrainGpu::set(DeviceBuffer<Float3>          d_x,
279                              DeviceBuffer<Float3>          d_v,
280                              const DeviceBuffer<Float3>    d_f,
281                              const InteractionDefinitions& idef,
282                              const t_mdatoms&              md)
283 {
284     impl_->set(d_x, d_v, d_f, idef, md);
285 }
286
287 void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
288 {
289     impl_->setPbc(pbcType, box);
290 }
291
292 GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
293 {
294     return impl_->getCoordinatesReadySync();
295 }
296
297 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
298 {
299     return LincsGpu::isNumCoupledConstraintsSupported(mtop);
300 }
301
302 bool UpdateConstrainGpu::areConstraintsSupported()
303 {
304     return sc_haveGpuConstraintSupport;
305 }
306
307 } // namespace gmx