2ad8969019c3839b1455cc29ce33e701d6bcfcfd
[alexxy/gromacs.git] / src / gromacs / listed_forces / gpubonded.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5  * Copyright (c) 2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \libinternal \file
37  *
38  * \brief This file contains declarations of high-level functions used
39  * by mdrun to compute energies and forces for listed interactions.
40  *
41  * Clients of libgromacs that want to evaluate listed interactions
42  * should call functions declared here.
43  *
44  * \author Mark Abraham <mark.j.abraham@gmail.com>
45  *
46  * \inlibraryapi
47  * \ingroup module_listed_forces
48  */
49 #ifndef GMX_LISTED_FORCES_GPUBONDED_H
50 #define GMX_LISTED_FORCES_GPUBONDED_H
51
52 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
53 #include "gromacs/math/vectypes.h"
54 #include "gromacs/pbcutil/pbc.h"
55 #include "gromacs/topology/idef.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/classhelpers.h"
58
59 class DeviceContext;
60 class DeviceStream;
61
62 struct gmx_enerdata_t;
63 struct gmx_ffparams_t;
64 struct gmx_mtop_t;
65 struct t_inputrec;
66 struct gmx_wallcycle;
67
68
69 namespace gmx
70 {
71
72 class StepWorkload;
73
74 /*! \brief The number on bonded function types supported on GPUs */
75 static constexpr int numFTypesOnGpu = 8;
76
77 /*! \brief List of all bonded function types supported on GPUs
78  *
79  * \note This list should be in sync with the actual GPU code.
80  * \note Perturbed interactions are not supported on GPUs.
81  * \note The function types in the list are ordered on increasing value.
82  * \note Currently bonded are only supported with CUDA, not with OpenCL.
83  */
84 constexpr std::array<int, numFTypesOnGpu> fTypesOnGpu = { F_BONDS,  F_ANGLES, F_UREY_BRADLEY,
85                                                           F_PDIHS,  F_RBDIHS, F_IDIHS,
86                                                           F_PIDIHS, F_LJ14 };
87
88 /*! \brief Checks whether the GROMACS build allows to compute bonded interactions on a GPU.
89  *
90  * \param[out] error  If non-null, the diagnostic message when bondeds cannot run on a GPU.
91  *
92  * \returns true when this build can run bonded interactions on a GPU, false otherwise.
93  *
94  * \throws std::bad_alloc when out of memory.
95  */
96 bool buildSupportsGpuBondeds(std::string* error);
97
98 /*! \brief Checks whether the input system allows to compute bonded interactions on a GPU.
99  *
100  * \param[in]  ir     Input system.
101  * \param[in]  mtop   Complete system topology to search for supported interactions.
102  * \param[out] error  If non-null, the error message if the input is not supported on GPU.
103  *
104  * \returns true if PME can run on GPU with this input, false otherwise.
105  */
106 bool inputSupportsGpuBondeds(const t_inputrec& ir, const gmx_mtop_t& mtop, std::string* error);
107
108 class GpuBonded
109 {
110 public:
111     /*! \brief Construct the manager with constant data and the stream to use.
112      *
113      * \param[in] ffparams                   Force-field parameters.
114      * \param[in] electrostaticsScaleFactor  Scaling factor for the electrostatic potential
115      *                                       (Coulomb constant, multiplied by the Fudge factor).
116      * \param[in] deviceContext              GPU device context (not used in CUDA).
117      * \param[in] deviceStream               GPU device stream.
118      * \param[in] wcycle                     The wallclock counter.
119      *
120      */
121     GpuBonded(const gmx_ffparams_t& ffparams,
122               float                 electrostaticsScaleFactor,
123               const DeviceContext&  deviceContext,
124               const DeviceStream&   deviceStream,
125               gmx_wallcycle*        wcycle);
126     //! Destructor
127     ~GpuBonded();
128
129     /*! \brief Update lists of interactions from idef suitable for the GPU,
130      * using the data structures prepared for PP work.
131      *
132      * Intended to be called after each neighbour search
133      * stage. Copies the bonded interactions assigned to the GPU
134      * to device data structures, and updates device buffers that
135      * may have been updated after search.
136      *
137      * \param[in]     nbnxnAtomOrder  Mapping between rvec and NBNXM formats.
138      * \param[in]     idef            List of interactions to compute.
139      * \param[in]     xqDevice        Device buffer with coordinates and charge in xyzq-format.
140      * \param[in,out] forceDevice     Device buffer with forces.
141      * \param[in,out] fshiftDevice    Device buffer with shift forces.
142      */
143     void updateInteractionListsAndDeviceBuffers(ArrayRef<const int>           nbnxnAtomOrder,
144                                                 const InteractionDefinitions& idef,
145                                                 void*                         xqDevice,
146                                                 DeviceBuffer<RVec>            forceDevice,
147                                                 DeviceBuffer<RVec>            fshiftDevice);
148     /*! \brief
149      * Update PBC data.
150      *
151      * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
152      *
153      * \param[in] pbcType The type of the periodic boundary.
154      * \param[in] box     The periodic boundary box matrix.
155      * \param[in] canMoleculeSpanPbc  Whether one molecule can have atoms in different PBC cells.
156      */
157     void setPbc(PbcType pbcType, const matrix box, bool canMoleculeSpanPbc);
158
159     /*! \brief Returns whether there are bonded interactions
160      * assigned to the GPU
161      *
162      * \returns If the list of interaction has elements.
163      */
164     bool haveInteractions() const;
165
166     /*! \brief Launches bonded kernel on a GPU
167      *
168      * \param[in]  stepWork  Simulation step work to determine if energy/virial are to be computed on this step.
169      */
170     void launchKernel(const gmx::StepWorkload& stepWork);
171
172     /*! \brief Sets the PBC and launches bonded kernel on a GPU
173      *
174      * \param[in] pbcType The type of the periodic boundary.
175      * \param[in] box     The periodic boundary box matrix.
176      * \param[in] canMoleculeSpanPbc  Whether one molecule can have atoms in different PBC cells.
177      * \param[in] stepWork  Simulation step work to determine if energy/virial are to be computed on this step.
178      */
179     void setPbcAndlaunchKernel(PbcType                  pbcType,
180                                const matrix             box,
181                                bool                     canMoleculeSpanPbc,
182                                const gmx::StepWorkload& stepWork);
183
184     /*! \brief Launches the transfer of computed bonded energies.
185      */
186     void launchEnergyTransfer();
187
188     /*! \brief Waits on the energy transfer, and accumulates bonded energies to \c enerd.
189      *
190      * \param[in,out] enerd The energy data object to add energy terms to.
191      */
192     void waitAccumulateEnergyTerms(gmx_enerdata_t* enerd);
193
194     /*! \brief Clears the device side energy buffer
195      */
196     void clearEnergies();
197
198 private:
199     class Impl;
200     PrivateImplPointer<Impl> impl_;
201 };
202
203 } // namespace gmx
204
205 #endif