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 /*! \defgroup module_listed_forces Interactions between lists of particles
37 * \ingroup group_mdrun
39 * \brief Handles computing energies and forces for listed
42 * Located here is the code for
43 * - computing energies and forces for interactions between a small
44 number of particles, e.g bonds, position restraints and listed
45 non-bonded interactions (e.g. 1-4).
46 * - high-level functions used by mdrun for computing a set of such
48 * - managing thread-wise decomposition, thread-local buffer output,
49 and reduction of output data across threads.
51 * \author Mark Abraham <mark.j.abraham@gmail.com>
52 * \author Berk Hess <hess@kth.se>
55 /*! \libinternal \file
57 * \brief This file contains declarations of high-level functions used
58 * by mdrun to compute energies and forces for listed interactions.
60 * Clients of libgromacs that want to evaluate listed interactions
61 * should call functions declared here.
63 * \author Mark Abraham <mark.j.abraham@gmail.com>
64 * \author Berk Hess <hess@kth.se>
67 * \ingroup module_listed_forces
69 #ifndef GMX_LISTED_FORCES_LISTED_FORCES_H
70 #define GMX_LISTED_FORCES_LISTED_FORCES_H
77 #include "gromacs/math/vectypes.h"
78 #include "gromacs/topology/idef.h"
79 #include "gromacs/topology/ifunc.h"
80 #include "gromacs/utility/classhelpers.h"
82 struct bonded_threading_t;
83 struct gmx_enerdata_t;
84 struct gmx_ffparams_t;
85 struct gmx_grppairener_t;
86 struct gmx_localtop_t;
87 struct gmx_multisim_t;
106 class ArrayRefWithPadding;
109 //! Type of CPU function to compute a bonded interaction.
110 using BondedFunction = real (*)(int nbonds,
111 const t_iatom iatoms[],
112 const t_iparams iparams[],
118 gmx::ArrayRef<real> dvdlambda,
121 t_disresdata* disresdata,
122 t_oriresdata* oriresdata,
125 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
126 BondedFunction bondedFunction(int ftype);
129 * \brief Class for calculating listed interactions, uses OpenMP parallelization
131 * Listed interactions can be divided over multiple instances of ListedForces
132 * using the selection flags passed to the constructor.
137 //! Enum for selecting groups of listed interaction types
138 enum class InteractionGroup : int
140 Pairs, //!< Pair interactions
141 Dihedrals, //!< Dihedrals, including cmap
143 Rest, //!< All listed interactions that are not any of the above
144 Count //!< The number of items above
147 //! Type for specifying selections of groups of interaction types
148 using InteractionSelection = std::bitset<static_cast<int>(InteractionGroup::Count)>;
150 //! Returns a selection with all listed interaction types selected
151 static InteractionSelection interactionSelectionAll()
153 InteractionSelection is;
157 /*! \brief Constructor
159 * \param[in] ffparams The force field parameters
160 * \param[in] numEnergyGroups The number of energy groups, used for storage of pair energies
161 * \param[in] numThreads The number of threads used for computed listed interactions
162 * \param[in] interactionSelection Select of interaction groups through bits set
163 * \param[in] fplog Log file for printing env.var. override, can be nullptr
165 ListedForces(const gmx_ffparams_t& ffparams,
168 InteractionSelection interactionSelection,
171 //! Move constructor, default, but in the source file to hide implementation classes
172 ListedForces(ListedForces&& o) noexcept;
174 //! Destructor which is actually default but in the source file to hide implementation classes
177 /*! \brief Copy the listed interactions from \p idef and set up the thread parallelization
179 * \param[in] domainIdef Interaction definitions for all listed interactions to be computed on this domain/rank
180 * \param[in] numAtomsForce Force are, potentially, computed for atoms 0 to \p numAtomsForce
181 * \param[in] useGpu Whether a GPU is used to compute (part of) the listed interactions
183 void setup(const InteractionDefinitions& domainIdef, int numAtomsForce, bool useGpu);
185 /*! \brief Do all aspects of energy and force calculations for mdrun
186 * on the set of listed interactions
188 * xWholeMolecules only needs to contain whole molecules when orientation
189 * restraints need to be computed and can be empty otherwise.
191 void calculate(struct gmx_wallcycle* wcycle,
193 const t_lambda* fepvals,
195 const gmx_multisim_t* ms,
196 gmx::ArrayRefWithPadding<const gmx::RVec> coordinates,
197 gmx::ArrayRef<const gmx::RVec> xWholeMolecules,
199 const history_t* hist,
200 gmx::ForceOutputs* forceOutputs,
201 const t_forcerec* fr,
202 const struct t_pbc* pbc,
203 gmx_enerdata_t* enerd,
205 gmx::ArrayRef<const real> lambda,
207 int* global_atom_index,
208 const gmx::StepWorkload& stepWork);
210 //! Returns whether bonded interactions are assigned to the CPU
211 bool haveCpuBondeds() const;
213 /*! \brief Returns whether listed forces are computed on the CPU
215 * NOTE: the current implementation returns true if there are position restraints
216 * or any bonded interactions computed on the CPU.
218 bool haveCpuListedForces(const t_fcdata& fcdata) const;
220 //! Returns true if there are position, distance or orientation restraints
221 bool haveRestraints(const t_fcdata& fcdata) const;
224 //! Pointer to the interaction definitions
225 InteractionDefinitions const* idef_ = nullptr;
226 //! Interaction defintions used for storing selections
227 InteractionDefinitions idefSelection_;
228 //! Thread parallelization setup, unique_ptr to avoid declaring bonded_threading_t
229 std::unique_ptr<bonded_threading_t> threading_;
230 //! Tells which interactions to select for computation
231 const InteractionSelection interactionSelection_;
232 //! Force buffer for free-energy forces
233 std::vector<real> forceBufferLambda_;
234 //! Shift force buffer for free-energy forces
235 std::vector<gmx::RVec> shiftForceBufferLambda_;
236 //! Temporary array for storing foreign lambda group pair energies
237 std::unique_ptr<gmx_grppairener_t> foreignEnergyGroups_;
239 GMX_DISALLOW_COPY_AND_ASSIGN(ListedForces);