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