SIMD accelerated Urey-Bradley
[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,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.
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 <cstdio>
51
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"
61
62
63 struct t_graph;
64 struct t_pbc;
65
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 */
71
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 */
76                int *t1, int *t2, int *t3);
77
78 /*! \brief Do an update of the forces for dihedral potentials */
79 void do_dih_fup(int i, int j, int k, int l, real ddphi,
80                 rvec r_ij, rvec r_kj, rvec r_kl,
81                 rvec m, rvec n, rvec4 f[], rvec fshift[],
82                 const struct t_pbc *pbc, const struct t_graph *g,
83                 const rvec *x, int t1, int t2, int t3);
84
85 /*! \brief Make a dihedral fall in the range (-pi,pi) */
86 void make_dp_periodic(real *dp);
87
88 /*! \brief Compute CMAP dihedral energies and forces */
89 real
90     cmap_dihs(int nbonds,
91               const t_iatom forceatoms[], const t_iparams forceparams[],
92               const gmx_cmap_t *cmap_grid,
93               const rvec x[], rvec4 f[], rvec fshift[],
94               const struct t_pbc *pbc, const struct t_graph *g,
95               real gmx_unused lambda, real gmx_unused *dvdlambda,
96               const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
97               int  gmx_unused *global_atom_index);
98
99 //! \cond
100 /*************************************************************************
101  *
102  *  Bonded force functions
103  *
104  *************************************************************************/
105 t_ifunc bonds, g96bonds, morse_bonds, cubic_bonds, FENE_bonds, restraint_bonds;
106 t_ifunc angles, g96angles, cross_bond_bond, cross_bond_angle, urey_bradley, quartic_angles, linear_angles;
107 t_ifunc restrangles;
108 t_ifunc pdihs, idihs, rbdihs;
109 t_ifunc restrdihs, cbtdihs;
110 t_ifunc tab_bonds, tab_angles, tab_dihs;
111 t_ifunc polarize, anharm_polarize, water_pol, thole_pol, angres, angresz, dihres, unimplemented;
112
113 /* As pdihs(), but without calculating energies and shift forces */
114 void
115     pdihs_noener(int nbonds,
116                  const t_iatom forceatoms[], const t_iparams forceparams[],
117                  const rvec x[], rvec4 f[],
118                  const struct t_pbc gmx_unused *pbc,
119                  const struct t_graph gmx_unused *g,
120                  real lambda,
121                  const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
122                  int gmx_unused *global_atom_index);
123
124 /* TODO these declarations should be internal to the module */
125
126 /* As angles(), but using SIMD to calculate many angles at once.
127  * This routines does not calculate energies and shift forces.
128  */
129 void
130     angles_noener_simd(int nbonds,
131                        const t_iatom forceatoms[], const t_iparams forceparams[],
132                        const rvec x[], rvec4 f[],
133                        const struct t_pbc *pbc,
134                        const struct t_graph gmx_unused *g,
135                        real gmx_unused lambda,
136                        const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
137                        int gmx_unused *global_atom_index);
138
139 /* As urey_bradley, but using SIMD to calculate many potentials at once.
140  * This routines does not calculate energies and shift forces.
141  */
142 void urey_bradley_noener_simd(int nbonds,
143                               const t_iatom forceatoms[], const t_iparams forceparams[],
144                               const rvec x[], rvec4 f[],
145                               const t_pbc *pbc, const t_graph gmx_unused *g,
146                               real gmx_unused lambda,
147                               const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
148                               int gmx_unused *global_atom_index);
149
150 /* As pdihs_noener(), but using SIMD to calculate many dihedrals at once. */
151 void
152     pdihs_noener_simd(int nbonds,
153                       const t_iatom forceatoms[], const t_iparams forceparams[],
154                       const rvec x[], rvec4 f[],
155                       const struct t_pbc *pbc,
156                       const struct t_graph gmx_unused *g,
157                       real gmx_unused lambda,
158                       const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
159                       int gmx_unused *global_atom_index);
160
161 /* As rbdihs(), when not needing energy or shift force, using SIMD to calculate many dihedrals at once. */
162 void
163     rbdihs_noener_simd(int nbonds,
164                        const t_iatom forceatoms[], const t_iparams forceparams[],
165                        const rvec x[], rvec4 f[],
166                        const struct t_pbc *pbc,
167                        const struct t_graph gmx_unused *g,
168                        real gmx_unused lambda,
169                        const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
170                        int gmx_unused *global_atom_index);
171
172 //! \endcond
173
174 #endif