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>
54 /*! \libinternal \file
56 * \brief This file contains declarations of high-level functions used
57 * by mdrun to compute energies and forces for listed interactions.
59 * Clients of libgromacs that want to evaluate listed interactions
60 * should call functions declared here.
62 * \author Mark Abraham <mark.j.abraham@gmail.com>
65 * \ingroup module_listed_forces
67 #ifndef GMX_LISTED_FORCES_LISTED_FORCES_H
68 #define GMX_LISTED_FORCES_LISTED_FORCES_H
70 #include "gromacs/math/vectypes.h"
71 #include "gromacs/topology/ifunc.h"
72 #include "gromacs/utility/basedefinitions.h"
74 struct gmx_enerdata_t
;
75 struct gmx_grppairener_t
;
76 struct gmx_multisim_t
;
96 //! Type of CPU function to compute a bonded interaction.
97 using BondedFunction
= real (*)(int nbonds
,
98 const t_iatom iatoms
[],
99 const t_iparams iparams
[],
111 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
112 BondedFunction
bondedFunction(int ftype
);
114 /*! \brief Calculates all listed force interactions.
116 * Note that pbc_full is used only for position restraints, and is
117 * not initialized if there are none. */
118 void calc_listed(const t_commrec
* cr
,
119 const gmx_multisim_t
* ms
,
120 struct gmx_wallcycle
* wcycle
,
124 gmx::ForceOutputs
* forceOutputs
,
125 const t_forcerec
* fr
,
126 const struct t_pbc
* pbc
,
127 const struct t_pbc
* pbc_full
,
128 const struct t_graph
* g
,
129 gmx_enerdata_t
* enerd
,
133 struct t_fcdata
* fcd
,
135 const gmx::StepWorkload
& stepWork
);
137 /*! \brief As calc_listed(), but only determines the potential energy
138 * for the perturbed interactions.
140 * The shift forces in fr are not affected. */
141 void calc_listed_lambda(const t_idef
* idef
,
143 const t_forcerec
* fr
,
144 const struct t_pbc
* pbc
,
145 const struct t_graph
* g
,
146 gmx_grppairener_t
* grpp
,
148 gmx::ArrayRef
<real
> dvdl
,
152 struct t_fcdata
* fcd
,
153 int* global_atom_index
);
155 /*! \brief Do all aspects of energy and force calculations for mdrun
156 * on the set of listed interactions */
157 void do_force_listed(struct gmx_wallcycle
* wcycle
,
159 const t_lambda
* fepvals
,
161 const gmx_multisim_t
* ms
,
165 gmx::ForceOutputs
* forceOutputs
,
166 const t_forcerec
* fr
,
167 const struct t_pbc
* pbc
,
168 const struct t_graph
* graph
,
169 gmx_enerdata_t
* enerd
,
173 struct t_fcdata
* fcd
,
174 int* global_atom_index
,
175 const gmx::StepWorkload
& stepWork
);
177 /*! \brief Returns true if there are position, distance or orientation restraints. */
178 bool haveRestraints(const t_idef
& idef
, const t_fcdata
& fcd
);
180 /*! \brief Returns true if there are CPU (i.e. not GPU-offloaded) bonded interactions to compute. */
181 bool haveCpuBondeds(const t_forcerec
& fr
);
183 /*! \brief Returns true if there are listed interactions to compute.
185 * NOTE: the current implementation returns true if there are position restraints
186 * or any bonded interactions computed on the CPU.
188 bool haveCpuListedForces(const t_forcerec
& fr
, const t_idef
& idef
, const t_fcdata
& fcd
);