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 Declares GPU implementation class for update and constraints.
39 * This header file is needed to include from both the device-side
40 * kernels file, and the host-side management code.
42 * \author Artem Zhmurov <zhmurov@gmail.com>
44 * \ingroup module_mdlib
46 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H
47 #define GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H
49 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
50 #include "gromacs/mdlib/update_constrain_gpu.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/pbcutil/pbc_aiuc.h"
55 # include "gromacs/gpu_utils/gputraits.cuh"
58 class GpuEventSynchronizer;
69 /*! \internal \brief Class with interfaces and data for GPU version of Update-Constraint. */
70 class UpdateConstrainGpu::Impl
74 /*! \brief Create Update-Constrain object.
76 * The constructor is given a non-nullptr \p deviceStream, in which all the update and constrain
77 * routines are executed. \p xUpdatedOnDevice should mark the completion of all kernels that
78 * modify coordinates. The event is maintained outside this class and also passed to all (if
79 * any) consumers of the updated coordinates. The \p xUpdatedOnDevice also can not be a nullptr
80 * because the markEvent(...) method is called unconditionally.
82 * \param[in] ir Input record data: LINCS takes number of iterations and order of
84 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
85 * and target O-H and H-H distances from this object.
86 * \param[in] numTempScaleValues Number of temperature scaling groups. Set zero for no temperature coupling.
87 * \param[in] deviceContext GPU device context.
88 * \param[in] deviceStream GPU stream to use.
89 * \param[in] xUpdatedOnDevice The event synchronizer to use to mark that
90 * update is done on the GPU.
91 * \param[in] wcycle The wallclock counter
93 Impl(const t_inputrec& ir,
94 const gmx_mtop_t& mtop,
95 int numTempScaleValues,
96 const DeviceContext& deviceContext,
97 const DeviceStream& deviceStream,
98 GpuEventSynchronizer* xUpdatedOnDevice,
99 gmx_wallcycle* wcycle);
105 * Integrates the equation of motion using Leap-Frog algorithm and applies
106 * LINCS and SETTLE constraints.
107 * If computeVirial is true, constraints virial is written at the provided pointer.
108 * doTempCouple should be true if:
109 * 1. The temperature coupling is enabled.
110 * 2. This is the temperature coupling step.
111 * Parameters virial/lambdas can be nullptr if computeVirial/doTempCouple are false.
113 * \param[in] fReadyOnDevice Event synchronizer indicating that the forces are ready in
115 * \param[in] dt Timestep.
116 * \param[in] updateVelocities If the velocities should be constrained.
117 * \param[in] computeVirial If virial should be updated.
118 * \param[out] virial Place to save virial tensor.
119 * \param[in] doTemperatureScaling If velocities should be scaled for temperature coupling.
120 * \param[in] tcstat Temperature coupling data.
121 * \param[in] doParrinelloRahman If current step is a Parrinello-Rahman pressure coupling step.
122 * \param[in] dtPressureCouple Period between pressure coupling steps.
123 * \param[in] prVelocityScalingMatrix Parrinello-Rahman velocity scaling matrix.
125 void integrate(GpuEventSynchronizer* fReadyOnDevice,
127 bool updateVelocities,
130 bool doTemperatureScaling,
131 gmx::ArrayRef<const t_grp_tcstat> tcstat,
132 bool doParrinelloRahman,
133 float dtPressureCouple,
134 const matrix prVelocityScalingMatrix);
136 /*! \brief Scale coordinates on the GPU for the pressure coupling.
138 * After pressure coupling step, the box size may change. Hence, the coordinates should be
139 * scaled so that all the particles fit in the new box.
141 * \param[in] scalingMatrix Coordinates scaling matrix.
143 void scaleCoordinates(const matrix scalingMatrix);
145 /*! \brief Scale velocities on the GPU for the pressure coupling.
147 * After pressure coupling step, the box size may change. In the C-Rescale algorithm, velocities should be scaled.
149 * \param[in] scalingMatrix Velocities scaling matrix.
151 void scaleVelocities(const matrix scalingMatrix);
153 /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
155 * \param[in,out] d_x Device buffer with coordinates.
156 * \param[in,out] d_v Device buffer with velocities.
157 * \param[in] d_f Device buffer with forces.
158 * \param[in] idef System topology
159 * \param[in] md Atoms data.
161 void set(DeviceBuffer<RVec> d_x,
162 DeviceBuffer<RVec> d_v,
163 const DeviceBuffer<RVec> d_f,
164 const InteractionDefinitions& idef,
165 const t_mdatoms& md);
170 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
172 * \param[in] pbcType The type of the periodic boundary.
173 * \param[in] box The periodic boundary box matrix.
175 void setPbc(PbcType pbcType, const matrix box);
177 /*! \brief Return the synchronizer associated with the event indicated that the coordinates are ready on the device.
179 GpuEventSynchronizer* getCoordinatesReadySync();
182 * Returns whether the maximum number of coupled constraints is supported
183 * by the GPU LINCS code.
185 * \param[in] mtop The molecular topology
187 static bool isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop);
190 //! GPU context object
191 const DeviceContext& deviceContext_;
193 const DeviceStream& deviceStream_;
194 //! GPU kernel launch config
195 KernelLaunchConfig coordinateScalingKernelLaunchConfig_;
197 //! Periodic boundary data
203 //! Local copy of the pointer to the device positions buffer
204 DeviceBuffer<float3> d_x_;
205 //! Local copy of the pointer to the device velocities buffer
206 DeviceBuffer<float3> d_v_;
207 //! Local copy of the pointer to the device forces buffer
208 DeviceBuffer<float3> d_f_;
210 //! Device buffer for intermediate positions (maintained internally)
211 DeviceBuffer<float3> d_xp_;
212 //! Number of elements in shifted coordinates buffer
214 //! Allocation size for the shifted coordinates buffer
215 int numXpAlloc_ = -1;
218 //! 1/mass for all atoms (GPU)
219 DeviceBuffer<real> d_inverseMasses_;
220 //! Number of elements in reciprocal masses buffer
221 int numInverseMasses_ = -1;
222 //! Allocation size for the reciprocal masses buffer
223 int numInverseMassesAlloc_ = -1;
225 //! Leap-Frog integrator
226 std::unique_ptr<LeapFrogGpu> integrator_;
227 //! LINCS GPU object to use for non-water constraints
228 std::unique_ptr<LincsGpu> lincsGpu_;
229 //! SETTLE GPU object for water constrains
230 std::unique_ptr<SettleGpu> settleGpu_;
232 //! An pointer to the event to indicate when the update of coordinates is complete
233 GpuEventSynchronizer* coordinatesReady_;
234 //! The wallclock counter
235 gmx_wallcycle* wcycle_ = nullptr;
240 #endif // GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H