2 * This file is part of the GROMACS molecular simulation package.
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.
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 using CUDA.
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/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"
73 /*!\brief Number of CUDA threads in a block
75 * \todo Check if using smaller block size will lead to better prformance.
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;
81 /*! \brief Scaling matrix struct.
83 * \todo Should be generalized.
87 float xx, yy, zz, yx, zx, zy;
90 __launch_bounds__(c_maxThreadsPerBlock) __global__
91 static void scaleCoordinates_kernel(const int numAtoms,
92 float3* __restrict__ gm_x,
93 const ScalingMatrix scalingMatrix)
95 int threadIndex = blockIdx.x * blockDim.x + threadIdx.x;
96 if (threadIndex < numAtoms)
98 float3 x = gm_x[threadIndex];
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;
104 gm_x[threadIndex] = x;
108 void UpdateConstrainGpu::Impl::integrate(GpuEventSynchronizer* fReadyOnDevice,
110 const bool updateVelocities,
111 const bool computeVirial,
113 const bool doTemperatureScaling,
114 gmx::ArrayRef<const t_grp_tcstat> tcstat,
115 const bool doParrinelloRahman,
116 const float dtPressureCouple,
117 const matrix prVelocityScalingMatrix)
119 // Clearing virial matrix
120 // TODO There is no point in having separate virial matrix for constraints
123 // Make sure that the forces are ready on device before proceeding with the update.
124 fReadyOnDevice->enqueueWaitEvent(deviceStream_);
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_);
136 // scaledVirial -> virial (methods above returns scaled values)
137 float scaleFactor = 0.5f / (dt * dt);
138 for (int i = 0; i < DIM; i++)
140 for (int j = 0; j < DIM; j++)
142 virial[i][j] = scaleFactor * virial[i][j];
146 coordinatesReady_->markEvent(deviceStream_);
151 void UpdateConstrainGpu::Impl::scaleCoordinates(const matrix scalingMatrix)
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];
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();
170 UpdateConstrainGpu::Impl::Impl(const t_inputrec& ir,
171 const gmx_mtop_t& mtop,
172 const DeviceContext& deviceContext,
173 const DeviceStream& deviceStream,
174 GpuEventSynchronizer* xUpdatedOnDevice) :
175 deviceContext_(deviceContext),
176 deviceStream_(deviceStream),
177 coordinatesReady_(xUpdatedOnDevice)
179 GMX_ASSERT(xUpdatedOnDevice != nullptr, "The event synchronizer can not be nullptr.");
182 integrator_ = std::make_unique<LeapFrogGpu>(deviceContext_, deviceStream_);
183 lincsGpu_ = std::make_unique<LincsGpu>(ir.nLincsIter, ir.nProjOrder, deviceContext_, deviceStream_);
184 settleGpu_ = std::make_unique<SettleGpu>(mtop, deviceContext_, deviceStream_);
186 coordinateScalingKernelLaunchConfig_.blockSize[0] = c_threadsPerBlock;
187 coordinateScalingKernelLaunchConfig_.blockSize[1] = 1;
188 coordinateScalingKernelLaunchConfig_.blockSize[2] = 1;
189 coordinateScalingKernelLaunchConfig_.sharedMemorySize = 0;
192 UpdateConstrainGpu::Impl::~Impl() {}
194 void UpdateConstrainGpu::Impl::set(DeviceBuffer<RVec> d_x,
195 DeviceBuffer<RVec> d_v,
196 const DeviceBuffer<RVec> d_f,
197 const InteractionDefinitions& idef,
199 const int numTempScaleValues)
201 GMX_ASSERT(d_x != nullptr, "Coordinates device buffer should not be null.");
202 GMX_ASSERT(d_v != nullptr, "Velocities device buffer should not be null.");
203 GMX_ASSERT(d_f != nullptr, "Forces device buffer should not be null.");
205 d_x_ = reinterpret_cast<float3*>(d_x);
206 d_v_ = reinterpret_cast<float3*>(d_v);
207 d_f_ = reinterpret_cast<float3*>(d_f);
211 reallocateDeviceBuffer(&d_xp_, numAtoms_, &numXp_, &numXpAlloc_, deviceContext_);
213 reallocateDeviceBuffer(&d_inverseMasses_, numAtoms_, &numInverseMasses_,
214 &numInverseMassesAlloc_, deviceContext_);
216 // Integrator should also update something, but it does not even have a method yet
217 integrator_->set(numAtoms_, md.invmass, numTempScaleValues, md.cTC);
218 lincsGpu_->set(idef, numAtoms_, md.invmass);
219 settleGpu_->set(idef);
221 coordinateScalingKernelLaunchConfig_.gridSize[0] =
222 (numAtoms_ + c_threadsPerBlock - 1) / c_threadsPerBlock;
225 void UpdateConstrainGpu::Impl::setPbc(const PbcType pbcType, const matrix box)
227 setPbcAiuc(numPbcDimensions(pbcType), box, &pbcAiuc_);
230 GpuEventSynchronizer* UpdateConstrainGpu::Impl::getCoordinatesReadySync()
232 return coordinatesReady_;
235 UpdateConstrainGpu::UpdateConstrainGpu(const t_inputrec& ir,
236 const gmx_mtop_t& mtop,
237 const DeviceContext& deviceContext,
238 const DeviceStream& deviceStream,
239 GpuEventSynchronizer* xUpdatedOnDevice) :
240 impl_(new Impl(ir, mtop, deviceContext, deviceStream, xUpdatedOnDevice))
244 UpdateConstrainGpu::~UpdateConstrainGpu() = default;
246 void UpdateConstrainGpu::integrate(GpuEventSynchronizer* fReadyOnDevice,
248 const bool updateVelocities,
249 const bool computeVirial,
251 const bool doTemperatureScaling,
252 gmx::ArrayRef<const t_grp_tcstat> tcstat,
253 const bool doParrinelloRahman,
254 const float dtPressureCouple,
255 const matrix prVelocityScalingMatrix)
257 impl_->integrate(fReadyOnDevice, dt, updateVelocities, computeVirial, virialScaled, doTemperatureScaling,
258 tcstat, doParrinelloRahman, dtPressureCouple, prVelocityScalingMatrix);
261 void UpdateConstrainGpu::scaleCoordinates(const matrix scalingMatrix)
263 impl_->scaleCoordinates(scalingMatrix);
266 void UpdateConstrainGpu::set(DeviceBuffer<RVec> d_x,
267 DeviceBuffer<RVec> d_v,
268 const DeviceBuffer<RVec> d_f,
269 const InteractionDefinitions& idef,
271 const int numTempScaleValues)
273 impl_->set(d_x, d_v, d_f, idef, md, numTempScaleValues);
276 void UpdateConstrainGpu::setPbc(const PbcType pbcType, const matrix box)
278 impl_->setPbc(pbcType, box);
281 GpuEventSynchronizer* UpdateConstrainGpu::getCoordinatesReadySync()
283 return impl_->getCoordinatesReadySync();
286 bool UpdateConstrainGpu::isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop)
288 return LincsGpu::isNumCoupledConstraintsSupported(mtop);