Apply clang-format-11
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_gpu_impl.h
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 Declares GPU implementation class for update and constraints.
38  *
39  * This header file is needed to include from both the device-side
40  * kernels file, and the host-side management code.
41  *
42  * \author Artem Zhmurov <zhmurov@gmail.com>
43  *
44  * \ingroup module_mdlib
45  */
46 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H
47 #define GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H
48
49 #include "gmxpre.h"
50
51 #include "gromacs/mdlib/leapfrog_gpu.h"
52 #include "gromacs/mdlib/lincs_gpu.h"
53 #include "gromacs/mdlib/settle_gpu.h"
54 #include "gromacs/mdlib/update_constrain_gpu.h"
55 #include "gromacs/mdtypes/inputrec.h"
56
57 class GpuEventSynchronizer;
58
59 namespace gmx
60 {
61
62 /*! \internal \brief Class with interfaces and data for GPU version of Update-Constraint. */
63 class UpdateConstrainGpu::Impl
64 {
65
66 public:
67     /*! \brief Create Update-Constrain object.
68      *
69      * The constructor is given a non-nullptr \p deviceStream, in which all the update and constrain
70      * routines are executed. \p xUpdatedOnDevice should mark the completion of all kernels that
71      * modify coordinates. The event is maintained outside this class and also passed to all (if
72      * any) consumers of the updated coordinates. The \p xUpdatedOnDevice also can not be a nullptr
73      * because the markEvent(...) method is called unconditionally.
74      *
75      * \param[in] ir                  Input record data: LINCS takes number of iterations and order of
76      *                                projection from it.
77      * \param[in] mtop                Topology of the system: SETTLE gets the masses for O and H atoms
78      *                                and target O-H and H-H distances from this object.
79      * \param[in] numTempScaleValues  Number of temperature scaling groups. Set zero for no temperature coupling.
80      * \param[in] deviceContext       GPU device context.
81      * \param[in] deviceStream        GPU stream to use.
82      * \param[in] xUpdatedOnDevice    The event synchronizer to use to mark that
83      *                                update is done on the GPU.
84      * \param[in] wcycle              The wallclock counter
85      */
86     Impl(const t_inputrec&     ir,
87          const gmx_mtop_t&     mtop,
88          int                   numTempScaleValues,
89          const DeviceContext&  deviceContext,
90          const DeviceStream&   deviceStream,
91          GpuEventSynchronizer* xUpdatedOnDevice,
92          gmx_wallcycle*        wcycle);
93
94     ~Impl();
95
96     /*! \brief Integrate
97      *
98      * Integrates the equation of motion using Leap-Frog algorithm and applies
99      * LINCS and SETTLE constraints.
100      * If computeVirial is true, constraints virial is written at the provided pointer.
101      * doTempCouple should be true if:
102      *   1. The temperature coupling is enabled.
103      *   2. This is the temperature coupling step.
104      * Parameters virial/lambdas can be nullptr if computeVirial/doTempCouple are false.
105      *
106      * \param[in]  fReadyOnDevice           Event synchronizer indicating that the forces are ready in
107      *                                      the device memory.
108      * \param[in]  dt                       Timestep.
109      * \param[in]  updateVelocities         If the velocities should be constrained.
110      * \param[in]  computeVirial            If virial should be updated.
111      * \param[out] virial                   Place to save virial tensor.
112      * \param[in]  doTemperatureScaling     If velocities should be scaled for temperature coupling.
113      * \param[in]  tcstat                   Temperature coupling data.
114      * \param[in]  doParrinelloRahman       If current step is a Parrinello-Rahman pressure coupling step.
115      * \param[in]  dtPressureCouple         Period between pressure coupling steps.
116      * \param[in]  prVelocityScalingMatrix  Parrinello-Rahman velocity scaling matrix.
117      */
118     void integrate(GpuEventSynchronizer*             fReadyOnDevice,
119                    real                              dt,
120                    bool                              updateVelocities,
121                    bool                              computeVirial,
122                    tensor                            virial,
123                    bool                              doTemperatureScaling,
124                    gmx::ArrayRef<const t_grp_tcstat> tcstat,
125                    bool                              doParrinelloRahman,
126                    float                             dtPressureCouple,
127                    const matrix                      prVelocityScalingMatrix);
128
129     /*! \brief Scale coordinates on the GPU for the pressure coupling.
130      *
131      * After pressure coupling step, the box size may change. Hence, the coordinates should be
132      * scaled so that all the particles fit in the new box.
133      *
134      * \param[in] scalingMatrix Coordinates scaling matrix.
135      */
136     void scaleCoordinates(const matrix scalingMatrix);
137
138     /*! \brief Scale velocities on the GPU for the pressure coupling.
139      *
140      * After pressure coupling step, the box size may change. In the C-Rescale algorithm, velocities should be scaled.
141      *
142      * \param[in] scalingMatrix Velocities scaling matrix.
143      */
144     void scaleVelocities(const matrix scalingMatrix);
145
146     /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
147      *
148      * \param[in,out]  d_x            Device buffer with coordinates.
149      * \param[in,out]  d_v            Device buffer with velocities.
150      * \param[in]      d_f            Device buffer with forces.
151      * \param[in] idef                System topology
152      * \param[in] md                  Atoms data.
153      */
154     void set(DeviceBuffer<Float3>          d_x,
155              DeviceBuffer<Float3>          d_v,
156              const DeviceBuffer<Float3>    d_f,
157              const InteractionDefinitions& idef,
158              const t_mdatoms&              md);
159
160     /*! \brief
161      * Update PBC data.
162      *
163      * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
164      *
165      * \param[in] pbcType The type of the periodic boundary.
166      * \param[in] box     The periodic boundary box matrix.
167      */
168     void setPbc(PbcType pbcType, const matrix box);
169
170     /*! \brief Return the synchronizer associated with the event indicated that the coordinates are ready on the device.
171      */
172     GpuEventSynchronizer* getCoordinatesReadySync();
173
174     /*! \brief
175      * Returns whether the maximum number of coupled constraints is supported
176      * by the GPU LINCS code.
177      *
178      * \param[in] mtop The molecular topology
179      */
180     static bool isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop);
181
182 private:
183     //! GPU context object
184     const DeviceContext& deviceContext_;
185     //! GPU stream
186     const DeviceStream& deviceStream_;
187     //! GPU kernel launch config
188     KernelLaunchConfig coordinateScalingKernelLaunchConfig_;
189
190     //! Periodic boundary data
191     PbcAiuc pbcAiuc_;
192
193     //! Number of atoms
194     int numAtoms_;
195
196     //! Local copy of the pointer to the device positions buffer
197     DeviceBuffer<Float3> d_x_;
198     //! Local copy of the pointer to the device velocities buffer
199     DeviceBuffer<Float3> d_v_;
200     //! Local copy of the pointer to the device forces buffer
201     DeviceBuffer<Float3> d_f_;
202
203     //! Device buffer for intermediate positions (maintained internally)
204     DeviceBuffer<Float3> d_xp_;
205     //! Number of elements in shifted coordinates buffer
206     int numXp_ = -1;
207     //! Allocation size for the shifted coordinates buffer
208     int numXpAlloc_ = -1;
209
210
211     //! 1/mass for all atoms (GPU)
212     DeviceBuffer<real> d_inverseMasses_;
213     //! Number of elements in reciprocal masses buffer
214     int numInverseMasses_ = -1;
215     //! Allocation size for the reciprocal masses buffer
216     int numInverseMassesAlloc_ = -1;
217
218     //! Leap-Frog integrator
219     std::unique_ptr<LeapFrogGpu> integrator_;
220     //! LINCS GPU object to use for non-water constraints
221     std::unique_ptr<LincsGpu> lincsGpu_;
222     //! SETTLE GPU object for water constrains
223     std::unique_ptr<SettleGpu> settleGpu_;
224
225     //! An pointer to the event to indicate when the update of coordinates is complete
226     GpuEventSynchronizer* coordinatesReady_;
227     //! The wallclock counter
228     gmx_wallcycle* wcycle_ = nullptr;
229 };
230
231 /*! \brief Scaling matrix struct.
232  *
233  * \todo Should be generalized.
234  */
235 struct ScalingMatrix
236 {
237     ScalingMatrix(const matrix m) :
238         xx(m[XX][XX]), yy(m[YY][YY]), zz(m[ZZ][ZZ]), yx(m[YY][XX]), zx(m[ZZ][XX]), zy(m[ZZ][YY])
239     {
240     }
241     float xx, yy, zz, yx, zx, zy;
242 };
243
244 } // namespace gmx
245
246 #endif // GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H