2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015,2016,2017,2018 by the GROMACS development team.
5 * Copyright (c) 2019,2020,2021, 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.
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.
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.
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.
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.
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.
36 /*! \libinternal \file
38 * \brief This file contains declarations of high-level functions used
39 * by mdrun to compute energies and forces for listed interactions.
41 * Clients of libgromacs that want to evaluate listed interactions
42 * should call functions declared here.
44 * \author Mark Abraham <mark.j.abraham@gmail.com>
47 * \ingroup module_listed_forces
49 #ifndef GMX_LISTED_FORCES_LISTED_FORCES_GPU_H
50 #define GMX_LISTED_FORCES_LISTED_FORCES_GPU_H
54 #include "gromacs/gpu_utils/devicebuffer_datatype.h"
55 #include "gromacs/math/vectypes.h"
56 #include "gromacs/pbcutil/pbc.h"
57 #include "gromacs/topology/idef.h"
62 struct gmx_enerdata_t;
63 struct gmx_ffparams_t;
76 /*! \brief The number on bonded function types supported on GPUs */
77 static constexpr int numFTypesOnGpu = 8;
79 /*! \brief List of all bonded function types supported on GPUs
81 * \note This list should be in sync with the actual GPU code.
82 * \note Perturbed interactions are not supported on GPUs.
83 * \note The function types in the list are ordered on increasing value.
84 * \note Currently bonded are only supported with CUDA, not with OpenCL.
86 constexpr std::array<int, numFTypesOnGpu> fTypesOnGpu = { F_BONDS, F_ANGLES, F_UREY_BRADLEY,
87 F_PDIHS, F_RBDIHS, F_IDIHS,
90 /*! \brief Checks whether the GROMACS build allows to compute bonded interactions on a GPU.
92 * \param[out] error If non-null, the diagnostic message when bondeds cannot run on a GPU.
94 * \returns true when this build can run bonded interactions on a GPU, false otherwise.
96 * \throws std::bad_alloc when out of memory.
98 bool buildSupportsListedForcesGpu(std::string* error);
100 /*! \brief Checks whether the input system allows to compute bonded interactions on a GPU.
102 * \param[in] ir Input system.
103 * \param[in] mtop Complete system topology to search for supported interactions.
104 * \param[out] error If non-null, the error message if the input is not supported on GPU.
106 * \returns true if PME can run on GPU with this input, false otherwise.
108 bool inputSupportsListedForcesGpu(const t_inputrec& ir, const gmx_mtop_t& mtop, std::string* error);
110 class ListedForcesGpu
113 /*! \brief Construct the manager with constant data and the stream to use.
115 * \param[in] ffparams Force-field parameters.
116 * \param[in] electrostaticsScaleFactor Scaling factor for the electrostatic potential
117 * (Coulomb constant, multiplied by the Fudge factor).
118 * \param[in] deviceContext GPU device context (not used in CUDA).
119 * \param[in] deviceStream GPU device stream.
120 * \param[in] wcycle The wallclock counter.
123 ListedForcesGpu(const gmx_ffparams_t& ffparams,
124 float electrostaticsScaleFactor,
125 const DeviceContext& deviceContext,
126 const DeviceStream& deviceStream,
127 gmx_wallcycle* wcycle);
131 /*! \brief Update lists of interactions from idef suitable for the GPU,
132 * using the data structures prepared for PP work.
134 * Intended to be called after each neighbour search
135 * stage. Copies the bonded interactions assigned to the GPU
136 * to device data structures, and updates device buffers that
137 * may have been updated after search.
139 * \param[in] nbnxnAtomOrder Mapping between rvec and NBNXM formats.
140 * \param[in] idef List of interactions to compute.
141 * \param[in] xqDevice Device buffer with coordinates and charge in xyzq-format.
142 * \param[in,out] forceDevice Device buffer with forces.
143 * \param[in,out] fshiftDevice Device buffer with shift forces.
145 void updateInteractionListsAndDeviceBuffers(ArrayRef<const int> nbnxnAtomOrder,
146 const InteractionDefinitions& idef,
148 DeviceBuffer<RVec> forceDevice,
149 DeviceBuffer<RVec> fshiftDevice);
153 * Converts PBC data from t_pbc into the PbcAiuc format and stores the latter.
155 * \param[in] pbcType The type of the periodic boundary.
156 * \param[in] box The periodic boundary box matrix.
157 * \param[in] canMoleculeSpanPbc Whether one molecule can have atoms in different PBC cells.
159 void setPbc(PbcType pbcType, const matrix box, bool canMoleculeSpanPbc);
161 /*! \brief Returns whether there are bonded interactions
162 * assigned to the GPU
164 * \returns If the list of interaction has elements.
166 bool haveInteractions() const;
168 /*! \brief Launches bonded kernel on a GPU
170 * \param[in] stepWork Simulation step work to determine if energy/virial are to be computed on this step.
172 void launchKernel(const gmx::StepWorkload& stepWork);
174 /*! \brief Sets the PBC and launches bonded kernel on a GPU
176 * \param[in] pbcType The type of the periodic boundary.
177 * \param[in] box The periodic boundary box matrix.
178 * \param[in] canMoleculeSpanPbc Whether one molecule can have atoms in different PBC cells.
179 * \param[in] stepWork Simulation step work to determine if energy/virial are to be computed on this step.
181 void setPbcAndlaunchKernel(PbcType pbcType,
183 bool canMoleculeSpanPbc,
184 const gmx::StepWorkload& stepWork);
186 /*! \brief Launches the transfer of computed bonded energies.
188 void launchEnergyTransfer();
190 /*! \brief Waits on the energy transfer, and accumulates bonded energies to \c enerd.
192 * \param[in,out] enerd The energy data object to add energy terms to.
194 void waitAccumulateEnergyTerms(gmx_enerdata_t* enerd);
196 /*! \brief Clears the device side energy buffer
198 void clearEnergies();
202 std::unique_ptr<Impl> impl_;
207 #endif // GMX_LISTED_FORCES_LISTED_FORCES_GPU_H