b8abfdd0c6fc1d53faad19424e4d28e488ba8509
[alexxy/gromacs.git] / src / gromacs / listed-forces / bonded.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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, 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.
10  *
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.
15  *
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.
20  *
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.
25  *
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.
33  *
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.
36  */
37 /*! \libinternal \file
38  *
39  * \brief This file contains declarations necessary for low-level
40  * functions for computing energies and forces for bonded interactions.
41  *
42  * \author Mark Abraham <mark.j.abraham@gmail.com>
43  * \inlibraryapi
44  * \ingroup module_listed-forces
45  */
46
47 #ifndef GMX_LISTED_FORCES_BONDED_H
48 #define GMX_LISTED_FORCES_BONDED_H
49
50 #include <stdio.h>
51
52 #include "gromacs/legacyheaders/types/fcdata.h"
53 #include "gromacs/legacyheaders/types/forcerec.h"
54 #include "gromacs/legacyheaders/types/ifunc.h"
55 #include "gromacs/legacyheaders/types/interaction_const.h"
56 #include "gromacs/legacyheaders/types/mdatom.h"
57 #include "gromacs/legacyheaders/types/nrnb.h"
58 #include "gromacs/math/vec.h"
59 #include "gromacs/topology/idef.h"
60 #include "gromacs/utility/basedefinitions.h"
61
62 #ifdef __cplusplus
63 extern "C" {
64 #endif
65
66 struct gmx_wallcycle;
67 struct t_graph;
68 struct t_pbc;
69
70 /*! \brief Calculate bond-angle. No PBC is taken into account (use mol-shift) */
71 real bond_angle(const rvec xi, const rvec xj, const rvec xk,
72                 const struct t_pbc *pbc,
73                 rvec r_ij, rvec r_kj, real *costh,
74                 int *t1, int *t2);  /* out */
75
76 /*! \brief Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
77 real dih_angle(const rvec xi, const rvec xj, const rvec xk, const rvec xl,
78                const struct t_pbc *pbc,
79                rvec r_ij, rvec r_kj, rvec r_kl, rvec m, rvec n, /* out */
80                real *sign,
81                int *t1, int *t2, int *t3);
82
83 /*! \brief Do an update of the forces for dihedral potentials */
84 void do_dih_fup(int i, int j, int k, int l, real ddphi,
85                 rvec r_ij, rvec r_kj, rvec r_kl,
86                 rvec m, rvec n, rvec f[], rvec fshift[],
87                 const struct t_pbc *pbc, const struct t_graph *g,
88                 const rvec *x, int t1, int t2, int t3);
89
90 /*! \brief Make a dihedral fall in the range (-pi,pi) */
91 void make_dp_periodic(real *dp);
92
93 /*! \brief Compute CMAP dihedral energies and forces */
94 real
95     cmap_dihs(int nbonds,
96               const t_iatom forceatoms[], const t_iparams forceparams[],
97               const gmx_cmap_t *cmap_grid,
98               const rvec x[], rvec f[], rvec fshift[],
99               const struct t_pbc *pbc, const struct t_graph *g,
100               real gmx_unused lambda, real gmx_unused *dvdlambda,
101               const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
102               int  gmx_unused *global_atom_index);
103
104 //! \cond
105 /*************************************************************************
106  *
107  *  Bonded force functions
108  *
109  *************************************************************************/
110 t_ifunc bonds, g96bonds, morse_bonds, cubic_bonds, FENE_bonds, restraint_bonds;
111 t_ifunc angles, g96angles, cross_bond_bond, cross_bond_angle, urey_bradley, quartic_angles, linear_angles;
112 t_ifunc restrangles;
113 t_ifunc pdihs, idihs, rbdihs;
114 t_ifunc restrdihs, cbtdihs;
115 t_ifunc tab_bonds, tab_angles, tab_dihs;
116 t_ifunc polarize, anharm_polarize, water_pol, thole_pol, angres, angresz, dihres, unimplemented;
117
118 /* As pdihs above, but without calculating energies and shift forces */
119 void
120     pdihs_noener(int nbonds,
121                  const t_iatom forceatoms[], const t_iparams forceparams[],
122                  const rvec x[], rvec f[],
123                  const struct t_pbc gmx_unused *pbc,
124                  const struct t_graph gmx_unused *g,
125                  real lambda,
126                  const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
127                  int gmx_unused *global_atom_index);
128 //! \endcond
129
130 #ifdef __cplusplus
131 }
132 #endif
133
134 #endif