2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team.
6 * Copyright (c) 2013,2014, by the GROMACS development team, led by
7 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8 * and including many others, as listed in the AUTHORS file in the
9 * top-level source directory and at http://www.gromacs.org.
11 * GROMACS is free software; you can redistribute it and/or
12 * modify it under the terms of the GNU Lesser General Public License
13 * as published by the Free Software Foundation; either version 2.1
14 * of the License, or (at your option) any later version.
16 * GROMACS is distributed in the hope that it will be useful,
17 * but WITHOUT ANY WARRANTY; without even the implied warranty of
18 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 * Lesser General Public License for more details.
21 * You should have received a copy of the GNU Lesser General Public
22 * License along with GROMACS; if not, see
23 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
26 * If you want to redistribute modifications to GROMACS, please
27 * consider that scientific software is very special. Version
28 * control is crucial - bugs must be traceable. We will be happy to
29 * consider code for inclusion in the official distribution, but
30 * derived work must not be called official GROMACS. Details are found
31 * in the README & COPYING files - if they are missing, get the
32 * official version at http://www.gromacs.org.
34 * To help us fund GROMACS development, we humbly ask that you cite
35 * the research papers on the package. Check out http://www.gromacs.org.
37 /*! \defgroup module_listed-forces Interactions between lists of particles
38 * \ingroup group_mdrun
40 * \brief Computes energies and forces for interactions between a
41 * small number of particles, e.g bonds.
43 * More functionality will move into this module shortly.
45 * \author Mark Abraham <mark.j.abraham@gmail.com>
52 * This file contains function declarations necessary for mdrun and tools
53 * to compute energies and forces for bonded interactions.
55 * \author Mark Abraham <mark.j.abraham@gmail.com>
57 * \ingroup module_listed-forces
60 #ifndef GMX_LISTED_FORCES_BONDED_H
61 #define GMX_LISTED_FORCES_BONDED_H
65 #include "gromacs/legacyheaders/nrnb.h"
66 #include "gromacs/legacyheaders/typedefs.h"
75 /*! \brief Return whether this is an interaction that actually
76 * calculates a potential and works on multiple atoms (not e.g. a
77 * connection or a position restraint).
79 * \todo This function could go away when idef is not a big bucket of
82 ftype_is_bonded_potential(int ftype);
84 /*! \brief Calculates all bonded force interactions. */
85 void calc_bonds(const gmx_multisim_t *ms,
87 const rvec x[], history_t *hist,
88 rvec f[], t_forcerec *fr,
89 const struct t_pbc *pbc, const struct t_graph *g,
90 gmx_enerdata_t *enerd, t_nrnb *nrnb, real *lambda,
92 t_fcdata *fcd, int *ddgatindex,
95 /*! \brief As calc_bonds, but only determines the potential energy
96 * for the perturbed interactions.
97 * The shift forces in fr are not affected. */
98 void calc_bonds_lambda(const t_idef *idef,
101 const struct t_pbc *pbc, const struct t_graph *g,
102 gmx_grppairener_t *grpp, real *epot, t_nrnb *nrnb,
105 t_fcdata *fcd, int *global_atom_index);
107 /*! \brief Position restraints require a different pbc treatment from other bondeds */
108 real posres(int nbonds,
109 const t_iatom forceatoms[], const t_iparams forceparams[],
110 const rvec x[], rvec f[], rvec vir_diag,
112 real lambda, real *dvdlambda,
113 int refcoord_scaling, int ePBC, rvec comA, rvec comB);
115 /*! \brief Flat-bottom posres. Same PBC treatment as in normal position restraints */
116 real fbposres(int nbonds,
117 const t_iatom forceatoms[], const t_iparams forceparams[],
118 const rvec x[], rvec f[], rvec vir_diag,
119 struct t_pbc *pbc, int refcoord_scaling, int ePBC, rvec com);
121 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
122 real bond_angle(const rvec xi, const rvec xj, const rvec xk,
123 const struct t_pbc *pbc,
124 rvec r_ij, rvec r_kj, real *costh,
125 int *t1, int *t2); /* out */
127 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
128 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
129 const struct t_pbc *pbc,
130 rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n, /* out */
132 int *t1, int *t2, int *t3);
134 /*! \brief Do an update of the forces for dihedral potentials */
135 void do_dih_fup(int i, int j, int k, int l, real ddphi,
136 rvec r_ij, rvec r_kj, rvec r_kl,
137 rvec m, rvec n, rvec f[], rvec fshift[],
138 const struct t_pbc *pbc, const struct t_graph *g,
139 const rvec *x, int t1, int t2, int t3);
141 /*! \brief Make a dihedral fall in the range (-pi,pi) */
142 void make_dp_periodic(real *dp);
145 /*************************************************************************
147 * Bonded force functions
149 *************************************************************************/
150 t_ifunc bonds, g96bonds, morse_bonds, cubic_bonds, FENE_bonds, restraint_bonds;
151 t_ifunc angles, g96angles, cross_bond_bond, cross_bond_angle, urey_bradley, quartic_angles, linear_angles;
153 t_ifunc pdihs, idihs, rbdihs;
154 t_ifunc restrdihs, cbtdihs;
155 t_ifunc tab_bonds, tab_angles, tab_dihs;
156 t_ifunc polarize, anharm_polarize, water_pol, thole_pol, angres, angresz, dihres, unimplemented;