Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / listed_forces / listed_forces.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2017,2018,2019, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \defgroup module_listed_forces Interactions between lists of particles
36  * \ingroup group_mdrun
37  *
38  * \brief Handles computing energies and forces for listed
39  * interactions.
40  *
41  * Located here is the code for
42  * - computing energies and forces for interactions between a small
43      number of particles, e.g bonds, position restraints and listed
44      non-bonded interactions (e.g. 1-4).
45  * - high-level functions used by mdrun for computing a set of such
46      quantities
47  * - managing thread-wise decomposition, thread-local buffer output,
48      and reduction of output data across threads.
49  *
50  * \author Mark Abraham <mark.j.abraham@gmail.com>
51  *
52  */
53 /*! \libinternal \file
54  *
55  * \brief This file contains declarations of high-level functions used
56  * by mdrun to compute energies and forces for listed interactions.
57  *
58  * Clients of libgromacs that want to evaluate listed interactions
59  * should call functions declared here.
60  *
61  * \author Mark Abraham <mark.j.abraham@gmail.com>
62  *
63  * \inlibraryapi
64  * \ingroup module_listed_forces
65  */
66 #ifndef GMX_LISTED_FORCES_LISTED_FORCES_H
67 #define GMX_LISTED_FORCES_LISTED_FORCES_H
68
69 #include "gromacs/math/vectypes.h"
70 #include "gromacs/topology/ifunc.h"
71 #include "gromacs/utility/basedefinitions.h"
72
73 struct gmx_enerdata_t;
74 struct gmx_grppairener_t;
75 struct gmx_multisim_t;
76 class history_t;
77 struct t_commrec;
78 struct t_fcdata;
79 struct t_forcerec;
80 struct t_idef;
81 struct t_graph;
82 struct t_lambda;
83 struct t_mdatoms;
84 struct t_nrnb;
85 class t_state;
86
87 namespace gmx
88 {
89 class ForceOutputs;
90 class StepWorkload;
91 } // namespace gmx
92
93 //! Type of CPU function to compute a bonded interaction.
94 using BondedFunction = real (*)(int              nbonds,
95                                 const t_iatom    iatoms[],
96                                 const t_iparams  iparams[],
97                                 const rvec       x[],
98                                 rvec4            f[],
99                                 rvec             fshift[],
100                                 const t_pbc*     pbc,
101                                 const t_graph*   g,
102                                 real             lambda,
103                                 real*            dvdlambda,
104                                 const t_mdatoms* md,
105                                 t_fcdata*        fcd,
106                                 int*             ddgatindex);
107
108 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
109 BondedFunction bondedFunction(int ftype);
110
111 /*! \brief Calculates all listed force interactions.
112  *
113  * Note that pbc_full is used only for position restraints, and is
114  * not initialized if there are none. */
115 void calc_listed(const t_commrec*         cr,
116                  const gmx_multisim_t*    ms,
117                  struct gmx_wallcycle*    wcycle,
118                  const t_idef*            idef,
119                  const rvec               x[],
120                  history_t*               hist,
121                  gmx::ForceOutputs*       forceOutputs,
122                  const t_forcerec*        fr,
123                  const struct t_pbc*      pbc,
124                  const struct t_pbc*      pbc_full,
125                  const struct t_graph*    g,
126                  gmx_enerdata_t*          enerd,
127                  t_nrnb*                  nrnb,
128                  const real*              lambda,
129                  const t_mdatoms*         md,
130                  struct t_fcdata*         fcd,
131                  int*                     ddgatindex,
132                  const gmx::StepWorkload& stepWork);
133
134 /*! \brief As calc_listed(), but only determines the potential energy
135  * for the perturbed interactions.
136  *
137  * The shift forces in fr are not affected. */
138 void calc_listed_lambda(const t_idef*         idef,
139                         const rvec            x[],
140                         const t_forcerec*     fr,
141                         const struct t_pbc*   pbc,
142                         const struct t_graph* g,
143                         gmx_grppairener_t*    grpp,
144                         real*                 epot,
145                         t_nrnb*               nrnb,
146                         const real*           lambda,
147                         const t_mdatoms*      md,
148                         struct t_fcdata*      fcd,
149                         int*                  global_atom_index);
150
151 /*! \brief Do all aspects of energy and force calculations for mdrun
152  * on the set of listed interactions */
153 void do_force_listed(struct gmx_wallcycle*    wcycle,
154                      const matrix             box,
155                      const t_lambda*          fepvals,
156                      const t_commrec*         cr,
157                      const gmx_multisim_t*    ms,
158                      const t_idef*            idef,
159                      const rvec               x[],
160                      history_t*               hist,
161                      gmx::ForceOutputs*       forceOutputs,
162                      const t_forcerec*        fr,
163                      const struct t_pbc*      pbc,
164                      const struct t_graph*    graph,
165                      gmx_enerdata_t*          enerd,
166                      t_nrnb*                  nrnb,
167                      const real*              lambda,
168                      const t_mdatoms*         md,
169                      struct t_fcdata*         fcd,
170                      int*                     global_atom_index,
171                      const gmx::StepWorkload& stepWork);
172
173 /*! \brief Returns true if there are position, distance or orientation restraints. */
174 bool haveRestraints(const t_idef& idef, const t_fcdata& fcd);
175
176 /*! \brief Returns true if there are CPU (i.e. not GPU-offloaded) bonded interactions to compute. */
177 bool haveCpuBondeds(const t_forcerec& fr);
178
179 /*! \brief Returns true if there are listed interactions to compute.
180  *
181  * NOTE: the current implementation returns true if there are position restraints
182  * or any bonded interactions computed on the CPU.
183  */
184 bool haveCpuListedForces(const t_forcerec& fr, const t_idef& idef, const t_fcdata& fcd);
185
186 #endif