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