4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
12 * Copyright (c) 1991-1999
13 * BIOSON Research Institute, Dept. of Biophysical Chemistry
14 * University of Groningen, The Netherlands
17 * GROMACS: A message-passing parallel molecular dynamics implementation
18 * H.J.C. Berendsen, D. van der Spoel and R. van Drunen
19 * Comp. Phys. Comm. 91, 43-56 (1995)
21 * Also check out our WWW page:
22 * http://md.chem.rug.nl/~gmx
27 * Green Red Orange Magenta Azure Cyan Skyblue
34 /* check kernel/toppush.c when you change these numbers */
36 #define MAXFORCEPARAM 6
39 typedef atom_id t_iatom
;
41 /* this MUST correspond to the
42 t_interaction_function[F_NRE] in gmxlib/ifunc.c */
81 F_NRE
/* This number is for the total number of energies */
86 /* Some parameters have A and B values for free energy calculations.
87 * The B values are not used for regular simulations of course.
88 * Free Energy for nonbondeds can be computed by changing the atom type.
89 * The harmonic type is used for all harmonic potentials:
90 * bonds, angles and improper dihedrals
92 struct {real a
,b
,c
; } bham
;
93 struct {real rA
,krA
,rB
,krB
; } harmonic
;
94 /* No free energy supported for WPOL */
95 struct {real kx
,ky
,kz
,rOH
,rHH
,rOD
; } wpol
;
96 struct {real c6
,c12
; } lj
;
97 struct {real c6A
,c12A
,c6B
,c12B
; } lj14
;
98 /* Proper dihedrals can not have different multiplicity when
99 * doing free energy calculations, because the potential would not
100 * be periodic anymore.
102 struct {real phiA
,cpA
;int mult
;real phiB
,cpB
; } pdihs
;
103 struct {real dA
,dB
; } shake
;
104 /* Settle can not be used for Free energy calculations.
105 * Use shake (or lincs) instead.
106 * The rest of the things cannot (yet) be used in FEP studies either.
108 struct {real doh
,dhh
; } settle
;
109 struct {real b0
,cb
,beta
; } morse
;
110 struct {real pos0
[DIM
],fc
[DIM
]; } posres
;
111 struct {real rbc
[NR_RBDIHS
]; } rbdihs
;
112 struct {real a
,b
,c
,d
,e
,f
; } dummy
;
113 struct {real low
,up1
,up2
,fac
;int type
,index
; } disres
;
114 struct {real buf
[MAXFORCEPARAM
]; } generic
; /* Conversion */
117 typedef int t_functype
;
122 int multinr
[MAXPROC
];
127 * The struct t_ilist defines a list of atoms with their interactions.
128 * General field description:
130 * the size (nr elements) of the interactions array (iatoms[]). This
131 * equals multinr[MAXPROC-1].
132 * int multinr[MAXPROC]
133 * specifies the number of type and atom id sequences in the iatoms[]
134 * array. Every element specifies the index of the first interaction
135 * for the next processor. The first processor starts at zero. So for
136 * n=0, the interactions run from 0 upto multinr[0]. The interactions
137 * for processor n (n>0) run from multinr[n-1] to index[n] (not including
140 * specifies which atoms are involved in an interaction of a certain
141 * type. The layout of this array is as follows:
143 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
144 * |type1|at1|at2|at3|type2|at1|at2|type1|at1|at2|at3|type3|at1|at2|
145 * +-----+---+---+---+-----+---+---+-----+---+---+---+-----+---+---+...
147 * So for interaction type type1 3 atoms are needed, and for type2 and
148 * type3 only 2. The type identifier is used to select the function to
149 * calculate the interaction and its actual parameters. This type
150 * identifier is an index in a params[] and functype[] array.
151 * The multinr[] array will be initialised for MAXPROC in such a way that up
152 * to the actual number of processors (during creation time), the array is
153 * filled with the indices, the remaining processors get empty parts by
154 * setting the indices to the largest value. In that way it is always possible
155 * to run this system on a larger multiprocessor ring however only the
156 * configured number of processors will we used. Running on less processors
157 * than configured is also possible by taking together adjacent
158 * processors. Note that in this case the load balance might get worse.
159 * The single processor version is implemented by simply using the complete
160 * configuration as one piece.
168 t_functype
*functype
;
175 * The struct t_idef defines all the interactions for the complete
176 * simulation. The structure is setup in such a way that the multiprocessor
177 * version of the program can use it as easy as the single processor version.
178 * General field description:
180 * defines the number of elements in functype[] and param[].
182 * the processor id (if parallel machines)
184 * the number of atomtypes
185 * t_functype *functype
186 * array of length ntypes, defines for every force type what type of
187 * function to use. Every "bond" with the same function but different
188 * force parameters is a different force type. The type identifier in the
189 * forceatoms[] array is an index in this array.
191 * array of length ntypes, defines the parameters for every interaction
192 * type. The type identifier in the actual interaction list
193 * (bondeds.iatoms[] or shakes.iatoms[]) is an index in this array.
195 * The list of interactions for each type. Note that some,
196 * such as LJ and COUL will have 0 entries.