Make common Update-Constraints code backend-agnostic
[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
74 namespace gmx
75 {
76
77 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer*             fReadyOnDevice,
78                                          const real                        dt,
79                                          const bool                        updateVelocities,
80                                          const bool                        computeVirial,
81                                          tensor                            virial,
82                                          const bool                        doTemperatureScaling,
83                                          gmx::ArrayRef<const t_grp_tcstat> tcstat,
84                                          const bool                        doParrinelloRahman,
85                                          const float                       dtPressureCouple,
86                                          const matrix                      prVelocityScalingMatrix)
87 {
88     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
89     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
90
91     // Clearing virial matrix
92     // TODO There is no point in having separate virial matrix for constraints
93     clear_mat(virial);
94
95     // Make sure that the forces are ready on device before proceeding with the update.
96     fReadyOnDevice->enqueueWaitEvent(deviceStream_);
97
98     // The integrate should save a copy of the current coordinates in d_xp_ and write updated
99     // once into d_x_. The d_xp_ is only needed by constraints.
100     integrator_->integrate(
101             d_x_, d_xp_, d_v_, d_f_, dt, doTemperatureScaling, tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
102     // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
103     // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
104     // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
105     lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
106     settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
107
108     // scaledVirial -> virial (methods above returns scaled values)
109     float scaleFactor = 0.5F / (dt * dt);
110     for (int i = 0; i < DIM; i++)
111     {
112         for (int j = 0; j < DIM; j++)
113         {
114             virial[i][j] = scaleFactor * virial[i][j];
115         }
116     }
117
118     coordinatesReady_->markEvent(deviceStream_);
119
120     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
121     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
122
123     return;
124 }
125
126 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
127 {
128     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
129     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
130
131     ScalingMatrix mu(scalingMatrix);
132
133     launchScaleCoordinatesKernel(numAtoms_, d_x_, mu, deviceStream_);
134
135     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
136     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
137 }
138
139 void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
140 {
141     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
142     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
143
144     ScalingMatrix mu(scalingMatrix);
145
146     launchScaleCoordinatesKernel(numAtoms_, d_v_, mu, deviceStream_);
147
148     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
149     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
150 }
151
152 UpdateConstrainGpu::Impl::Impl(const t_inputrec&     ir,
153                                const gmx_mtop_t&     mtop,
154                                const int             numTempScaleValues,
155                                const DeviceContext&  deviceContext,
156                                const DeviceStream&   deviceStream,
157                                GpuEventSynchronizer* xUpdatedOnDevice,
158                                gmx_wallcycle*        wcycle) :
159     deviceContext_(deviceContext),
160     deviceStream_(deviceStream),
161     coordinatesReady_(xUpdatedOnDevice),
162     wcycle_(wcycle)
163 {
164     GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
165
166     integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_, numTempScaleValues);
167     lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
168     settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
169 }
170
171 UpdateConstrainGpu::Impl::~Impl() {}
172
173 void UpdateConstrainGpu::Impl::set(DeviceBuffer<Float3>          d_x,
174                                    DeviceBuffer<Float3>          d_v,
175                                    const DeviceBuffer<Float3>    d_f,
176                                    const InteractionDefinitions& idef,
177                                    const t_mdatoms&              md)
178 {
179     wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
180     wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
181
182     GMX_ASSERT(d_x, "Coordinates device buffer should not be null.");
183     GMX_ASSERT(d_v, "Velocities device buffer should not be null.");
184     GMX_ASSERT(d_f, "Forces device buffer should not be null.");
185
186     d_x_ = d_x;
187     d_v_ = d_v;
188     d_f_ = d_f;
189
190     numAtoms_ = md.nr;
191
192     reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
193
194     reallocateDeviceBuffer(
195             &d_inverseMasses_, numAtoms_, &numInverseMasses_, &numInverseMassesAlloc_, deviceContext_);
196
197     // Integrator should also update something, but it does not even have a method yet
198     integrator_->set(numAtoms_, md.invmass, md.cTC);
199     lincsGpu_->set(idef, numAtoms_, md.invmass);
200     settleGpu_->set(idef);
201
202     wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
203     wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
204 }
205
206 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
207 {
208     // TODO wallcycle
209     setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
210 }
211
212 GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
213 {
214     return coordinatesReady_;
215 }
216
217 UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec&     ir,
218                                        const gmx_mtop_t&     mtop,
219                                        const int             numTempScaleValues,
220                                        const DeviceContext&  deviceContext,
221                                        const DeviceStream&   deviceStream,
222                                        GpuEventSynchronizer* xUpdatedOnDevice,
223                                        gmx_wallcycle*        wcycle) :
224     impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
225 {
226 }
227
228 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
229
230 void UpdateConstrainGpu::integrate(GpuEventSynchronizer*             fReadyOnDevice,
231                                    const real                        dt,
232                                    const bool                        updateVelocities,
233                                    const bool                        computeVirial,
234                                    tensor                            virialScaled,
235                                    const bool                        doTemperatureScaling,
236                                    gmx::ArrayRef<const t_grp_tcstat> tcstat,
237                                    const bool                        doParrinelloRahman,
238                                    const float                       dtPressureCouple,
239                                    const matrix                      prVelocityScalingMatrix)
240 {
241     impl_->integrate(fReadyOnDevice,
242                      dt,
243                      updateVelocities,
244                      computeVirial,
245                      virialScaled,
246                      doTemperatureScaling,
247                      tcstat,
248                      doParrinelloRahman,
249                      dtPressureCouple,
250                      prVelocityScalingMatrix);
251 }
252
253 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
254 {
255     impl_->scaleCoordinates(scalingMatrix);
256 }
257
258 void UpdateConstrainGpu::scaleVelocities(const matrix scalingMatrix)
259 {
260     impl_->scaleVelocities(scalingMatrix);
261 }
262
263 void UpdateConstrainGpu::set(DeviceBuffer<Float3>          d_x,
264                              DeviceBuffer<Float3>          d_v,
265                              const DeviceBuffer<Float3>    d_f,
266                              const InteractionDefinitions& idef,
267                              const t_mdatoms&              md)
268 {
269     impl_->set(d_x, d_v, d_f, idef, md);
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