Remove use of graph in do_force()
[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 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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36 /*! \defgroup module_listed_forces Interactions between lists of particles
37  * \ingroup group_mdrun
38  *
39  * \brief Handles computing energies and forces for listed
40  * interactions.
41  *
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
47      quantities
48  * - managing thread-wise decomposition, thread-local buffer output,
49      and reduction of output data across threads.
50  *
51  * \author Mark Abraham <mark.j.abraham@gmail.com>
52  *
53  */
54 /*! \libinternal \file
55  *
56  * \brief This file contains declarations of high-level functions used
57  * by mdrun to compute energies and forces for listed interactions.
58  *
59  * Clients of libgromacs that want to evaluate listed interactions
60  * should call functions declared here.
61  *
62  * \author Mark Abraham <mark.j.abraham@gmail.com>
63  *
64  * \inlibraryapi
65  * \ingroup module_listed_forces
66  */
67 #ifndef GMX_LISTED_FORCES_LISTED_FORCES_H
68 #define GMX_LISTED_FORCES_LISTED_FORCES_H
69
70 #include "gromacs/math/vectypes.h"
71 #include "gromacs/topology/ifunc.h"
72 #include "gromacs/utility/basedefinitions.h"
73
74 struct gmx_enerdata_t;
75 struct gmx_grppairener_t;
76 struct gmx_multisim_t;
77 class history_t;
78 class InteractionDefinitions;
79 struct t_commrec;
80 struct t_fcdata;
81 struct t_forcerec;
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 template<typename>
92 class ArrayRef;
93 } // namespace gmx
94
95 //! Type of CPU function to compute a bonded interaction.
96 using BondedFunction = real (*)(int              nbonds,
97                                 const t_iatom    iatoms[],
98                                 const t_iparams  iparams[],
99                                 const rvec       x[],
100                                 rvec4            f[],
101                                 rvec             fshift[],
102                                 const t_pbc*     pbc,
103                                 real             lambda,
104                                 real*            dvdlambda,
105                                 const t_mdatoms* md,
106                                 t_fcdata*        fcd,
107                                 int*             ddgatindex);
108
109 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
110 BondedFunction bondedFunction(int ftype);
111
112 /*! \brief Do all aspects of energy and force calculations for mdrun
113  * on the set of listed interactions
114  *
115  * xWholeMolecules only needs to contain whole molecules when orientation
116  * restraints need to be computed and can be empty otherwise.
117  */
118 void do_force_listed(struct gmx_wallcycle*          wcycle,
119                      const matrix                   box,
120                      const t_lambda*                fepvals,
121                      const t_commrec*               cr,
122                      const gmx_multisim_t*          ms,
123                      const InteractionDefinitions&  idef,
124                      const rvec                     x[],
125                      gmx::ArrayRef<const gmx::RVec> xWholeMolecules,
126                      history_t*                     hist,
127                      gmx::ForceOutputs*             forceOutputs,
128                      const t_forcerec*              fr,
129                      const struct t_pbc*            pbc,
130                      gmx_enerdata_t*                enerd,
131                      t_nrnb*                        nrnb,
132                      const real*                    lambda,
133                      const t_mdatoms*               md,
134                      struct t_fcdata*               fcd,
135                      int*                           global_atom_index,
136                      const gmx::StepWorkload&       stepWork);
137
138 /*! \brief Returns true if there are position, distance or orientation restraints. */
139 bool haveRestraints(const InteractionDefinitions& idef, const t_fcdata& fcd);
140
141 /*! \brief Returns true if there are CPU (i.e. not GPU-offloaded) bonded interactions to compute. */
142 bool haveCpuBondeds(const t_forcerec& fr);
143
144 /*! \brief Returns true if there are listed interactions to compute.
145  *
146  * NOTE: the current implementation returns true if there are position restraints
147  * or any bonded interactions computed on the CPU.
148  */
149 bool haveCpuListedForces(const t_forcerec& fr, const InteractionDefinitions& idef, const t_fcdata& fcd);
150
151 #endif