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.
35 /*! \libinternal \file
37 * \brief Declaration of high-level functions of GPU implementation of update and constrain class.
39 * \author Artem Zhmurov <zhmurov@gmail.com>
41 * \ingroup module_mdlib
44 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_GPU_H
45 #define GMX_MDLIB_UPDATE_CONSTRAIN_GPU_H
49 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
50 #include "gromacs/mdtypes/group.h"
51 #include "gromacs/timing/wallcycle.h"
52 #include "gromacs/utility/arrayref.h"
56 class GpuEventSynchronizer;
58 enum class PbcType : int;
59 class InteractionDefinitions;
67 class UpdateConstrainGpu
71 /*! \brief Create Update-Constrain object.
73 * The constructor is given a non-nullptr \p deviceStream, in which all the update and constrain
74 * routines are executed.
76 * \param[in] ir Input record data: LINCS takes number of iterations and order of
78 * \param[in] mtop Topology of the system: SETTLE gets the masses for O and H atoms
79 * and target O-H and H-H distances from this object.
80 * \param[in] numTempScaleValues Number of temperature scaling groups. Zero for no temperature scaling.
81 * \param[in] deviceContext GPU device context.
82 * \param[in] deviceStream GPU stream to use.
83 * \param[in] wcycle The wallclock counter
85 UpdateConstrainGpu(const t_inputrec& ir,
86 const gmx_mtop_t& mtop,
87 int numTempScaleValues,
88 const DeviceContext& deviceContext,
89 const DeviceStream& deviceStream,
90 gmx_wallcycle* wcycle);
92 ~UpdateConstrainGpu();
96 * This will extract temperature scaling factors from tcstat, transform them into the plain
97 * array and call the normal integrate method.
99 * \param[in] fReadyOnDevice Event synchronizer indicating that the forces are
100 * ready in the device memory.
101 * \param[in] dt Timestep.
102 * \param[in] updateVelocities If the velocities should be constrained.
103 * \param[in] computeVirial If virial should be updated.
104 * \param[out] virial Place to save virial tensor.
105 * \param[in] doTemperatureScaling If velocities should be scaled for temperature coupling.
106 * \param[in] tcstat Temperature coupling data.
107 * \param[in] doParrinelloRahman If current step is a Parrinello-Rahman pressure coupling step.
108 * \param[in] dtPressureCouple Period between pressure coupling steps.
109 * \param[in] prVelocityScalingMatrix Parrinello-Rahman velocity scaling matrix.
111 void integrate(GpuEventSynchronizer* fReadyOnDevice,
113 bool updateVelocities,
116 bool doTemperatureScaling,
117 gmx::ArrayRef<const t_grp_tcstat> tcstat,
118 bool doParrinelloRahman,
119 float dtPressureCouple,
120 const matrix prVelocityScalingMatrix);
122 /*! \brief Scale coordinates on the GPU for the pressure coupling.
124 * After pressure coupling step, the box size may change. Hence, the coordinates should be
125 * scaled so that all the particles fit in the new box.
127 * \param[in] scalingMatrix Coordinates scaling matrix.
129 void scaleCoordinates(const matrix scalingMatrix);
131 /*! \brief Scale velocities on the GPU for the pressure coupling.
133 * After pressure coupling step, the box size may change. In the C-Rescale algorithm, velocities should be scaled.
135 * \param[in] scalingMatrix Velocities scaling matrix.
137 void scaleVelocities(const matrix scalingMatrix);
139 /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
141 * \param[in,out] d_x Device buffer with coordinates.
142 * \param[in,out] d_v Device buffer with velocities.
143 * \param[in] d_f Device buffer with forces.
144 * \param[in] idef System topology
145 * \param[in] md Atoms data.
147 void set(DeviceBuffer<RVec> d_x,
148 DeviceBuffer<RVec> d_v,
149 DeviceBuffer<RVec> d_f,
150 const InteractionDefinitions& idef,
151 const t_mdatoms& md);
156 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
158 * \param[in] pbcType The type of the periodic boundary (Xyz, NO, XY or Screw).
159 * \param[in] box The periodic boundary box matrix.
161 void setPbc(PbcType pbcType, const matrix box);
163 /*! \brief Return the synchronizer associated with the event that indicates
164 * that the coordinates are ready on the device.
166 GpuEventSynchronizer* xUpdatedOnDeviceEvent();
169 * Returns whether the maximum number of coupled constraints is supported
170 * by the GPU LINCS code.
172 * \param[in] mtop The molecular topology
174 static bool isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop);
177 * Returns whether the constraints are supported by the GPU code.
179 * Currently true for CUDA, false for others.
181 static bool areConstraintsSupported();
185 std::unique_ptr<Impl> impl_;
190 #endif // GMX_MDLIB_UPDATE_CONSTRAIN_GPU_H