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 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
44 #include "visibility.h"
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.
62 void calc_bonds(FILE *fplog
, const gmx_multisim_t
*ms
,
64 rvec x
[], history_t
*hist
,
65 rvec f
[], t_forcerec
*fr
,
66 const t_pbc
*pbc
, const t_graph
*g
,
67 gmx_enerdata_t
*enerd
, t_nrnb
*nrnb
, real
*lambda
,
69 t_fcdata
*fcd
, int *ddgatindex
,
70 t_atomtypes
*atype
, gmx_genborn_t
*born
,
72 gmx_bool bPrintSepPot
, gmx_large_int_t step
);
74 * The function calc_bonds() calculates all bonded force interactions.
75 * The "bonds" are specified as follows:
77 * the total number of bonded interactions.
79 * specifies which atoms are involved in a bond of a certain
80 * type, see also struct t_idef.
81 * t_functype *functype
82 * defines for every bonded force type what type of function to
83 * use, see also struct t_idef.
84 * t_iparams *forceparams
85 * defines the parameters for every bond type, see also struct
88 * total potential energy split up over the function types.
90 * global atom number indices, should be NULL when not using DD.
91 * gmx_bool bPrintSepPot
92 * if TRUE print local potential and dVdlambda for each bonded type.
94 * used with bPrintSepPot
96 * the total potential energy (sum over epot).
100 void calc_bonds_lambda(FILE *fplog
,
104 const t_pbc
*pbc
, const t_graph
*g
,
105 gmx_grppairener_t
*grpp
, real
*epot
, t_nrnb
*nrnb
,
108 t_fcdata
*fcd
, int *global_atom_index
);
109 /* As calc_bonds, but only determines the potential energy
110 * for the perturbed interactions.
111 * The shift forces in fr are not affected.
115 real
posres(int nbonds
,
116 const t_iatom forceatoms
[], const t_iparams forceparams
[],
117 const rvec x
[], rvec f
[], rvec vir_diag
,
119 real lambda
, real
*dvdlambda
,
120 int refcoord_scaling
, int ePBC
, rvec comA
, rvec comB
);
121 /* Position restraints require a different pbc treatment from other bondeds */
124 real
bond_angle(const rvec xi
, const rvec xj
, const rvec xk
,
126 rvec r_ij
, rvec r_kj
, real
*costh
,
127 int *t1
, int *t2
); /* out */
128 /* Calculate bond-angle. No PBC is taken into account (use mol-shift) */
131 real
dih_angle(const rvec xi
, const rvec xj
, const rvec xk
, const rvec xl
,
133 rvec r_ij
, rvec r_kj
, rvec r_kl
, rvec m
, rvec n
, /* out */
135 int *t1
, int *t2
, int *t3
);
136 /* Calculate dihedral-angle. No PBC is taken into account (use mol-shift) */
138 void do_dih_fup(int i
, int j
, int k
, int l
, real ddphi
,
139 rvec r_ij
, rvec r_kj
, rvec r_kl
,
140 rvec m
, rvec n
, rvec f
[], rvec fshift
[],
141 const t_pbc
*pbc
, const t_graph
*g
,
142 const rvec
*x
, int t1
, int t2
, int t3
);
143 /* Do an update of the forces for dihedral potentials */
145 void make_dp_periodic(real
*dp
);
146 /* make a dihedral fall in the range (-pi,pi) */
148 /*************************************************************************
150 * Bonded force functions
152 *************************************************************************/
153 t_ifunc bonds
, g96bonds
, morse_bonds
, cubic_bonds
, FENE_bonds
, restraint_bonds
;
154 t_ifunc angles
, g96angles
, cross_bond_bond
, cross_bond_angle
, urey_bradley
, quartic_angles
, linear_angles
;
155 t_ifunc pdihs
, idihs
, rbdihs
;
156 t_ifunc tab_bonds
, tab_angles
, tab_dihs
;
157 t_ifunc polarize
, anharm_polarize
, water_pol
, thole_pol
, angres
, angresz
, dihres
, unimplemented
;
160 /* Divided the bonded interactions over the threads, count=fr->nthreads
161 * and set up the bonded thread-force buffer reduction.
162 * This should be called each time the bonded setup changes;
163 * i.e. at start-up without domain decomposition and at DD.
166 void setup_bonded_threading(t_forcerec
*fr
, t_idef
*idef
);
172 #endif /* _bondf_h */