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, 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/basedefinitions.h"
81 #include "gromacs/utility/classhelpers.h"
83 struct bonded_threading_t;
84 struct gmx_enerdata_t;
85 struct gmx_ffparams_t;
86 struct gmx_grppairener_t;
87 struct gmx_localtop_t;
88 struct gmx_multisim_t;
105 class ArrayRefWithPadding;
108 //! Type of CPU function to compute a bonded interaction.
109 using BondedFunction = real (*)(int nbonds,
110 const t_iatom iatoms[],
111 const t_iparams iparams[],
122 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
123 BondedFunction bondedFunction(int ftype);
126 * \brief Class for calculating listed interactions, uses OpenMP parallelization
128 * Listed interactions can be divided over multiple instances of ListedForces
129 * using the selection flags passed to the constructor.
134 //! Enum for selecting groups of listed interaction types
135 enum class InteractionGroup : int
137 Pairs, //!< Pair interactions
138 Dihedrals, //!< Dihedrals, including cmap
140 Rest, //!< All listed interactions that are not any of the above
141 Count //!< The number of items above
144 //! Type for specifying selections of groups of interaction types
145 using InteractionSelection = std::bitset<static_cast<int>(InteractionGroup::Count)>;
147 //! Returns a selection with all listed interaction types selected
148 static InteractionSelection interactionSelectionAll()
150 InteractionSelection is;
154 /*! \brief Constructor
156 * \param[in] ffparams The force field parameters
157 * \param[in] numEnergyGroups The number of energy groups, used for storage of pair energies
158 * \param[in] numThreads The number of threads used for computed listed interactions
159 * \param[in] interactionSelection Select of interaction groups through bits set
160 * \param[in] fplog Log file for printing env.var. override, can be nullptr
162 ListedForces(const gmx_ffparams_t& ffparams,
165 InteractionSelection interactionSelection,
168 //! Move constructor, default, but in the source file to hide implementation classes
169 ListedForces(ListedForces&& o) noexcept;
171 //! Destructor which is actually default but in the source file to hide implementation classes
174 /*! \brief Copy the listed interactions from \p idef and set up the thread parallelization
176 * \param[in] domainIdef Interaction definitions for all listed interactions to be computed on this domain/rank
177 * \param[in] numAtomsForce Force are, potentially, computed for atoms 0 to \p numAtomsForce
178 * \param[in] useGpu Whether a GPU is used to compute (part of) the listed interactions
180 void setup(const InteractionDefinitions& domainIdef, int numAtomsForce, bool useGpu);
182 /*! \brief Do all aspects of energy and force calculations for mdrun
183 * on the set of listed interactions
185 * xWholeMolecules only needs to contain whole molecules when orientation
186 * restraints need to be computed and can be empty otherwise.
188 void calculate(struct gmx_wallcycle* wcycle,
190 const t_lambda* fepvals,
192 const gmx_multisim_t* ms,
193 gmx::ArrayRefWithPadding<const gmx::RVec> coordinates,
194 gmx::ArrayRef<const gmx::RVec> xWholeMolecules,
197 gmx::ForceOutputs* forceOutputs,
198 const t_forcerec* fr,
199 const struct t_pbc* pbc,
200 gmx_enerdata_t* enerd,
204 int* global_atom_index,
205 const gmx::StepWorkload& stepWork);
207 //! Returns whether bonded interactions are assigned to the CPU
208 bool haveCpuBondeds() const;
210 /*! \brief Returns whether listed forces are computed on the CPU
212 * NOTE: the current implementation returns true if there are position restraints
213 * or any bonded interactions computed on the CPU.
215 bool haveCpuListedForces(const t_fcdata& fcdata) const;
217 //! Returns true if there are position, distance or orientation restraints
218 bool haveRestraints(const t_fcdata& fcdata) const;
221 //! Pointer to the interaction definitions
222 InteractionDefinitions const* idef_ = nullptr;
223 //! Interaction defintions used for storing selections
224 InteractionDefinitions idefSelection_;
225 //! Thread parallelization setup, unique_ptr to avoid declaring bonded_threading_t
226 std::unique_ptr<bonded_threading_t> threading_;
227 //! Tells which interactions to select for computation
228 const InteractionSelection interactionSelection_;
229 //! Force buffer for free-energy forces
230 std::vector<real> forceBufferLambda_;
231 //! Shift force buffer for free-energy forces
232 std::vector<gmx::RVec> shiftForceBufferLambda_;
234 GMX_DISALLOW_COPY_AND_ASSIGN(ListedForces);