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