2 * This file is part of the GROMACS molecular simulation package.
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.
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.
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.
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.
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.
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.
37 * \brief Implements update and constraints class.
39 * The class combines Leap-Frog integrator with LINCS and SETTLE constraints.
41 * \todo The computational procedures in members should be integrated to improve
42 * computational performance.
44 * \author Artem Zhmurov <zhmurov@gmail.com>
46 * \ingroup module_mdlib
50 #include "update_constrain_gpu_impl.h"
59 #include "gromacs/gpu_utils/device_context.h"
60 #include "gromacs/gpu_utils/device_stream.h"
61 #include "gromacs/gpu_utils/devicebuffer.h"
62 #include "gromacs/gpu_utils/gpueventsynchronizer.h"
63 #include "gromacs/gpu_utils/gputraits.h"
64 #include "gromacs/mdlib/leapfrog_gpu.h"
65 #include "gromacs/mdlib/update_constrain_gpu.h"
66 #include "gromacs/mdlib/update_constrain_gpu_internal.h"
67 #include "gromacs/mdtypes/mdatom.h"
68 #include "gromacs/timing/wallcycle.h"
69 #include "gromacs/topology/mtop_util.h"
71 static constexpr bool sc_haveGpuConstraintSupport = GMX_GPU_CUDA;
76 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer* fReadyOnDevice,
78 const bool updateVelocities,
79 const bool computeVirial,
81 const bool doTemperatureScaling,
82 gmx::ArrayRef<const t_grp_tcstat> tcstat,
83 const bool doParrinelloRahman,
84 const float dtPressureCouple,
85 const matrix prVelocityScalingMatrix)
87 wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
88 wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
90 // Clearing virial matrix
91 // TODO There is no point in having separate virial matrix for constraints
94 // Make sure that the forces are ready on device before proceeding with the update.
95 fReadyOnDevice->enqueueWaitEvent(deviceStream_);
97 // The integrate should save a copy of the current coordinates in d_xp_ and write updated
98 // once into d_x_. The d_xp_ is only needed by constraints.
99 integrator_->integrate(
100 d_x_, d_xp_, d_v_, d_f_, dt, doTemperatureScaling, tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
101 // Constraints need both coordinates before (d_x_) and after (d_xp_) update. However, after constraints
102 // are applied, the d_x_ can be discarded. So we intentionally swap the d_x_ and d_xp_ here to avoid the
103 // d_xp_ -> d_x_ copy after constraints. Note that the integrate saves them in the wrong order as well.
104 if (sc_haveGpuConstraintSupport)
106 lincsGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
107 settleGpu_->apply(d_xp_, d_x_, updateVelocities, d_v_, 1.0 / dt, computeVirial, virial, pbcAiuc_);
110 // scaledVirial -> virial (methods above returns scaled values)
111 float scaleFactor = 0.5F / (dt * dt);
112 for (int i = 0; i < DIM; i++)
114 for (int j = 0; j < DIM; j++)
116 virial[i][j] = scaleFactor * virial[i][j];
120 coordinatesReady_->markEvent(deviceStream_);
122 wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
123 wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
126 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
128 wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
129 wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
131 ScalingMatrix mu(scalingMatrix);
133 launchScaleCoordinatesKernel(numAtoms_, d_x_, mu, deviceStream_);
135 wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
136 wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
139 void UpdateConstrainGpu::Impl::scaleVelocities(const matrix scalingMatrix)
141 wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
142 wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
144 ScalingMatrix mu(scalingMatrix);
146 launchScaleCoordinatesKernel(numAtoms_, d_v_, mu, deviceStream_);
148 wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
149 wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
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),
164 GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
166 integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_, numTempScaleValues);
167 if (sc_haveGpuConstraintSupport)
169 lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
170 settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
174 UpdateConstrainGpu::Impl::~Impl() {}
176 void UpdateConstrainGpu::Impl::set(DeviceBuffer<Float3> d_x,
177 DeviceBuffer<Float3> d_v,
178 const DeviceBuffer<Float3> d_f,
179 const InteractionDefinitions& idef,
182 wallcycle_start_nocount(wcycle_, WallCycleCounter::LaunchGpu);
183 wallcycle_sub_start(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
185 GMX_ASSERT(d_x, "Coordinates device buffer should not be null.");
186 GMX_ASSERT(d_v, "Velocities device buffer should not be null.");
187 GMX_ASSERT(d_f, "Forces device buffer should not be null.");
195 reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
197 reallocateDeviceBuffer(
198 &d_inverseMasses_, numAtoms_, &numInverseMasses_, &numInverseMassesAlloc_, deviceContext_);
200 // Integrator should also update something, but it does not even have a method yet
201 integrator_->set(numAtoms_, md.invmass, md.cTC);
202 if (sc_haveGpuConstraintSupport)
204 lincsGpu_->set(idef, numAtoms_, md.invmass);
205 settleGpu_->set(idef);
209 GMX_ASSERT(idef.il[F_SETTLE].empty(), "SETTLE not supported");
210 GMX_ASSERT(idef.il[F_CONSTR].empty(), "LINCS not supported");
213 wallcycle_sub_stop(wcycle_, WallCycleSubCounter::LaunchGpuUpdateConstrain);
214 wallcycle_stop(wcycle_, WallCycleCounter::LaunchGpu);
217 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
220 setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
223 GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
225 return coordinatesReady_;
228 UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec& ir,
229 const gmx_mtop_t& mtop,
230 const int numTempScaleValues,
231 const DeviceContext& deviceContext,
232 const DeviceStream& deviceStream,
233 GpuEventSynchronizer* xUpdatedOnDevice,
234 gmx_wallcycle* wcycle) :
235 impl_(new Impl(ir, mtop, numTempScaleValues, deviceContext, deviceStream, xUpdatedOnDevice, wcycle))
239 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
241 void UpdateConstrainGpu::integrate(GpuEventSynchronizer* fReadyOnDevice,
243 const bool updateVelocities,
244 const bool computeVirial,
246 const bool doTemperatureScaling,
247 gmx::ArrayRef<const t_grp_tcstat> tcstat,
248 const bool doParrinelloRahman,
249 const float dtPressureCouple,
250 const matrix prVelocityScalingMatrix)
252 impl_->integrate(fReadyOnDevice,
257 doTemperatureScaling,
261 prVelocityScalingMatrix);
264 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
266 impl_->scaleCoordinates(scalingMatrix);
269 void UpdateConstrainGpu::scaleVelocities(const matrix scalingMatrix)
271 impl_->scaleVelocities(scalingMatrix);
274 void UpdateConstrainGpu::set(DeviceBuffer<Float3> d_x,
275 DeviceBuffer<Float3> d_v,
276 const DeviceBuffer<Float3> d_f,
277 const InteractionDefinitions& idef,
280 impl_->set(d_x, d_v, d_f, idef, md);
283 void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
285 impl_->setPbc(pbcType, box);
288 GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
290 return impl_->getCoordinatesReadySync();
293 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
295 return LincsGpu::isNumCoupledConstraintsSupported(mtop);
298 bool UpdateConstrainGpu::areConstraintsSupported()
300 return sc_haveGpuConstraintSupport;