3 * This source code is part of
7 * GROningen MAchine for Chemical Simulations
10 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
11 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
12 * Copyright (c) 2001-2004, The GROMACS development team,
13 * check out http://www.gromacs.org for more information.
15 * This program is free software; you can redistribute it and/or
16 * modify it under the terms of the GNU General Public License
17 * as published by the Free Software Foundation; either version 2
18 * of the License, or (at your option) any later version.
20 * If you want to redistribute modifications, please consider that
21 * scientific software is very special. Version control is crucial -
22 * bugs must be traceable. We will be happy to consider code for
23 * inclusion in the official distribution, but derived work must not
24 * be called official GROMACS. Details are found in the README & COPYING
25 * files - if they are missing, get the official version at www.gromacs.org.
27 * To help us fund GROMACS development, we humbly ask that you cite
28 * the papers on the package - you can find them in the top README file.
30 * For more info, check our website at http://www.gromacs.org
33 * GRoups of Organic Molecules in ACtion for Science
47 /* check kernel/toppush.c when you change these numbers */
49 #define MAXFORCEPARAM 12
53 typedef atom_id t_iatom
;
55 /* this MUST correspond to the
56 t_interaction_function[F_NRE] in gmxlib/ifunc.c */
135 F_NRE
/* This number is for the total number of energies */
140 /* Some parameters have A and B values for free energy calculations.
141 * The B values are not used for regular simulations of course.
142 * Free Energy for nonbondeds can be computed by changing the atom type.
143 * The harmonic type is used for all harmonic potentials:
144 * bonds, angles and improper dihedrals
146 struct {real a
,b
,c
; } bham
;
147 struct {real rA
,krA
,rB
,krB
; } harmonic
;
148 struct {real lowA
,up1A
,up2A
,kA
,lowB
,up1B
,up2B
,kB
; } restraint
;
149 /* No free energy supported for cubic bonds, FENE, WPOL or cross terms */
150 struct {real b0
,kb
,kcub
; } cubic
;
151 struct {real bm
,kb
; } fene
;
152 struct {real r1e
,r2e
,krr
; } cross_bb
;
153 struct {real r1e
,r2e
,r3e
,krt
; } cross_ba
;
154 struct {real theta
,ktheta
,r13
,kUB
; } u_b
;
155 struct {real theta
,c
[5]; } qangle
;
156 struct {real alpha
; } polarize
;
157 struct {real al_x
,al_y
,al_z
,rOH
,rHH
,rOD
; } wpol
;
158 struct {real a
,alpha1
,alpha2
,rfac
; } thole
;
159 struct {real c6
,c12
; } lj
;
160 struct {real c6A
,c12A
,c6B
,c12B
; } lj14
;
161 struct {real fqq
,qi
,qj
,c6
,c12
; } ljc14
;
162 struct {real qi
,qj
,c6
,c12
; } ljcnb
;
163 /* Proper dihedrals can not have different multiplicity when
164 * doing free energy calculations, because the potential would not
165 * be periodic anymore.
167 struct {real phiA
,cpA
;int mult
;real phiB
,cpB
; } pdihs
;
168 struct {real dA
,dB
; } constr
;
169 /* Settle can not be used for Free energy calculations of water bond geometry.
170 * Use shake (or lincs) instead if you have to change the water bonds.
172 struct {real doh
,dhh
; } settle
;
173 /* No free energy supported for morse bonds */
174 struct {real b0
,cb
,beta
; } morse
;
175 struct {real pos0A
[DIM
],fcA
[DIM
],pos0B
[DIM
],fcB
[DIM
]; } posres
;
176 struct {real rbcA
[NR_RBDIHS
], rbcB
[NR_RBDIHS
]; } rbdihs
;
177 struct {real a
,b
,c
,d
,e
,f
; } vsite
;
178 struct {int n
; real a
; } vsiten
;
179 /* NOTE: npair is only set after reading the tpx file */
180 struct {real low
,up1
,up2
,kfac
;int type
,label
,npair
; } disres
;
181 struct {real phi
,dphi
,kfac
;int label
,power
; } dihres
;
182 struct {int ex
,power
,label
; real c
,obs
,kfac
; } orires
;
183 struct {int table
;real kA
;real kB
; } tab
;
184 struct {real sar
,st
,pi
,gbr
,bmlt
; } gb
;
185 struct {int cmapA
,cmapB
; } cmap
;
186 struct {real buf
[MAXFORCEPARAM
]; } generic
; /* Conversion */
189 typedef int t_functype
;
192 * The nonperturbed/perturbed interactions are now separated (sorted) in the
193 * ilist, such that the first 0..(nr_nonperturbed-1) ones are exactly that, and
194 * the remaining ones from nr_nonperturbed..(nr-1) are perturbed bonded
206 * The struct t_ilist defines a list of atoms with their interactions.
207 * General field description:
209 * the size (nr elements) of the interactions array (iatoms[]).
211 * specifies which atoms are involved in an interaction of a certain
212 * type. The layout of this array is as follows:
214 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
215 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
216 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
218 * So for interaction type type1 3 atoms are needed, and for type2 and
219 * type3 only 2. The type identifier is used to select the function to
220 * calculate the interaction and its actual parameters. This type
221 * identifier is an index in a params[] and functype[] array.
226 real
*cmap
; /* Has length 4*grid_spacing*grid_spacing, */
227 /* there are 4 entries for each cmap type (V,dVdx,dVdy,d2dVdxdy) */
232 int ngrid
; /* Number of allocated cmap (cmapdata_t ) grids */
233 int grid_spacing
; /* Grid spacing */
234 cmapdata_t
*cmapdata
; /* Pointer to grid with actual, pre-interpolated data */
242 t_functype
*functype
;
244 double reppow
; /* The repulsion power for VdW: C12*r^-reppow */
245 real fudgeQQ
; /* The scaling factor for Coulomb 1-4: f*q1*q2 */
246 gmx_cmap_t cmap_grid
; /* The dihedral correction maps */
250 ilsortUNKNOWN
, ilsortNO_FE
, ilsortFE_UNSORTED
, ilsortFE_SORTED
257 t_functype
*functype
;
260 gmx_cmap_t cmap_grid
;
261 t_iparams
*iparams_posres
;
262 int iparams_posres_nalloc
;
269 * The struct t_idef defines all the interactions for the complete
270 * simulation. The structure is setup in such a way that the multinode
271 * version of the program can use it as easy as the single node version.
272 * General field description:
274 * defines the number of elements in functype[] and param[].
276 * the node id (if parallel machines)
278 * the number of atomtypes
279 * t_functype *functype
280 * array of length ntypes, defines for every force type what type of
281 * function to use. Every "bond" with the same function but different
282 * force parameters is a different force type. The type identifier in the
283 * forceatoms[] array is an index in this array.
285 * array of length ntypes, defines the parameters for every interaction
286 * type. The type identifier in the actual interaction list
287 * (bondeds.iatoms[] or shakes.iatoms[]) is an index in this array.
289 * The list of interactions for each type. Note that some,
290 * such as LJ and COUL will have 0 entries.
294 int n
; /* n+1 is the number of points */
295 real scale
; /* distance between two points */
296 real
*tab
; /* the actual tables, per point there are 4 numbers */