1cfeed429d2e92e84a4019d9366ab9361cc81fd7
[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/classhelpers.h"
57
58 class DeviceContext;
59 class DeviceStream;
60
61 struct gmx_enerdata_t;
62 struct gmx_ffparams_t;
63 struct gmx_mtop_t;
64 struct t_inputrec;
65 struct gmx_wallcycle;
66
67
68 namespace gmx
69 {
70
71 template<typename>
72 class ArrayRef;
73 class StepWorkload;
74
75 /*! \brief The number on bonded function types supported on GPUs */
76 static constexpr int numFTypesOnGpu = 8;
77
78 /*! \brief List of all bonded function types supported on GPUs
79  *
80  * \note This list should be in sync with the actual GPU code.
81  * \note Perturbed interactions are not supported on GPUs.
82  * \note The function types in the list are ordered on increasing value.
83  * \note Currently bonded are only supported with CUDA, not with OpenCL.
84  */
85 constexpr std::array<int, numFTypesOnGpu> fTypesOnGpu = { F_BONDS,  F_ANGLES, F_UREY_BRADLEY,
86                                                           F_PDIHS,  F_RBDIHS, F_IDIHS,
87                                                           F_PIDIHS, F_LJ14 };
88
89 /*! \brief Checks whether the GROMACS build allows to compute bonded interactions on a GPU.
90  *
91  * \param[out] error  If non-null, the diagnostic message when bondeds cannot run on a GPU.
92  *
93  * \returns true when this build can run bonded interactions on a GPU, false otherwise.
94  *
95  * \throws std::bad_alloc when out of memory.
96  */
97 bool buildSupportsGpuBondeds(std::string* error);
98
99 /*! \brief Checks whether the input system allows to compute bonded interactions on a GPU.
100  *
101  * \param[in]  ir     Input system.
102  * \param[in]  mtop   Complete system topology to search for supported interactions.
103  * \param[out] error  If non-null, the error message if the input is not supported on GPU.
104  *
105  * \returns true if PME can run on GPU with this input, false otherwise.
106  */
107 bool inputSupportsGpuBondeds(const t_inputrec& ir, const gmx_mtop_t& mtop, std::string* error);
108
109 class GpuBonded
110 {
111 public:
112     /*! \brief Construct the manager with constant data and the stream to use.
113      *
114      * \param[in] ffparams                   Force-field parameters.
115      * \param[in] electrostaticsScaleFactor  Scaling factor for the electrostatic potential
116      *                                       (Coulomb constant, multiplied by the Fudge factor).
117      * \param[in] deviceContext              GPU device context (not used in CUDA).
118      * \param[in] deviceStream               GPU device stream.
119      * \param[in] wcycle                     The wallclock counter.
120      *
121      */
122     GpuBonded(const gmx_ffparams_t& ffparams,
123               float                 electrostaticsScaleFactor,
124               const DeviceContext&  deviceContext,
125               const DeviceStream&   deviceStream,
126               gmx_wallcycle*        wcycle);
127     //! Destructor
128     ~GpuBonded();
129
130     /*! \brief Update lists of interactions from idef suitable for the GPU,
131      * using the data structures prepared for PP work.
132      *
133      * Intended to be called after each neighbour search
134      * stage. Copies the bonded interactions assigned to the GPU
135      * to device data structures, and updates device buffers that
136      * may have been updated after search.
137      *
138      * \param[in]     nbnxnAtomOrder  Mapping between rvec and NBNXM formats.
139      * \param[in]     idef            List of interactions to compute.
140      * \param[in]     xqDevice        Device buffer with coordinates and charge in xyzq-format.
141      * \param[in,out] forceDevice     Device buffer with forces.
142      * \param[in,out] fshiftDevice    Device buffer with shift forces.
143      */
144     void updateInteractionListsAndDeviceBuffers(ArrayRef<const int>           nbnxnAtomOrder,
145                                                 const InteractionDefinitions& idef,
146                                                 void*                         xqDevice,
147                                                 DeviceBuffer<RVec>            forceDevice,
148                                                 DeviceBuffer<RVec>            fshiftDevice);
149     /*! \brief
150      * Update PBC data.
151      *
152      * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
153      *
154      * \param[in] pbcType The type of the periodic boundary.
155      * \param[in] box     The periodic boundary box matrix.
156      * \param[in] canMoleculeSpanPbc  Whether one molecule can have atoms in different PBC cells.
157      */
158     void setPbc(PbcType pbcType, const matrix box, bool canMoleculeSpanPbc);
159
160     /*! \brief Returns whether there are bonded interactions
161      * assigned to the GPU
162      *
163      * \returns If the list of interaction has elements.
164      */
165     bool haveInteractions() const;
166
167     /*! \brief Launches bonded kernel on a GPU
168      *
169      * \param[in]  stepWork  Simulation step work to determine if energy/virial are to be computed on this step.
170      */
171     void launchKernel(const gmx::StepWorkload& stepWork);
172
173     /*! \brief Sets the PBC and launches bonded kernel on a GPU
174      *
175      * \param[in] pbcType The type of the periodic boundary.
176      * \param[in] box     The periodic boundary box matrix.
177      * \param[in] canMoleculeSpanPbc  Whether one molecule can have atoms in different PBC cells.
178      * \param[in] stepWork  Simulation step work to determine if energy/virial are to be computed on this step.
179      */
180     void setPbcAndlaunchKernel(PbcType                  pbcType,
181                                const matrix             box,
182                                bool                     canMoleculeSpanPbc,
183                                const gmx::StepWorkload& stepWork);
184
185     /*! \brief Launches the transfer of computed bonded energies.
186      */
187     void launchEnergyTransfer();
188
189     /*! \brief Waits on the energy transfer, and accumulates bonded energies to \c enerd.
190      *
191      * \param[in,out] enerd The energy data object to add energy terms to.
192      */
193     void waitAccumulateEnergyTerms(gmx_enerdata_t* enerd);
194
195     /*! \brief Clears the device side energy buffer
196      */
197     void clearEnergies();
198
199 private:
200     class Impl;
201     PrivateImplPointer<Impl> impl_;
202 };
203
204 } // namespace gmx
205
206 #endif