9050039b13b7b52846746a1858d34a5be7c3320d
[alexxy/gromacs.git] / src / gromacs / mdlib / update_constrain_gpu.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 /*! \libinternal \file
36  *
37  * \brief Declaration of high-level functions of GPU implementation of update and constrain class.
38  *
39  * \author Artem Zhmurov <zhmurov@gmail.com>
40  *
41  * \ingroup module_mdlib
42  * \inlibraryapi
43  */
44 #ifndef GMX_MDLIB_UPDATE_CONSTRAIN_GPU_H
45 #define GMX_MDLIB_UPDATE_CONSTRAIN_GPU_H
46
47 #include <memory>
48
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"
53
54 class DeviceContext;
55 class DeviceStream;
56 class GpuEventSynchronizer;
57 struct gmx_mtop_t;
58 enum class PbcType : int;
59 class InteractionDefinitions;
60 struct t_inputrec;
61 struct t_mdatoms;
62 struct t_pbc;
63
64 namespace gmx
65 {
66
67 class UpdateConstrainGpu
68 {
69
70 public:
71     /*! \brief Create Update-Constrain object.
72      *
73      * The constructor is given a non-nullptr \p deviceStream, in which all the update and constrain
74      * routines are executed. \p xUpdatedOnDevice should mark the completion of all kernels that
75      * modify coordinates. The event is maintained outside this class and also passed to all (if
76      * any) consumers of the updated coordinates. The \p xUpdatedOnDevice also can not be a nullptr
77      * because the markEvent(...) method is called unconditionally.
78      *
79      * \param[in] ir                  Input record data: LINCS takes number of iterations and order of
80      *                                projection from it.
81      * \param[in] mtop                Topology of the system: SETTLE gets the masses for O and H atoms
82      *                                and target O-H and H-H distances from this object.
83      * \param[in] numTempScaleValues  Number of temperature scaling groups. Zero for no temperature scaling.
84      * \param[in] deviceContext       GPU device context.
85      * \param[in] deviceStream        GPU stream to use.
86      * \param[in] xUpdatedOnDevice    The event synchronizer to use to mark that update is done
87      *                                on the GPU.
88      * \param[in] wcycle              The wallclock counter
89      */
90     UpdateConstrainGpu(const t_inputrec&     ir,
91                        const gmx_mtop_t&     mtop,
92                        int                   numTempScaleValues,
93                        const DeviceContext&  deviceContext,
94                        const DeviceStream&   deviceStream,
95                        GpuEventSynchronizer* xUpdatedOnDevice,
96                        gmx_wallcycle*        wcycle);
97
98     ~UpdateConstrainGpu();
99
100     /*! \brief Integrate
101      *
102      * This will extract temperature scaling factors from tcstat, transform them into the plain
103      * array and call the normal integrate method.
104      *
105      * \param[in]  fReadyOnDevice           Event synchronizer indicating that the forces are
106      *                                      ready in the device memory.
107      * \param[in]  dt                       Timestep.
108      * \param[in]  updateVelocities         If the velocities should be constrained.
109      * \param[in]  computeVirial            If virial should be updated.
110      * \param[out] virial                   Place to save virial tensor.
111      * \param[in]  doTemperatureScaling     If velocities should be scaled for temperature coupling.
112      * \param[in]  tcstat                   Temperature coupling data.
113      * \param[in]  doParrinelloRahman       If current step is a Parrinello-Rahman pressure coupling step.
114      * \param[in]  dtPressureCouple         Period between pressure coupling steps.
115      * \param[in]  prVelocityScalingMatrix  Parrinello-Rahman velocity scaling matrix.
116      */
117     void integrate(GpuEventSynchronizer*             fReadyOnDevice,
118                    real                              dt,
119                    bool                              updateVelocities,
120                    bool                              computeVirial,
121                    tensor                            virial,
122                    bool                              doTemperatureScaling,
123                    gmx::ArrayRef<const t_grp_tcstat> tcstat,
124                    bool                              doParrinelloRahman,
125                    float                             dtPressureCouple,
126                    const matrix                      prVelocityScalingMatrix);
127
128     /*! \brief Scale coordinates on the GPU for the pressure coupling.
129      *
130      * After pressure coupling step, the box size may change. Hence, the coordinates should be
131      * scaled so that all the particles fit in the new box.
132      *
133      * \param[in] scalingMatrix Coordinates scaling matrix.
134      */
135     void scaleCoordinates(const matrix scalingMatrix);
136
137     /*! \brief Scale velocities on the GPU for the pressure coupling.
138      *
139      * After pressure coupling step, the box size may change. In the C-Rescale algorithm, velocities should be scaled.
140      *
141      * \param[in] scalingMatrix Velocities scaling matrix.
142      */
143     void scaleVelocities(const matrix scalingMatrix);
144
145     /*! \brief Set the pointers and update data-structures (e.g. after NB search step).
146      *
147      * \param[in,out]  d_x                 Device buffer with coordinates.
148      * \param[in,out]  d_v                 Device buffer with velocities.
149      * \param[in]      d_f                 Device buffer with forces.
150      * \param[in]      idef                System topology
151      * \param[in]      md                  Atoms data.
152      */
153     void set(DeviceBuffer<RVec>            d_x,
154              DeviceBuffer<RVec>            d_v,
155              DeviceBuffer<RVec>            d_f,
156              const InteractionDefinitions& idef,
157              const t_mdatoms&              md);
158
159     /*! \brief
160      * Update PBC data.
161      *
162      * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
163      *
164      * \param[in] pbcType The type of the periodic boundary (Xyz, NO, XY or Screw).
165      * \param[in] box     The periodic boundary box matrix.
166      */
167     void setPbc(PbcType pbcType, const matrix box);
168
169     /*! \brief Return the synchronizer associated with the event indicated that the coordinates are ready on the device.
170      */
171     GpuEventSynchronizer* getCoordinatesReadySync();
172
173     /*! \brief
174      * Returns whether the maximum number of coupled constraints is supported
175      * by the GPU LINCS code.
176      *
177      * \param[in] mtop The molecular topology
178      */
179     static bool isNumCoupledConstraintsSupported(const gmx_mtop_t& mtop);
180
181     /*! \brief
182      * Returns whether the constraints are supported by the GPU code.
183      *
184      * Currently true for CUDA, false for others.
185      */
186     static bool areConstraintsSupported();
187
188 private:
189     class Impl;
190     std::unique_ptr<Impl> impl_;
191 };
192
193 } // namespace gmx
194
195 #endif // GMX_MDLIB_UPDATE_CONSTRAIN_GPU_H