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.
54 int glatnr(int *global_atom_index
, int i
);
55 /* Returns the global topology atom number belonging to local atom index i.
56 * This function is intended for writing ascii output
57 * and returns atom numbers starting at 1.
58 * When global_atom_index=NULL returns i+1.
61 void calc_bonds(const gmx_multisim_t
*ms
,
63 rvec x
[], history_t
*hist
,
64 rvec f
[], t_forcerec
*fr
,
65 const struct t_pbc
*pbc
, const struct t_graph
*g
,
66 gmx_enerdata_t
*enerd
, t_nrnb
*nrnb
, real
*lambda
,
68 t_fcdata
*fcd
, int *ddgatindex
,
69 t_atomtypes
*atype
, gmx_genborn_t
*born
,
72 * The function calc_bonds() calculates all bonded force interactions.
73 * The "bonds" are specified as follows:
75 * the total number of bonded interactions.
77 * specifies which atoms are involved in a bond of a certain
78 * type, see also struct t_idef.
79 * t_functype *functype
80 * defines for every bonded force type what type of function to
81 * use, see also struct t_idef.
82 * t_iparams *forceparams
83 * defines the parameters for every bond type, see also struct
86 * total potential energy split up over the function types.
88 * global atom number indices, should be NULL when not using DD.
90 * the total potential energy (sum over epot).
93 void calc_bonds_lambda(const t_idef
*idef
,
96 const struct t_pbc
*pbc
, const struct t_graph
*g
,
97 gmx_grppairener_t
*grpp
, real
*epot
, t_nrnb
*nrnb
,
100 t_fcdata
*fcd
, int *global_atom_index
);
101 /* As calc_bonds, but only determines the potential energy
102 * for the perturbed interactions.
103 * The shift forces in fr are not affected.
106 real
posres(int nbonds
,
107 const t_iatom forceatoms
[], const t_iparams forceparams
[],
108 const rvec x
[], rvec f
[], rvec vir_diag
,
110 real lambda
, real
*dvdlambda
,
111 int refcoord_scaling
, int ePBC
, rvec comA
, rvec comB
);
112 /* Position restraints require a different pbc treatment from other bondeds */
114 real
fbposres(int nbonds
,
115 const t_iatom forceatoms
[], const t_iparams forceparams
[],
116 const rvec x
[], rvec f
[], rvec vir_diag
,
117 struct t_pbc
*pbc
, int refcoord_scaling
, int ePBC
, rvec com
);
118 /* Flat-bottom posres. Same PBC treatment as in normal position restraints */
120 real
bond_angle(const rvec xi
, const rvec xj
, const rvec xk
,
121 const struct t_pbc
*pbc
,
122 rvec r_ij
, rvec r_kj
, real
*costh
,
123 int *t1
, int *t2
); /* out */
124 /* Calculate bond-angle. No PBC is taken into account (use mol-shift) */
126 real
dih_angle(const rvec xi
, const rvec xj
, const rvec xk
, const rvec xl
,
127 const struct t_pbc
*pbc
,
128 rvec r_ij
, rvec r_kj
, rvec r_kl
, rvec m
, rvec n
, /* out */
130 int *t1
, int *t2
, int *t3
);
131 /* Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
133 void do_dih_fup(int i
, int j
, int k
, int l
, real ddphi
,
134 rvec r_ij
, rvec r_kj
, rvec r_kl
,
135 rvec m
, rvec n
, rvec f
[], rvec fshift
[],
136 const struct t_pbc
*pbc
, const struct t_graph
*g
,
137 const rvec
*x
, int t1
, int t2
, int t3
);
138 /* Do an update of the forces for dihedral potentials */
140 void make_dp_periodic(real
*dp
);
141 /* make a dihedral fall in the range (-pi,pi) */
143 /*************************************************************************
145 * Bonded force functions
147 *************************************************************************/
148 t_ifunc bonds
, g96bonds
, morse_bonds
, cubic_bonds
, FENE_bonds
, restraint_bonds
;
149 t_ifunc angles
, g96angles
, cross_bond_bond
, cross_bond_angle
, urey_bradley
, quartic_angles
, linear_angles
;
151 t_ifunc pdihs
, idihs
, rbdihs
;
152 t_ifunc restrdihs
, cbtdihs
;
153 t_ifunc tab_bonds
, tab_angles
, tab_dihs
;
154 t_ifunc polarize
, anharm_polarize
, water_pol
, thole_pol
, angres
, angresz
, dihres
, unimplemented
;
157 /* Divided the bonded interactions over the threads, count=fr->nthreads
158 * and set up the bonded thread-force buffer reduction.
159 * This should be called each time the bonded setup changes;
160 * i.e. at start-up without domain decomposition and at DD.
162 void setup_bonded_threading(t_forcerec
*fr
, t_idef
*idef
);
168 #endif /* _bondf_h */