Set temperature scaling groups upon construction of integrator
[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,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 /*! \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
70      * modify coordinates. The event is maintained outside this class and also passed to all (if
71      * any) consumers of the updated coordinates. The \p xUpdatedOnDevice also can not be a nullptr
72      * because the 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] numTempScaleValues  Number of temperature scaling groups. Set zero for no temperature coupling.
79      * \param[in] deviceContext       GPU device context.
80      * \param[in] deviceStream        GPU stream to use.
81      * \param[in] xUpdatedOnDevice    The event synchronizer to use to mark that
82      *                                update is done on the GPU.
83      * \param[in] wcycle              The wallclock counter
84      */
85     Impl(const t_inputrec&     ir,
86          const gmx_mtop_t&     mtop,
87          int                   numTempScaleValues,
88          const DeviceContext&  deviceContext,
89          const DeviceStream&   deviceStream,
90          GpuEventSynchronizer* xUpdatedOnDevice,
91          gmx_wallcycle*        wcycle);
92
93     ~Impl();
94
95     /*! \brief Integrate
96      *
97      * Integrates the equation of motion using Leap-Frog algorithm and applies
98      * LINCS and SETTLE constraints.
99      * If computeVirial is true, constraints virial is written at the provided pointer.
100      * doTempCouple should be true if:
101      *   1. The temperature coupling is enabled.
102      *   2. This is the temperature coupling step.
103      * Parameters virial/lambdas can be nullptr if computeVirial/doTempCouple are false.
104      *
105      * \param[in]  fReadyOnDevice           Event synchronizer indicating that the forces are ready in
106      *                                      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              const 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.
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 private:
182     //! GPU context object
183     const DeviceContext& deviceContext_;
184     //! GPU stream
185     const DeviceStream& deviceStream_;
186     //! GPU kernel launch config
187     KernelLaunchConfig coordinateScalingKernelLaunchConfig_;
188
189     //! Periodic boundary data
190     PbcAiuc pbcAiuc_;
191
192     //! Number of atoms
193     int numAtoms_;
194
195     //! Local copy of the pointer to the device positions buffer
196     float3* d_x_;
197     //! Local copy of the pointer to the device velocities buffer
198     float3* d_v_;
199     //! Local copy of the pointer to the device forces buffer
200     float3* d_f_;
201
202     //! Device buffer for intermediate positions (maintained internally)
203     float3* d_xp_;
204     //! Number of elements in shifted coordinates buffer
205     int numXp_ = -1;
206     //! Allocation size for the shifted coordinates buffer
207     int numXpAlloc_ = -1;
208
209
210     //! 1/mass for all atoms (GPU)
211     real* d_inverseMasses_;
212     //! Number of elements in reciprocal masses buffer
213     int numInverseMasses_ = -1;
214     //! Allocation size for the reciprocal masses buffer
215     int numInverseMassesAlloc_ = -1;
216
217     //! Leap-Frog integrator
218     std::unique_ptr<LeapFrogGpu> integrator_;
219     //! LINCS GPU object to use for non-water constraints
220     std::unique_ptr<LincsGpu> lincsGpu_;
221     //! SETTLE GPU object for water constrains
222     std::unique_ptr<SettleGpu> settleGpu_;
223
224     //! An pointer to the event to indicate when the update of coordinates is complete
225     GpuEventSynchronizer* coordinatesReady_;
226     //! The wallclock counter
227     gmx_wallcycle* wcycle_ = nullptr;
228 };
229
230 } // namespace gmx
231
232 #endif // GMX_MDLIB_UPDATE_CONSTRAIN_GPU_IMPL_H