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,2015,2016,2017, 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 /*! \libinternal \file
39 * \brief This file contains declarations necessary for low-level
40 * functions for computing energies and forces for bonded interactions.
42 * \author Mark Abraham <mark.j.abraham@gmail.com>
44 * \ingroup module_listed-forces
47 #ifndef GMX_LISTED_FORCES_BONDED_H
48 #define GMX_LISTED_FORCES_BONDED_H
52 #include "gromacs/gmxlib/nrnb.h"
53 #include "gromacs/math/vec.h"
54 #include "gromacs/mdtypes/fcdata.h"
55 #include "gromacs/mdtypes/forcerec.h"
56 #include "gromacs/mdtypes/interaction_const.h"
57 #include "gromacs/mdtypes/mdatom.h"
58 #include "gromacs/topology/idef.h"
59 #include "gromacs/topology/ifunc.h"
60 #include "gromacs/utility/basedefinitions.h"
66 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
67 real
bond_angle(const rvec xi
, const rvec xj
, const rvec xk
,
68 const struct t_pbc
*pbc
,
69 rvec r_ij
, rvec r_kj
, real
*costh
,
70 int *t1
, int *t2
); /* out */
72 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
73 real
dih_angle(const rvec xi
, const rvec xj
, const rvec xk
, const rvec xl
,
74 const struct t_pbc
*pbc
,
75 rvec r_ij
, rvec r_kj
, rvec r_kl
, rvec m
, rvec n
, /* out */
77 int *t1
, int *t2
, int *t3
);
79 /*! \brief Do an update of the forces for dihedral potentials */
80 void do_dih_fup(int i
, int j
, int k
, int l
, real ddphi
,
81 rvec r_ij
, rvec r_kj
, rvec r_kl
,
82 rvec m
, rvec n
, rvec4 f
[], rvec fshift
[],
83 const struct t_pbc
*pbc
, const struct t_graph
*g
,
84 const rvec
*x
, int t1
, int t2
, int t3
);
86 /*! \brief Make a dihedral fall in the range (-pi,pi) */
87 void make_dp_periodic(real
*dp
);
89 /*! \brief Compute CMAP dihedral energies and forces */
92 const t_iatom forceatoms
[], const t_iparams forceparams
[],
93 const gmx_cmap_t
*cmap_grid
,
94 const rvec x
[], rvec4 f
[], rvec fshift
[],
95 const struct t_pbc
*pbc
, const struct t_graph
*g
,
96 real gmx_unused lambda
, real gmx_unused
*dvdlambda
,
97 const t_mdatoms gmx_unused
*md
, t_fcdata gmx_unused
*fcd
,
98 int gmx_unused
*global_atom_index
);
101 /*************************************************************************
103 * Bonded force functions
105 *************************************************************************/
106 t_ifunc bonds
, g96bonds
, morse_bonds
, cubic_bonds
, FENE_bonds
, restraint_bonds
;
107 t_ifunc angles
, g96angles
, cross_bond_bond
, cross_bond_angle
, urey_bradley
, quartic_angles
, linear_angles
;
109 t_ifunc pdihs
, idihs
, rbdihs
;
110 t_ifunc restrdihs
, cbtdihs
;
111 t_ifunc tab_bonds
, tab_angles
, tab_dihs
;
112 t_ifunc polarize
, anharm_polarize
, water_pol
, thole_pol
, angres
, angresz
, dihres
, unimplemented
;
114 /* As pdihs(), but without calculating energies and shift forces */
116 pdihs_noener(int nbonds
,
117 const t_iatom forceatoms
[], const t_iparams forceparams
[],
118 const rvec x
[], rvec4 f
[],
119 const struct t_pbc gmx_unused
*pbc
,
120 const struct t_graph gmx_unused
*g
,
122 const t_mdatoms gmx_unused
*md
, t_fcdata gmx_unused
*fcd
,
123 int gmx_unused
*global_atom_index
);
125 /* TODO these declarations should be internal to the module */
127 /* As angles(), but using SIMD to calculate many angles at once.
128 * This routines does not calculate energies and shift forces.
131 angles_noener_simd(int nbonds
,
132 const t_iatom forceatoms
[], const t_iparams forceparams
[],
133 const rvec x
[], rvec4 f
[],
134 const struct t_pbc
*pbc
,
135 const struct t_graph gmx_unused
*g
,
136 real gmx_unused lambda
,
137 const t_mdatoms gmx_unused
*md
, t_fcdata gmx_unused
*fcd
,
138 int gmx_unused
*global_atom_index
);
140 /* As pdihs_noener(), but using SIMD to calculate many dihedrals at once. */
142 pdihs_noener_simd(int nbonds
,
143 const t_iatom forceatoms
[], const t_iparams forceparams
[],
144 const rvec x
[], rvec4 f
[],
145 const struct t_pbc
*pbc
,
146 const struct t_graph gmx_unused
*g
,
147 real gmx_unused lambda
,
148 const t_mdatoms gmx_unused
*md
, t_fcdata gmx_unused
*fcd
,
149 int gmx_unused
*global_atom_index
);
151 /* As rbdihs(), when not needing energy or shift force, using SIMD to calculate many dihedrals at once. */
153 rbdihs_noener_simd(int nbonds
,
154 const t_iatom forceatoms
[], const t_iparams forceparams
[],
155 const rvec x
[], rvec4 f
[],
156 const struct t_pbc
*pbc
,
157 const struct t_graph gmx_unused
*g
,
158 real gmx_unused lambda
,
159 const t_mdatoms gmx_unused
*md
, t_fcdata gmx_unused
*fcd
,
160 int gmx_unused
*global_atom_index
);