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