9c1d4526162d69ac75d24e04efd8371124f7006a
[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, 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 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 ForceWithVirial;
90 }
91
92 //! Type of CPU function to compute a bonded interaction.
93 using BondedFunction = real(*)(int nbonds, const t_iatom iatoms[],
94                                const t_iparams iparams[],
95                                const rvec x[], rvec4 f[], rvec fshift[],
96                                const t_pbc *pbc, const t_graph *g,
97                                real lambda, real *dvdlambda,
98                                const t_mdatoms *md, t_fcdata *fcd,
99                                int *ddgatindex);
100
101 //! Getter for finding a callable CPU function to compute an \c ftype interaction.
102 BondedFunction bondedFunction(int ftype);
103
104 /*! \brief Calculates all listed force interactions.
105  *
106  * Note that pbc_full is used only for position restraints, and is
107  * not initialized if there are none. */
108 void calc_listed(const t_commrec *cr,
109                  const gmx_multisim_t *ms,
110                  struct gmx_wallcycle *wcycle,
111                  const t_idef *idef,
112                  const rvec x[], history_t *hist,
113                  rvec f[],
114                  gmx::ForceWithVirial *forceWithVirial,
115                  const t_forcerec *fr,
116                  const struct t_pbc *pbc, const struct t_pbc *pbc_full,
117                  const struct t_graph *g,
118                  gmx_enerdata_t *enerd, t_nrnb *nrnb, const real *lambda,
119                  const t_mdatoms *md,
120                  struct t_fcdata *fcd, int *ddgatindex,
121                  int force_flags);
122
123 /*! \brief As calc_listed(), but only determines the potential energy
124  * for the perturbed interactions.
125  *
126  * The shift forces in fr are not affected. */
127 void calc_listed_lambda(const t_idef *idef,
128                         const rvec x[],
129                         const t_forcerec *fr,
130                         const struct t_pbc *pbc, const struct t_graph *g,
131                         gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
132                         const real *lambda,
133                         const t_mdatoms *md,
134                         struct t_fcdata *fcd, int *global_atom_index);
135
136 /*! \brief Do all aspects of energy and force calculations for mdrun
137  * on the set of listed interactions */
138 void
139 do_force_listed(struct gmx_wallcycle           *wcycle,
140                 matrix                          box,
141                 const t_lambda                 *fepvals,
142                 const t_commrec                *cr,
143                 const gmx_multisim_t           *ms,
144                 const t_idef                   *idef,
145                 const rvec                      x[],
146                 history_t                      *hist,
147                 rvec                           *forceForUseWithShiftForces,
148                 gmx::ForceWithVirial           *forceWithVirial,
149                 const t_forcerec               *fr,
150                 const struct t_pbc             *pbc,
151                 const struct t_graph           *graph,
152                 gmx_enerdata_t                 *enerd,
153                 t_nrnb                         *nrnb,
154                 const real                     *lambda,
155                 const t_mdatoms                *md,
156                 struct t_fcdata                *fcd,
157                 int                            *global_atom_index,
158                 int                             flags);
159
160 #endif