9336933c400ec67f71e522edd07ff2031e786a74
[alexxy/gromacs.git] / docs / manual / topology.tex
1 %
2 % This file is part of the GROMACS molecular simulation package.
3 %
4 % Copyright (c) 2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
5 % Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 % and including many others, as listed in the AUTHORS file in the
7 % top-level source directory and at http://www.gromacs.org.
8 %
9 % GROMACS is free software; you can redistribute it and/or
10 % modify it under the terms of the GNU Lesser General Public License
11 % as published by the Free Software Foundation; either version 2.1
12 % of the License, or (at your option) any later version.
13 %
14 % GROMACS is distributed in the hope that it will be useful,
15 % but WITHOUT ANY WARRANTY; without even the implied warranty of
16 % MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17 % Lesser General Public License for more details.
18 %
19 % You should have received a copy of the GNU Lesser General Public
20 % License along with GROMACS; if not, see
21 % http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 % Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23 %
24 % If you want to redistribute modifications to GROMACS, please
25 % consider that scientific software is very special. Version
26 % control is crucial - bugs must be traceable. We will be happy to
27 % consider code for inclusion in the official distribution, but
28 % derived work must not be called official GROMACS. Details are found
29 % in the README & COPYING files - if they are missing, get the
30 % official version at http://www.gromacs.org.
31 %
32 % To help us fund GROMACS development, we humbly ask that you cite
33 % the research papers on the package. Check out http://www.gromacs.org.
34
35 \chapter{Topologies}
36 \label{ch:top}
37 \section{Introduction}
38 {\gromacs} must know on which atoms and combinations of atoms the
39 various contributions to the potential functions (see
40 \chref{ff}) must act. It must
41 also know what \normindex{parameter}s must be applied to the various
42 functions. All this is described in the {\em \normindex{topology}} file
43 {\tt *.top}, which lists the {\em constant attributes} of each atom.
44 There are many more atom types than elements, but only atom types
45 present in biological systems are parameterized in the force field,
46 plus some metals, ions and silicon. The bonded and special
47 interactions are determined by fixed lists that are included in the
48 topology file. Certain non-bonded interactions must be excluded (first
49 and second neighbors), as these are already treated in bonded
50 interactions.  In addition, there are {\em dynamic attributes} of
51 atoms - their positions, velocities and forces. These do not
52 strictly belong to the molecular topology, and are stored in the
53 coordinate file {\tt *.gro} (positions and velocities), or trajectory
54 file {\tt *.trr} (positions, velocities, forces).
55
56 This chapter describes the setup of the topology file, the
57 {\tt *.top} file and the database files: what the parameters
58 stand for and how/where to change them if needed.
59 First, all file formats are explained.
60 Section \ssecref{fffiles} describes the organization of
61 the files in each force field.
62
63 {\bf Note:} if you construct your own topologies, we encourage you
64 to upload them to our topology archive at {\wwwpage}! Just imagine
65 how thankful you'd have been if your topology had been available
66 there before you started. The same goes for new force fields or
67 modified versions of the standard force fields - contribute them
68 to the force field archive!
69
70 \section{Particle type}
71 \label{sec:parttype}
72
73 In {\gromacs}, there are three types of \normindex{particle}s, see
74 \tabref{ptype}. Only regular atoms and virtual interaction sites are used
75 in {\gromacs}; shells are necessary for
76 polarizable models like the Shell-Water models~\cite{Maaren2001a}.
77
78 \begin{table}
79 \centerline{
80 \begin{tabular}{|l|c|}
81 \dline
82 Particle                        & Symbol        \\
83 \hline
84 \seeindex{atom}{particle}s      & A   \\
85 \seeindex{shell}{particle}s     & S   \\
86 \normindex{virtual interaction sites}   & V (or D)   \\
87 \dline
88 \end{tabular}
89 }
90 \caption{Particle types in {\gromacs}}
91 \label{tab:ptype}
92 \end{table}
93
94 \subsection{Atom types}
95 \label{subsec:atomtype}
96
97 Each force field defines a set of \swapindex{atom}{type}s,
98 which have a characteristic name or number, and mass (in
99 a.m.u.). These listings are found in the {\tt atomtypes.atp}
100 file (.atp = {\bf a}tom {\bf t}ype {\bf p}arameter file).
101 Therefore, it is in this file that you can begin to change
102 and/or add an atom type. A sample from the {\tt gromos43a1.ff} 
103 force field is listed below.
104
105 {\small
106 \begin{verbatim}
107     O  15.99940 ;     carbonyl oxygen (C=O)
108    OM  15.99940 ;     carboxyl oxygen (CO-)
109    OA  15.99940 ;     hydroxyl, sugar or ester oxygen
110    OW  15.99940 ;     water oxygen
111     N  14.00670 ;     peptide nitrogen (N or NH)
112    NT  14.00670 ;     terminal nitrogen (NH2)
113    NL  14.00670 ;     terminal nitrogen (NH3)
114    NR  14.00670 ;     aromatic nitrogen
115    NZ  14.00670 ;     Arg NH (NH2)
116    NE  14.00670 ;     Arg NE (NH)
117     C  12.01100 ;     bare carbon
118   CH1  13.01900 ;     aliphatic or sugar CH-group
119   CH2  14.02700 ;     aliphatic or sugar CH2-group
120   CH3  15.03500 ;     aliphatic CH3-group
121 \end{verbatim}}
122
123 {\bf Note:} {\gromacs} makes use of the atom types as a name, {\em
124 not} as a number (as {\eg} in {\gromos}).
125
126 %Atomic detail is used except for hydrogen atoms bound to (aliphatic)
127 %carbon atoms, which are treated as {\em \swapindex{united}{atom}s}. No
128 %special \normindex{hydrogen-bond} term is included. {\bf Note} that other force field
129 %like OPLS/AA, CHARMM, and AMBER use all atoms.
130
131 %\subsection{Nucleus}
132 %{\em Necessary for \normindex{polarisability}, not implemented yet.}
133 %
134 %\subsection{Shell}
135 %{\em Necessary for polarisability, not implemented yet.}
136 %
137 %\subsection{Bond shell}
138 %{\em Necessary for polarisability, not implemented yet.}
139
140 \subsection{Virtual sites}
141 \label{sec:vsitetop}
142 Some force fields use \normindex{virtual interaction sites}
143 (interaction sites that are constructed from other particle positions)
144 on which certain interactions are located
145 ({\eg} on benzene rings, to reproduce the correct
146 \normindex{quadrupole}). This is described in~\secref{virtual_sites}.
147
148 To make virtual sites in your system, you should include a section
149 {\tt [~virtual_sites?~]} (for backward compatibility the old name
150 {\tt [~dummies?~]} can also be used) in your topology file,
151 where the `{\tt ?}' stands
152 for the number constructing particles for the virtual site. This will be
153 `{\tt 2}' for type 2, `{\tt 3}' for types 3, 3fd, 3fad and 3out and
154 `{\tt 4}' for type 4fdn. The last of these replace an older 4fd type (with the `type' value 1) 
155 that could occasionally be unstable; while it is still supported internally
156 in the code, the old 4fd type should not be used in new input files.
157  The different types are explained
158 in~\secref{virtual_sites}.
159
160 Parameters for type 2 should look like this:
161 {\small
162 \begin{verbatim}
163 [ virtual_sites2 ]
164 ; Site  from        funct  a
165 5       1     2     1      0.7439756
166 \end{verbatim}}
167
168 for type 3 like this:
169 {\small
170 \begin{verbatim}
171 [ virtual_sites3 ]
172 ; Site  from               funct   a          b
173 5       1     2     3      1       0.7439756  0.128012
174 \end{verbatim}}
175
176 for type 3fd like this:
177 {\small
178 \begin{verbatim}
179 [ virtual_sites3 ]
180 ; Site  from               funct   a          d
181 5       1     2     3      2       0.5        -0.105
182 \end{verbatim}}
183
184 for type 3fad like this:
185 {\small
186 \begin{verbatim}
187 [ virtual_sites3 ]
188 ; Site  from               funct   theta      d
189 5       1     2     3      3       120        0.5
190 \end{verbatim}}
191
192 for type 3out like this:
193 {\small
194 \begin{verbatim}
195 [ virtual_sites3 ]
196 ; Site  from               funct   a          b          c
197 5       1     2     3      4       -0.4       -0.4       6.9281
198 \end{verbatim}}
199
200 for type 4fdn like this:
201 {\small
202 \begin{verbatim}
203 [ virtual_sites4 ]
204 ; Site  from                      funct   a          b          c
205 5       1     2     3     4       2       1.0        0.9       0.105
206 \end{verbatim}}
207
208 This will result in the construction of a virtual site, number 5
209 (first column `{\tt Site}'), based on the positions of the atoms
210 whose indices are 1 and 2 or 1, 2 and 3 or 1, 2, 3 and 4 (next two,
211 three or four columns `{\tt from}') following the rules determined by the function number
212 (next column `{\tt funct}') with the parameters specified (last one,
213 two or three columns `{\tt a b} . .'). Obviously, the atom numbers
214 (including virtual site number) depend
215 on the molecule. It may be instructive to study the topologies for
216 TIP4P or TIP5P water models that are included with the {\gromacs} distribution.
217
218 {\bf Note} that if any constant bonded interactions are defined between
219 virtual sites and/or normal atoms, they will be removed by {\tt grompp}
220 (unless the option {tt -normvsbds} is used).
221 This removal of bonded interactions is done after generating exclusions,
222 as the generation of exclusions is based on ``chemically'' bonded interactions.
223
224 Virtual sites can be constructed in a more generic way using basic geometric
225 parameters.  The directive that can be used is {\tt [ virtual_sitesn ]}.  Required
226 parameters are listed in~\tabref{topfile2}.  An example entry for defining a virtual
227 site at the center of geometry of a given set of atoms might be:
228
229 {\small
230 \begin{verbatim}
231 [ virtual_sitesn ]
232 ; Site   funct    from
233 5        1        1     2     3     4
234 \end{verbatim}}
235
236 \section{Parameter files}
237
238 \label{sec:paramfiles}
239
240 \subsection{Atoms}
241 The {\em static} properties (see \tabref{statprop} assigned to the
242 atom types are assigned based on data in several places.
243 The mass is listed in {\tt atomtypes.atp}
244 (see~\ssecref{atomtype}), whereas the charge is listed in {\tt *.rtp}
245 (.rtp = {\bf r}esidue {\bf t}opology {\bf p}arameter file,
246 see~\ssecref{rtp}).  This implies that the charges are only defined in
247 the \normindex{building block}s of amino acids, nucleic acids or
248 otherwise, as defined by the user. When generating a topology
249 ({\tt *.top}) using the {\tt \normindex{pdb2gmx}} program, the
250 information from these files is combined.
251  
252 \begin{table}
253 \centerline{
254 \begin{tabular}{|l|c|c|}
255 \dline
256 Property        & Symbol        & Unit          \\
257 \hline
258 Type            & -             & -             \\
259 Mass            & m             & a.m.u.        \\
260 Charge          & q             & electron      \\
261 epsilon         & $\epsilon$    & kJ/mol        \\
262 sigma           & $\sigma$      & nm            \\
263 \dline
264 \end{tabular}
265 }
266 \caption{Static atom type properties in {\gromacs}}
267 \label{tab:statprop}
268 \end{table}
269
270 %The following {\em dynamic} quantities are associated with an atom
271 %\begin{itemize}
272 %\item   Position {\bf x}
273 %\item   Velocity {\bf v}
274 %\end{itemize}
275 %These quantities are listed in the coordinate file, {\tt *.gro}
276 %(see section~\ssecref{grofile}).
277
278 \subsection{Non-bonded parameters}
279 \label{subsec:nbpar}
280 The \swapindex{non-bonded}{parameter}s consist of the van der Waals 
281 parameters V ({\tt c6} or $\sigma$, depending on the combination rule) and W 
282 ({\tt c12} or $\epsilon$), as listed in the file {\tt ffnonbonded.itp}, where 
283 {\tt ptype} is the particle type (see \tabref{ptype}). As with the bonded parameters, entries in {\tt [~*type~]} directives 
284 are applied to their counterparts in the topology file.  Missing parameters 
285 generate warnings, except as noted below in section~\ssecref{pairinteractions}.
286
287 {\small
288 \begin{verbatim}
289 [ atomtypes ]
290 ;name   at.num      mass      charge   ptype         V(c6)        W(c12)
291     O        8  15.99940       0.000       A   0.22617E-02   0.74158E-06
292    OM        8  15.99940       0.000       A   0.22617E-02   0.74158E-06
293    .....
294
295 [ nonbond_params ]
296   ; i    j func       V(c6)        W(c12)
297     O    O    1 0.22617E-02   0.74158E-06
298     O   OA    1 0.22617E-02   0.13807E-05
299     .....
300 \end{verbatim}}
301
302 {\bf Note} that most of the included force fields also include the {\tt at.num.} column, 
303 but this same information is implied in the OPLS-AA {\tt bond_type} column.
304 The interpretation of the parameters V and W depends on the combination rule 
305 that was chosen in the {\tt [~defaults~]} section of the topology file 
306 (see~\ssecref{topfile}):
307 \begin{eqnarray}
308 \mbox{for combination rule 1}: & &
309 \begin{array}{llllll}
310   \mbox{V}_{ii} & = & C^{(6)}_{i}  & = & 4\,\epsilon_i\sigma_i^{6} &
311   \mbox{[ kJ mol$^{-1}$ nm$^{6}$ ]}\\
312   \mbox{W}_{ii} & = & C^{(12)}_{i} & = & 4\,\epsilon_i\sigma_i^{12} &
313   \mbox{[ kJ mol$^{-1}$ nm$^{12}$ ]}\\
314 \end{array}
315 \\
316 \mbox{for combination rules 2 and 3}: & &
317 \begin{array}{llll}
318   \mbox{V}_{ii} & = & \sigma_i   & \mbox{[ nm ]} \\
319   \mbox{W}_{ii} & = & \epsilon_i & \mbox{[ kJ mol$^{-1}$ ]}
320 \end{array}
321 \end{eqnarray}
322 Some or all combinations for different atom types can be given in the 
323 {\tt [~nonbond_params~]} section, again with parameters V and W as defined 
324 above. Any combination that is not given will be computed from the parameters 
325 for the corresponding atom types, according to the \normindex{combination rule}:
326 \begin{eqnarray}
327 \mbox{for combination rules 1 and 3}: & &
328 \begin{array}{lll}
329   C^{(6)}_{ij}  & = & \left(C^{(6)}_i\,C^{(6)}_j\right)^{\frac{1}{2}} \\
330   C^{(12)}_{ij} & = & \left(C^{(12)}_i\,C^{(12)}_j\right)^{\frac{1}{2}}
331 \end{array}
332 \\
333 \mbox{for combination rule 2}: & &
334 \begin{array}{lll}
335   \sigma_{ij}   & = & \frac{1}{2}(\sigma_i+\sigma_j) \\
336   \epsilon_{ij} & = & \sqrt{\epsilon_i\,\epsilon_j}
337 \end{array}
338 \end{eqnarray}
339 When $\sigma$ and $\epsilon$ need to be supplied (rules 2 and 3),
340 it would seem it is impossible to have a non-zero $C^{12}$ combined
341 with a zero $C^6$ parameter. However, providing a negative $\sigma$
342 will do exactly that, such that $C^6$ is set to zero and $C^{12}$ is
343 calculated normally. This situation represents a special case in reading
344 the value of $\sigma$, and nothing more.
345
346 There is only one set of \normindex{combination rule}s
347 for Buckingham potentials:
348 \beq
349 \begin{array}{rcl}
350 A_{ij}   &=& \left(A_{ii} \, A_{jj}\right)^{1/2}    \\
351 B_{ij}   &=& 2 / \left(\frac{1}{B_{ii}} + \frac{1}{B_{jj}}\right)        \\
352 C_{ij}   &=& \left(C_{ii} \, C_{jj}\right)^{1/2}
353 \end{array}
354 \eeq
355
356 \subsection{Bonded parameters}
357 \label{subsec:bondparam}
358 The \swapindex{bonded}{parameter}s ({\ie} bonds, bond angles, improper and proper
359 dihedrals) are listed in {\tt ffbonded.itp}.~
360 % The term {\tt func} is 1 for
361 % harmonic and 2 for \gromosv{96} bond and angle potentials.
362 % For the dihedral, this is explained after this listing.
363 The entries in this database describe, respectively, the atom types
364 in the interactions, the type of the interaction, and the parameters
365 associated with that interaction. These parameters are then read
366 by {\tt \normindex{grompp}} when processing a topology and applied
367 to the relevant bonded parameters, {\ie} {\tt bondtypes} are applied to
368 entries in the {\tt [~bonds~]} directive, etc. Any bonded parameter that is
369 missing from the relevant {\tt [~*type~]} directive generates a fatal error.
370 The types of interactions are listed in \tabref{topfile2}.
371 Example excerpts from such files follow:
372
373 {\small 
374 \begin{verbatim}
375 [ bondtypes ]
376   ; i    j func        b0          kb
377     C    O    1   0.12300     502080.
378     C   OM    1   0.12500     418400.
379     ......
380
381 [ angletypes ]
382   ; i    j    k func       th0         cth
383    HO   OA    C    1   109.500     397.480
384    HO   OA  CH1    1   109.500     397.480
385    ......
386
387 [ dihedraltypes ]
388   ; i    l func        q0          cq
389  NR5*  NR5    2     0.000     167.360
390  NR5* NR5*    2     0.000     167.360
391  ......
392
393 [ dihedraltypes ]
394   ; j    k func      phi0          cp   mult
395     C   OA    1   180.000      16.736      2
396     C    N    1   180.000      33.472      2
397     ......
398
399 [ dihedraltypes ]
400 ;
401 ; Ryckaert-Bellemans Dihedrals
402 ;
403 ; aj    ak      funct
404 CP2     CP2     3       9.2789  12.156  -13.120 -3.0597 26.240  -31.495
405 \end{verbatim}}
406
407 In the {\tt ffbonded.itp} file, you can add bonded parameters. If you
408 want to include parameters for new atom types, make sure you define
409 them in {\tt atomtypes.atp} as well.
410
411 For most interaction types, bonded parameters are searched and assigned
412 using an exact match for all type names and allowing only a single set
413 of parameters. The exception to this rule are \normindex{dihedral} parameters.
414 For {\tt [ dihedraltypes ]} wildcard atom type names can be specified
415 with the letter {\tt X} in one or more of the four positions. Thus one
416 can for example assign proper dihedral parameters based on the types of
417 the middle two atoms. The parameters for the entry with the most exact matches,
418 i.e. the least wildcard matches, will be used. Note that {\gromacs} versions
419 older than 5.1.3 used the first match, which means that a full match would be
420 ignored if it is preceded by an entry that matches on wildcards. Thus it
421 is suggested to put wildcard entries at the end, in case someone might
422 use a forcefield with older versions of {\gromacs}.
423 In addition there is a dihedral type 9 which adds the possibility of
424 assigning multiple dihedral potentials, useful for combining terms with
425 different multiplicities. The different dihedral potential parameter
426 sets should be on directly adjacent lines in the {\tt [ dihedraltypes ]}
427 section.
428
429 \section{Molecule definition\index{molecule definition}}
430
431 \subsection{Moleculetype entries}
432 An organizational structure that usually corresponds to molecules is
433 the {\tt [ moleculetype ]} entry. This entry serves two main purposes. One is
434 to give structure to the topology file(s), usually corresponding to real
435 molecules. This makes the topology easier to read and writing it less labor
436 intensive. A second purpose is computational efficiency. The system definition
437 that is kept in memory is proportional in size of the {\tt moleculetype}
438 definitions. If a molecule is present in 100000 copies, this saves a factor
439 of 100000 in memory, which means the system usually fits in cache, which
440 can improve performance tremendously. Interactions that correspond to chemical
441 bonds, that generate exclusions, can only be defined between atoms within
442 a {\tt moleculetype}. It is allowed to have multiple molecules which are
443 not covalently bonded in one {\tt moleculetype} definition. Molecules can
444 be made infinitely long by connecting to themselves over periodic boundaries.
445 When such periodic molecules are present, an option in the {\tt mdp} file
446 needs to be set to tell {\gromacs} not to attempt to make molecules
447 that are broken over periodic boundaries whole again.
448
449 \subsection{Intermolecular interactions\index{intermolecular interaction}}
450 In some cases, one would like atoms in different molecules to also interact
451 with other interactions than the usual non-bonded interactions. This is often
452 the case in binding studies. When the molecules are covalently bound, e.g.
453 a ligand binding covalently to a protein, they are effectively one molecule
454 and they should be defined in one {\tt [ moleculetype ]} entry. Note that
455 {\tt pdb2gmx} has an option to put two or more molecules in one
456 {\tt [ moleculetype ]} entry. When molecules are not covalently bound,
457 it is much more convenient to use separate {\tt moleculetype} definitions
458 and specify the intermolecular interactions in the
459 {\tt [ intermolecular\_interactions] } section. In this section, which is
460 placed at the end of the topology (see \tabref{topfile1}), normal bonded
461 interactions can be specified using global atom indices. The only restrictions
462 are that no interactions can be used that generates exclusions and no
463 constraints can be used.
464
465 \subsection{Intramolecular pair interactions\index{intramolecular pair interaction}}
466 \label{subsec:pairinteractions}
467 Extra Lennard-Jones and electrostatic interactions between pairs
468 of atoms in a molecule can be added in the {\tt [~pairs~]} section of
469 a molecule definition. The parameters for these interactions can
470 be set independently from the non-bonded interaction parameters.
471 In the {\gromos} force fields, pairs are only used
472 to modify the \normindex{1-4 interaction}s (interactions of atoms
473 separated by three bonds). In these force fields the 1-4 interactions
474 are excluded from the non-bonded interactions (see \secref{excl}).
475
476 {\small
477 \begin{verbatim}
478
479 [ pairtypes ]
480   ; i    j func         cs6          cs12 ; THESE ARE 1-4 INTERACTIONS
481     O    O    1 0.22617E-02   0.74158E-06
482     O   OM    1 0.22617E-02   0.74158E-06
483     .....
484 \end{verbatim}}
485
486 The pair interaction parameters for the atom types
487 in {\tt ffnonbonded.itp} are listed in the {\tt [~pairtypes~]} section.
488 The {\gromos} force fields list all these interaction parameters
489 explicitly, but this section might be empty for force fields like
490 OPLS that calculate the \normindex{1-4 interaction}s by uniformly scaling the parameters.
491 Pair parameters that are not present in the {\tt [~pairtypes~]} section
492 are only generated when {\tt gen-pairs} is set to ``yes'' in the {\tt [~defaults~]}
493 directive of {\tt forcefield.itp} (see \ssecref{topfile}). 
494 When {\tt gen-pairs} is set to ``no,'' {\tt \normindex{grompp}}
495 will give a warning for each pair type for which no parameters are given.
496
497 The normal pair interactions, intended for \normindex{1-4 interaction}s,
498 have function type 1. Function type 2 and the {\tt [~pairs_nb~]} are intended
499 for free-energy simulations. When determining hydration
500 free energies, the solute needs to be decoupled from the solvent.
501 This can be done by adding a B-state topology (see \secref{fecalc})
502 that uses zero for all solute non-bonded parameters, {\ie} charges and LJ parameters.
503 However, the free energy difference between the A and
504 B states is not the total hydration free energy.  One has to
505 add the free energy for reintroducing the internal Coulomb and 
506 LJ interactions in the solute when in vacuum. This second step can be combined with
507 the first step when the Coulomb and LJ interactions within
508 the solute are not modified. For this purpose, there is a pairs
509 function type 2, which is identical to function type 1, except
510 that the B-state parameters are always identical to the A-state
511 parameters. For searching the parameters in the {\tt [~pairtypes~]} section,
512 no distinction is made between function type 1 and 2.
513 The pairs section {\tt [~pairs_nb~]} is intended to replace the non-bonded interaction.
514 It uses the unscaled charges and the non-bonded LJ parameters;
515 it also only uses the A-state parameters. {\bf Note} that
516 one should add exclusions for all atom pairs listed in {\tt [~pairs_nb~]},
517 otherwise such pairs will also end up in the normal neighbor lists.
518
519 Alternatively, this same behavior can be achieved without ever
520 touching the topology, by using the {\tt couple-moltype}, {\tt
521   couple-lambda0}, {\tt couple-lambda1}, and {\tt couple-intramol}
522 keywords.  See sections \secref{fecalc} and \secref{dgimplement} for
523 more information.
524
525 All three pair types always use plain Coulomb interactions,
526 even when Reaction-field, PME, Ewald or shifted Coulomb interactions
527 are selected for the non-bonded interactions.
528 Energies for types 1 and 2 are written to the energy and log file
529 in separate ``LJ-14'' and ``Coulomb-14'' entries per energy group pair.
530 Energies for {\tt [~pairs_nb~]} are added to the ``LJ-(SR)'' and ``Coulomb-(SR)'' terms.
531
532
533 \subsection{Exclusions}
534 \label{sec:excl}
535 The \normindex{exclusions} for non-bonded interactions are generated by {\tt
536 grompp} for neighboring atoms up to a certain number of bonds away, as
537 defined in the {\tt [~moleculetype~]} section in the topology file
538 (see \ssecref{topfile}). Particles are considered bonded when they are
539 connected by ``chemical'' bonds ({\tt [~bonds~]} types 1 to 5, 7 or 8)
540 or constraints ({\tt [~constraints~]} type 1).
541 Type 5 {\tt [~bonds~]} can be used to create a \normindex{connection}
542 between two atoms without creating an interaction.
543 There is a \normindex{harmonic interaction}
544 ({\tt [~bonds~]} type 6) that does not connect the atoms by a chemical bond.
545 There is also a second constraint type ({\tt [~constraints~]} type 2)
546 that fixes the distance, but does not connect
547 the atoms by a chemical bond.
548 For a complete list of all these interactions, see \tabref{topfile2}.
549
550 Extra exclusions within a molecule can be added manually
551 in a {\tt [~exclusions~]} section. Each line should start with one
552 atom index, followed by one or more atom indices. All non-bonded
553 interactions between the first atom and the other atoms will be excluded.
554
555 When all non-bonded interactions within or between groups of atoms need
556 to be excluded, is it more convenient and much more efficient to use
557 energy monitor group exclusions (see \secref{groupconcept}).
558
559
560 \section{Implicit solvation parameters\index{implicit solvation parameters}}
561 Starting with {\gromacs} 4.5, implicit solvent is supported. A section in the
562 topology has been introduced to list those parameters:
563
564 {\small
565 \begin{verbatim}
566 [ implicit_genborn_params ]
567 ; Atomtype  sar     st   pi      gbr      hct
568 NH1         0.155   1    1.028   0.17063  0.79 ; N
569 N           0.155   1    1       0.155    0.79 ; Proline backbone N
570 H           0.1     1    1       0.115    0.85 ; H
571 CT1         0.180   1    1.276   0.190    0.72 ; C
572 \end{verbatim}}
573
574 In this example the atom type is listed first, followed by five
575 numbers, and a comment (following a semicolon).
576
577 Values in columns 1-3 are not currently used. They pertain to more
578 elaborate surface area algorithms, the one from Qiu {\em et al.}~\cite{Still97} in
579 particular.  Column 4 contains the atomic van der Waals radii, which are used
580 in computing the Born radii. The dielectric offset is specified in
581 the {\tt *.mdp} file, and gets subtracted from the input van der Waals radii for the different
582 Born radii methods, as described by Onufriev {\em et al.}~\cite{Case04}.  Column 5 is the 
583 scale factor for the HCT and OBC models. The values are taken from the Tinker implementation of 
584 the HCT pairwise scaling method~\cite{Truhlar96}.  This method has been modified such that the
585 scaling factors have been adjusted to minimize differences between analytical surface areas and
586 GB using the HCT algorithm.  The scaling is further modified in that it is not applied pairwise
587 as proposed by Hawkins {\em et al.}~\cite{Truhlar96}, but on a per-atom (rather than a per-pair) 
588 basis.
589
590
591 \section{Constraint algorithms\index{constraint algorithms}}
592 \label{sec:constraints}
593 Constraints are defined in the {\tt [~constraints~]} section.
594 The format is two atom numbers followed by the function type,
595 which can be 1 or 2, and the constraint distance.
596 The only difference between the two types is that type 1 is used
597 for generating exclusions and type 2 is not (see \secref{excl}).
598 The distances are constrained using the LINCS or the SHAKE algorithm,
599 which can be selected in the {\tt *.mdp} file.
600 Both types of constraints can be perturbed in free-energy calculations
601 by adding a second constraint distance (see \ssecref{constraintforce}).
602 Several types of bonds and angles (see \tabref{topfile2}) can
603 be converted automatically to constraints by {\tt grompp}.
604 There are several options for this in the {\tt *.mdp} file.
605
606 We have also implemented the \normindex{SETTLE} algorithm~\cite{Miyamoto92},
607 which is an analytical solution of SHAKE, specifically for water. 
608 SETTLE can be selected in the topology file. See, for instance, the
609 SPC molecule definition:
610
611 {\small
612 \begin{verbatim}
613 [ moleculetype ]
614 ; molname       nrexcl
615 SOL             1
616
617 [ atoms ]
618 ; nr    at type res nr  ren nm  at nm   cg nr   charge
619 1       OW      1       SOL     OW1     1       -0.82
620 2       HW      1       SOL     HW2     1        0.41
621 3       HW      1       SOL     HW3     1        0.41
622
623 [ settles ]
624 ; OW    funct   doh     dhh
625 1       1       0.1     0.16333
626
627 [ exclusions ]
628 1       2       3
629 2       1       3
630 3       1       2
631 \end{verbatim}}
632
633 The {\tt [~settles~]} directive defines the first atom of the water molecule.
634 The settle funct is always 1, and the distance between O-H and H-H distances
635 must be given. {\bf Note} that the algorithm can also be used
636 for TIP3P and TIP4P~\cite{Jorgensen83}.
637 TIP3P just has another geometry. TIP4P has a virtual site, but since 
638 that is generated it does not need to be shaken (nor stirred).
639
640 \section{\normindex{pdb2gmx} input files}
641 \label{sec:pdb2gmxfiles}
642 The {\gromacs} program {\tt pdb2gmx} generates a topology for
643 the input coordinate file. Several formats are supported for
644 that coordinate file, but {\tt *.pdb} is the most commonly-used format
645 (hence the name {\tt pdb2gmx}).
646 {\tt pdb2gmx} searches for force fields in sub-directories of the {\gromacs} {\tt share/top}
647 directory and your working directory. Force fields are recognized from
648 the file {\tt forcefield.itp} in a directory with the extension {\tt .ff}.
649 The file {\tt forcefield.doc} may be present, and if so, its first line
650 will be used by {\tt pdb2gmx} to present a short description to the
651 user to help in choosing a force field. Otherwise, the user can
652 choose a force field with the {\tt -ff xxx} command-line argument
653 to {\tt pdb2gmx}, which indicates that a force field in a
654 {\tt xxx.ff} directory is desired. {\tt pdb2gmx} will search first in the
655 working directory, then in the {\gromacs} {\tt share/top} directory, and
656 use the first matching {\tt xxx.ff} directory found.
657
658 Two general files are read by {\tt pdb2gmx}: an atom type file
659 (extension {\tt .atp}, see~\ssecref{atomtype}) from the force-field directory,
660 and a file called {\tt residuetypes.dat} from either the working directory, or
661 the {\gromacs} {\tt share/top} directory. {\tt residuetypes.dat}
662 determines which residue names are considered protein, DNA, RNA,
663 water, and ions.
664
665 {\tt pdb2gmx} can read one or multiple databases with topological information
666 for different types of molecules. A set of files belonging to one database
667 should have the same basename, preferably telling something about the type
668 of molecules ({\eg} aminoacids, rna, dna). The possible files are:
669 \begin{itemize}
670 \item {\tt <basename>.rtp}
671 \item {\tt <basename>.r2b} (optional)
672 \item {\tt <basename>.arn} (optional)
673 \item {\tt <basename>.hdb} (optional)
674 \item {\tt <basename>.n.tdb} (optional)
675 \item {\tt <basename>.c.tdb} (optional)
676 \end{itemize}
677 Only the {\tt .rtp} file, which contains the topologies of the building
678 blocks, is mandatory. Information from other files will only be used 
679 for building blocks that come from an {\tt .rtp} file with the same base name.
680 The user can add building blocks to a force field by having additional
681 files with the same base name in their working directory. By default, only
682 extra building blocks can be defined, but calling {\tt pdb2gmx} with
683 the {\tt -rtpo} option will allow building blocks in a local file
684 to replace the default ones in the force field.
685
686
687 \subsection{Residue database}
688 \label{subsec:rtp}
689 The files holding the residue databases have the extension {\tt .rtp}.
690 Originally this file contained building blocks (amino acids) for proteins,
691 and is the {\gromacs} interpretation of the {\tt rt37c4.dat} file of {\gromos}.
692 So the residue database file contains information (bonds, charges, charge groups,
693 and improper dihedrals) for a frequently-used building block. It is
694 better {\em not} to change this file because it is standard input for
695 {\tt pdb2gmx}, but if changes are needed make them in the
696 {\tt *.top} file (see~\ssecref{topfile}), or in a {\tt .rtp} file
697 in the working directory as explained in \secref{pdb2gmxfiles}.
698 Defining topologies of new small molecules is probably easier
699 by writing an include topology file {\tt *.itp} directly.
700 This will be discussed in section~\ssecref{molitp}.
701 When adding a new protein residue to the database, don't forget to
702 add the residue name to the {\tt \normindex{residuetypes.dat}} file,
703 so that {\tt grompp}, {\tt make_ndx} and analysis tools can recognize
704 the residue as a protein residue (see \ssecref{defaultgroups}).
705
706 The {\tt .rtp} files are only used by {\tt pdb2gmx}.
707 As mentioned before, the only extra information this
708 program needs from the {\tt .rtp} database is bonds, charges of atoms,
709 charge groups, and improper dihedrals, because the rest is read from
710 the coordinate input file.
711 Some proteins contain residues that are not standard, but are
712 listed in the coordinate file. You have to construct a building block
713 for this ``strange'' residue, otherwise you will not obtain a
714 {\tt *.top} file. This also holds for molecules in the
715 coordinate file such as ligands, polyatomic ions, crystallization co-solvents, etc.
716 The residue database is constructed in the following way:
717
718 {\small
719 \begin{verbatim}
720 [ bondedtypes ]  ; mandatory
721 ; bonds  angles  dihedrals  impropers
722      1       1          1          2  ; mandatory
723
724 [ GLY ]  ; mandatory
725
726  [ atoms ]  ; mandatory 
727 ; name  type  charge  chargegroup 
728      N     N  -0.280     0
729      H     H   0.280     0
730     CA   CH2   0.000     1
731      C     C   0.380     2
732      O     O  -0.380     2
733
734  [ bonds ]  ; optional
735 ;atom1 atom2      b0      kb
736      N     H
737      N    CA
738     CA     C
739      C     O
740     -C     N
741
742  [ exclusions ]  ; optional
743 ;atom1 atom2
744
745  [ angles ]  ; optional
746 ;atom1 atom2 atom3    th0    cth
747
748  [ dihedrals ]  ; optional
749 ;atom1 atom2 atom3 atom4   phi0     cp   mult
750
751  [ impropers ]  ; optional
752 ;atom1 atom2 atom3 atom4     q0     cq
753      N    -C    CA     H
754     -C   -CA     N    -O
755
756 [ ZN ]
757
758  [ atoms ]
759     ZN    ZN   2.000     0
760 \end{verbatim}}
761
762 The file is free format; the only restriction is that there can be at
763 most one entry on a line.  The first field in the file is the
764 {\tt [~bondedtypes~]} field, which is followed by four numbers,
765 indicating the interaction type for bonds, angles, dihedrals, and
766 improper dihedrals.  The file contains residue entries, which consist
767 of atoms and (optionally) bonds, angles, dihedrals, and impropers.  The
768 charge group codes denote the charge group numbers. Atoms in the same
769 charge group should always be ordered consecutively. When using the
770 hydrogen database with {\tt pdb2gmx} for adding missing hydrogens
771 (see~\ssecref{hdb}), the atom names defined in the {\tt .rtp} entry
772 should correspond exactly to the naming convention used in the
773 hydrogen database. The atom names in the bonded interaction can be
774 preceded by a minus or a plus, indicating that the atom is in the
775 preceding or following residue respectively.  Explicit parameters added
776 to bonds, angles, dihedrals, and impropers override
777 the standard parameters in the {\tt .itp} files.  This should only be
778 used in special cases. Instead of parameters, a string can be added
779 for each bonded interaction.  This is used in \gromosv{96} {\tt .rtp}
780 files. These strings are copied to the topology file and can be
781 replaced by force-field parameters by the C-preprocessor in {\tt grompp}
782 using {\tt \#define} statements.
783
784 {\tt pdb2gmx} automatically generates all angles. This means that for 
785 most force fields the {\tt [~angles~]} field is only useful for overriding 
786 {\tt .itp} parameters. For the \gromosv{96} force field the interaction 
787 number of all angles needs to be specified.
788
789 {\tt pdb2gmx} automatically generates one proper dihedral for every rotatable
790 bond, preferably on heavy atoms. When the {\tt [~dihedrals~]} field is used,
791 no other dihedrals will be generated for the bonds corresponding to the
792 specified  dihedrals. It is possible to put more than one dihedral
793 function on a rotatable bond. In the case of CHARMM27 FF {\tt pdb2gmx}
794 can add correction maps to the dihedrals using the default {\tt -cmap} option.
795 Please refer to \ssecref{charmmff} for more information.
796
797 {\tt pdb2gmx} sets the number of exclusions to 3, which
798 means that interactions between atoms connected by at most 3 bonds are
799 excluded. Pair interactions are generated for all pairs of atoms that are
800 separated by 3 bonds (except pairs of hydrogens).
801 When more interactions need to be excluded, or some pair interactions should
802 not be generated, an {\tt [~exclusions~]} field can be added, followed by
803 pairs of atom names on separate lines. All non-bonded and pair interactions
804 between these atoms will be excluded.
805
806 \subsection{Residue to building block database}
807 Each force field has its own naming convention for residues.
808 Most residues have consistent naming, but some, especially those
809 with different protonation states, can have many different names.
810 The {\tt .r2b} files are used to convert standard residue names to
811 the force-field build block names. If no {\tt .r2b} is present
812 in the force-field directory or a residue is not listed, the building
813 block name is assumed to be identical to the residue name.
814 The {\tt .r2b} can contain 2 or 5 columns. The 2-column format
815 has the residue name in the first column and the building block name
816 in the second. The 5-column format has 3 additional columns with
817 the building block for the residue occurring in the N-terminus, C-terminus
818 and both termini at the same time (single residue molecule).
819 This is useful for, for instance, the AMBER force fields.
820 If one or more of the terminal versions are not present, a dash should be entered
821 in the corresponding column.
822
823 There is a {\gromacs} naming convention for residues which is only
824 apparent (except for the {\tt pdb2gmx} code) through the {\tt .r2b} file
825 and {\tt specbond.dat} files.
826 This convention is only of importance when you are adding residue types
827 to an {\tt .rtp} file. The convention is listed in \tabref{r2b}.
828 For special bonds with, for instance, a heme group, the {\gromacs} naming
829 convention is introduced through {\tt specbond.dat} (see~\ssecref{specbond}), which can
830 subsequently be translated by the {\tt .r2b} file, if required.
831
832 \begin{table}
833 \centerline{
834 \begin{tabular}{|ll|}
835 \dline
836 ARG  & protonated arginine \\
837 ARGN & neutral arginine \\
838 ASP  & negatively charged aspartic acid \\
839 ASPH & neutral aspartic acid \\
840 CYS  & neutral cysteine \\
841 CYS2 & cysteine with sulfur bound to another cysteine or a heme \\
842 GLU  & negatively charged glutamic acid \\
843 GLUH & neutral glutamic acid \\
844 HISD & neutral histidine with N$_\delta$ protonated \\
845 HISE & neutral histidine with N$_\epsilon$ protonated \\
846 HISH & positive histidine with both N$_\delta$ and N$_\epsilon$ protonated \\
847 HIS1 & histidine bound to a heme \\
848 LYSN & neutral lysine \\
849 LYS  & protonated lysine \\
850 HEME & heme \\
851 \dline
852 \end{tabular}
853 }
854 \caption{Internal {\gromacs} residue naming convention.}
855 \label{tab:r2b}
856 \end{table}
857
858 \subsection{Atom renaming database}
859 Force fields often use atom names that do not follow IUPAC or PDB convention.
860 The {\tt .arn} database is used to translate the atom names in the coordinate
861 file to the force-field names. Atoms that are not listed keep their names.
862 The file has three columns: the building block name,
863 the old atom name, and the new atom name, respectively. The residue name
864 supports question-mark wildcards that match a single character.
865
866 An additional general atom renaming file called {\tt xlateat.dat} is present
867 in the {\tt share/top} directory, which translates common non-standard
868 atom names in the coordinate file to IUPAC/PDB convention. Thus, when writing
869 force-field files, you can assume standard atom names and no further
870 atom name translation is required, except for translating from standard atom names
871 to the force-field ones.
872
873 \subsection{Hydrogen database}
874 \label{subsec:hdb}
875 The \swapindex{hydrogen}{database} is stored in {\tt .hdb} files. It
876 contains information for the {\tt pdb2gmx} program on how to connect
877 hydrogen atoms to existing atoms. In versions of the database before
878 {\gromacs} 3.3, hydrogen atoms were named after the atom they are
879 connected to: the first letter of the atom name was replaced by an
880 `H.' In the versions from 3.3 onwards, the H atom has to be listed explicitly,
881 because the old behavior was protein-specific and hence could not
882 be generalized to other molecules.
883 If more than one hydrogen atom is connected to the same atom, a
884 number will be added to the end of the hydrogen atom name. For
885 example, adding two hydrogen atoms to \texttt{ND2} (in asparagine), the
886 hydrogen atoms will be named \texttt{HD21} and \texttt{HD22}. This is
887 important since atom naming in the \texttt{.rtp} file (see~\ssecref{rtp})
888 must be the same. The format of the hydrogen database is as follows:
889
890 {\small
891 \begin{verbatim}
892 ; res   # additions
893         # H add type    H       i       j       k
894 ALA     1
895         1       1       H       N       -C      CA
896 ARG     4
897         1       2       H       N       CA      C
898         1       1       HE      NE      CD      CZ
899         2       3       HH1     NH1     CZ      NE
900         2       3       HH2     NH2     CZ      NE
901 \end{verbatim}}
902
903 On the first line we see the residue name (ALA or ARG) and the number
904 of kinds of hydrogen atoms that may be added to this residue by the
905 hydrogen database. After that follows one line for each addition, on which
906 we see:
907 \begin{itemize}
908 \item The number of H atoms added
909 \item The method for adding H atoms, which can be any of:
910 \begin{enumerate}
911 \item[1]{\em one planar hydrogen, {\eg} rings or peptide bond}\\
912 One hydrogen atom (n) is generated, lying in the plane of atoms
913 (i,j,k) on the plane bisecting angle (j-i-k) at a distance of 0.1 nm
914 from atom i, such that the angles (n-i-j) and (n-i-k) are $>$ 90$^{\rm o}$.
915
916 \item[2]{\em one single hydrogen, {\eg} hydroxyl}\\
917 One hydrogen atom (n) is generated at a distance of 0.1 nm from atom
918 i, such that angle (n-i-j)=109.5 degrees and dihedral (n-i-j-k)=trans.
919
920 \item[3]{\em two planar hydrogens, {\eg} ethylene -C=CH{$_2$}, or amide -C(=O)NH{$_2$}}\\
921 Two hydrogens (n1,n2) are generated at a distance of 0.1 nm from atom
922 i, such that angle (n1-i-j)=(n2-i-j)=120 degrees and dihedral
923 (n1-i-j-k)=cis and (n2-i-j-k)=trans, such that names are according to
924 IUPAC standards~\cite{iupac70}.
925
926 \item[4]{\em two or three tetrahedral hydrogens, {\eg} -CH{$_3$}}\\
927 Three (n1,n2,n3) or two (n1,n2) hydrogens are generated at a distance
928 of 0.1 nm from atom i, such that angle
929 (n1-i-j)=(n2-i-j)=(n3-i-j)=109.47$^{\rm o}$, dihedral (n1-i-j-k)=trans,
930 (n2-i-j-k)=trans+120 and (n3-i-j-k)=trans+240$^{\rm o}$.
931
932 \item[5]{\em one tetrahedral hydrogen, {\eg} C{$_3$}CH}\\
933 One hydrogen atom (n$^\prime$) is generated at a distance of 0.1 nm from atom
934 i in tetrahedral conformation such that angle
935 (n$^\prime$-i-j)=(n$^\prime$-i-k)=(n$^\prime$-i-l)=109.47$^{\rm o}$.
936
937 \item[6]{\em two tetrahedral hydrogens, {\eg} C-CH{$_2$}-C}\\
938 Two hydrogen atoms (n1,n2) are generated at a distance of 0.1 nm from
939 atom i in tetrahedral conformation on the plane bisecting angle j-i-k
940 with angle (n1-i-n2)=(n1-i-j)=(n1-i-k)=109.47$^{\rm o}$.
941
942 \item[7]{\em two water hydrogens}\\
943 Two hydrogens are generated around atom i according to
944 SPC~\cite{Berendsen81} water geometry. The symmetry axis will
945 alternate between three coordinate axes in both directions.
946
947 \item[10]{\em three water ``hydrogens''}\\
948 Two hydrogens are generated around atom i according to
949 SPC~\cite{Berendsen81} water geometry. The symmetry axis will
950 alternate between three coordinate axes in both directions. In addition,
951 an extra particle is generated on the position of the oxygen with
952 the first letter of the name replaced by `M'. This is for
953 use with four-atom water models such as TIP4P~\cite{Jorgensen83}.
954
955 \item[11]{\em four water ``hydrogens''}\\
956 Same as above, except that two additional
957 particles are generated on the position of the oxygen, with names
958 `LP1' and `LP2.' This is for
959 use with five-atom water models such as TIP5P~\cite{Mahoney2000a}.
960 \end{enumerate}
961
962 \item
963 The name of the new H atom (or its prefix, {\eg} {\tt HD2} for
964 the asparagine example given earlier).
965
966 \item
967 Three or four control atoms (i,j,k,l), where the first always is the
968 atom to which the H atoms are connected. The other two or three depend
969 on the code selected. For water, there is only one control atom.
970 \end{itemize}
971
972 Some more exotic cases can be approximately constructed from the above tools,
973 and with suitable use of energy minimization are good enough for beginning
974 MD simulations. For example secondary amine hydrogen, nitrenyl hydrogen
975 (C\nolinebreak[4]=\nolinebreak[4]NH) and even ethynyl hydrogen could be
976 approximately constructed using method 2 above for hydroxyl hydrogen.
977
978 \subsection{Termini database}
979 \label{subsec:tdb}
980 The \swapindex{termini}{database}s are stored in {\tt aminoacids.n.tdb} and
981 {\tt aminoacids.c.tdb} for the N- and C-termini respectively. They contain
982 information for the {\tt pdb2gmx} program on how to connect new atoms
983 to existing ones, which atoms should be removed or changed, and which
984 bonded interactions should be added. Their format is as follows
985 (from {\tt gromos43a1.ff/aminoacids.c.tdb}):
986
987 {\small
988 \begin{verbatim}
989 [ None ]
990
991 [ COO- ]
992 [ replace ]
993 C       C       C       12.011  0.27
994 O       O1      OM      15.9994 -0.635
995 OXT     O2      OM      15.9994 -0.635
996 [ add ]
997 2       8       O       C       CA      N
998         OM      15.9994 -0.635
999 [ bonds ]
1000 C       O1      gb_5
1001 C       O2      gb_5
1002 [ angles ]
1003 O1      C       O2      ga_37
1004 CA      C       O1      ga_21
1005 CA      C       O2      ga_21
1006 [ dihedrals ]
1007 N       CA      C       O2      gd_20
1008 [ impropers ]
1009 C       CA      O2      O1      gi_1
1010 \end{verbatim}}
1011
1012 The file is organized in blocks, each with a header specifying the
1013 name of the block. These blocks correspond to different types of
1014 termini that can be added to a molecule. In this example {\tt [~COO-~]}
1015 is the first block, corresponding to changing the terminal carbon
1016 atom into a deprotonated carboxyl group. {\tt [~None~]} is the
1017 second terminus type, corresponding to a terminus that leaves
1018 the molecule as it is. Block names cannot be any of the following:
1019 {\tt replace}, {\tt add}, {\tt delete}, {\tt bonds}, {\tt angles},
1020 {\tt dihedrals}, {\tt impropers}.  Doing so would interfere with
1021 the parameters of the block, and would probably also be very confusing
1022 to human readers.
1023
1024 For each block the following options are present:
1025 \begin{itemize}
1026 \item {\tt [~replace~]} \\
1027 Replace an existing atom by one with a different atom type, atom name,
1028 charge, and/or mass. This entry can be used to replace an atom that is
1029 present both in the input coordinates and in the {\tt .rtp} database,
1030 but also to only rename an atom in the input coordinates such that
1031 it matches the name in the force field. In the latter case, there
1032 should also be a corresponding {\tt [~add~]} section present that
1033 gives instructions to add the same atom, such that the position in the sequence
1034 and the bonding is known. Such an atom can be present in the input
1035 coordinates and kept, or not present and constructed by {\tt pdb2gmx}.
1036 For each atom to be replaced on line should be
1037 entered with the following fields:
1038 \begin{itemize}
1039 \item name of the atom to be replaced
1040 \item new atom name (optional)
1041 \item new atom type
1042 \item new mass
1043 \item new charge
1044 \end{itemize}
1045 \item {\tt [~add~]} \\
1046 Add new atoms. For each (group of) added atom(s), a two-line entry is
1047 necessary. The first line contains the same fields as an entry in the
1048 hydrogen database (name of the new atom, 
1049 number of atoms, type of addition, control atoms,
1050 see~\ssecref{hdb}), but the possible types of addition are extended
1051 by two more, specifically for C-terminal additions:
1052 \begin{enumerate}
1053 \item[8]{\em two carboxyl oxygens, -COO{$^-$}}\\
1054 Two oxygens (n1,n2) are generated according to rule 3, at a distance
1055 of 0.136 nm from atom i and an angle (n1-i-j)=(n2-i-j)=117 degrees
1056 \item[9]{\em carboxyl oxygens and hydrogen, -COOH}\\
1057 Two oxygens (n1,n2) are generated according to rule 3, at distances of
1058 0.123 nm and 0.125 nm from atom i for n1 and n2, respectively, and angles
1059 (n1-i-j)=121 and (n2-i-j)=115 degrees. One hydrogen (n$^\prime$) is generated
1060 around n2 according to rule 2, where n-i-j and n-i-j-k should be read
1061 as n$^\prime$-n2-i and n$^\prime$-n2-i-j, respectively.
1062 \end{enumerate}
1063 After this line, another line follows that specifies the details of
1064 the added atom(s), in the same way as for replacing atoms, {\ie}: 
1065 \begin{itemize}
1066 \item atom type
1067 \item mass
1068 \item charge
1069 \item charge group (optional)
1070 \end{itemize}
1071 Like in the hydrogen database (see~\ssecref{rtp}), when more than
1072 one atom is connected to an existing one, a number will be appended to
1073 the end of the atom name. {\bf Note} that, like in the hydrogen database, the
1074 atom name is now on the same line as the control atoms, whereas it was
1075 at the beginning of the second line prior to {\gromacs} version 3.3.
1076 When the charge group field is left out, the added atom will have
1077 the same charge group number as the atom that it is bonded to.
1078 \item {\tt [~delete~]}\\
1079 Delete existing atoms. One atom name per line.
1080 \item {\tt [~bonds~]}, {\tt [~angles~]}, {\tt [~dihedrals~]} and {\tt [~impropers~]}\\
1081 Add additional bonded parameters. The format is identical to that used
1082 in the {\tt *.rtp} file, see~\ssecref{rtp}.
1083 \end{itemize}
1084
1085 \subsection{Virtual site database}
1086 Since we cannot rely on the positions of hydrogens in input files, we need a special
1087 input file to decide the geometries and parameters with which to add virtual site
1088 hydrogens. For more complex virtual site constructs ({\eg} when entire aromatic side chains
1089 are made rigid) we also need information about the equilibrium bond lengths and angles
1090 for all atoms in the side chain. This information is specified in the {\tt .vsd} file for each force 
1091 field. Just as for the termini, there is one such file for each class of residues in 
1092 the {\tt .rtp} file. 
1093
1094 The virtual site database is not really a very simple list of information. The first couple of sections
1095 specify which mass centers (typically called MCH$_3$/MNH$_3$) to use for CH$_3$, NH$_3$, 
1096 and NH$_2$ groups. Depending on the 
1097 equilibrium bond lengths and angles between the hydrogens and heavy atoms we need to apply
1098 slightly different constraint distances between these mass centers. {\bf Note} that we do {\em not} have to
1099 specify the actual parameters (that is automatic), just the type of mass center to use. To accomplish this,
1100 there are three sections names \verb+[ CH3 ]+, \verb+[ NH3 ]+, and \verb+[ NH2 ]+. For each of these we
1101 expect three columns. The first column is the atom type bound to the 2/3 hydrogens, the second column
1102 is the next heavy atom type which this is bound, and the third column the type of mass center to use.
1103 As a special case, in the  \verb+[ NH2 ]+ section it is also possible to specify \verb+planar+ in the second
1104 column, which will use a different construction without mass center. There are currently different opinions
1105 in some force fields whether an NH$_2$ group should be planar or not, but we try hard to stick to the
1106 default equilibrium parameters of the force field.
1107
1108 The second part of the virtual site database contains explicit equilibrium bond lengths and angles
1109 for pairs/triplets of atoms in aromatic side chains. These entries are currently read by specific routines
1110 in the virtual site generation code, so if you would like to extend it {\eg} to nucleic acids you would also
1111 need to write new code there. These sections are named after the short amino acid names
1112 (\verb+[ PHE ]+, \verb+[ TYR ]+, \verb+[ TRP ]+, \verb+[ HID ]+, \verb+[ HIE ]+, \verb+[ HIP ]+), and simply
1113 contain 2 or 3 columns with atom names, followed by a number specifying the bond length (in nm) or angle
1114 (in degrees). {\bf Note} that these are approximations of the equilibrated geometry for the entire molecule, 
1115 which might not be identical to the equilibrium value for a single bond/angle if the molecule is strained.
1116
1117 \subsection{Special bonds}
1118 \label{subsec:specbond}
1119 The primary mechanism used by {\tt \normindex{pdb2gmx}} to generate
1120 inter-residue bonds relies on head-to-tail linking of backbone atoms
1121 in different residues to build a macromolecule. In some cases ({\eg}
1122 \normindex{disulfide bonds}, a \normindex{heme group},
1123 \normindex{branched polymers}), it is necessary to create
1124 inter-residue bonds that do not lie on the backbone. The file {\tt
1125   \normindex{specbond.dat}} takes care of this function. It is
1126 necessary that the residues belong to the same {\tt [~moleculetype~]}.
1127 The {\tt -merge} and {\tt -chainsep} functions of {\tt pdb2gmx} can be
1128 useful when managing special inter-residue bonds between different
1129 chains.
1130
1131 The first line of {\tt specbond.dat} indicates the number of entries that are in the file. If you
1132 add a new entry, be sure to increment this number. The remaining lines in the file provide the
1133 specifications for creating bonds. The format of the lines is as follows:
1134
1135 {\tt resA  atomA  nbondsA  resB  atomB  nbondsB  length  newresA  newresB }
1136
1137 The columns indicate:
1138 \begin{enumerate}
1139 \item {\tt resA} The name of residue A that participates in the bond.
1140 \item {\tt atomA} The name of the atom in residue A that forms the bond.
1141 \item {\tt nbondsA} The total number of bonds {\tt atomA} can form.
1142 \item {\tt resB} The name of residue B that participates in the bond.
1143 \item {\tt atomB} The name of the atom in residue B that forms the bond.
1144 \item {\tt nbondsB} The total number of bonds {\tt atomB} can form.
1145 \item {\tt length} The reference length for the bond. If {\tt atomA} and {\tt atomB} are not within
1146 {\tt length} $\pm$ 10\% in the coordinate file supplied to {\tt pdb2gmx}, no bond will be formed.
1147 \item {\tt newresA} The new name of residue A, if necessary. Some force fields use {\eg} CYS2 for 
1148 a cysteine in a disulfide or heme linkage.
1149 \item {\tt newresB} The new name of residue B, likewise.
1150 \end{enumerate}
1151
1152
1153 \section{File formats}
1154 \subsection{Topology file\swapindexquiet{topology}{file}}
1155 \label{subsec:topfile}
1156 The topology file is built following the {\gromacs} specification for a
1157 molecular topology.  A {\tt *.top} file can be generated by
1158 {\tt pdb2gmx}.
1159 All possible entries in the topology file are listed in
1160 Tables \ref{tab:topfile1} and \ref{tab:topfile2}.
1161 Also tabulated are: all the units
1162 of the parameters, which interactions can be perturbed for free energy
1163 calculations, which bonded interactions are used by {\tt grompp}
1164 for generating exclusions, and which bonded interactions can be converted
1165 to constraints by {\tt grompp}.
1166
1167 %\renewcommand\floatpagefraction{.2}
1168
1169 \newcommand{\tts}{\tt \small}
1170
1171 % move these figures so they end up on facing pages 
1172 % (first figure on even page)
1173 \newcommand{\kJmol}{kJ~mol$^{-1}$}
1174 \newcommand{\kJmolnm}[1]{\kJmol~nm$^{#1}$}
1175 \newcommand{\kJmolrad}[1]{\kJmol~rad$^{#1}$}
1176 \newcommand{\kJmoldeg}[1]{\kJmol~deg$^{#1}$}
1177
1178 \begin{table}[p]
1179 \centering{
1180 \begin{tabular}{|l|llllc|}
1181 \multicolumn{6}{c}{\bf \large Parameters} \\
1182 \dline
1183 interaction     & directive           & \#  & f. & parameters                           & F. E. \\
1184 type    &                             & at. & tp &                                      &       \\
1185 \dline
1186 {\em mandatory} & {\tts defaults}       & & &   non-bonded function type; & \\
1187                 &                       & & &   combination rule$^{(cr)}$; &\\
1188                 &                       & & &   generate pairs (no/yes); & \\
1189                 &                       & & &   fudge LJ (); fudge QQ () & \\
1190 \hline
1191 {\em mandatory} & {\tts atomtypes}      &   &   & atom type; m (u); q (e); particle type; & \\
1192                 &                       &   &   & V$^{(cr)}$; W$^{(cr)}$ & \\
1193 %\hline
1194                 & {\tts bondtypes}      & \multicolumn{3}{l}{(see \tabref{topfile2}, directive {\tts bonds})}           & \\
1195                 & {\tts pairtypes}      & \multicolumn{3}{l}{(see \tabref{topfile2}, directive {\tts pairs})}           & \\
1196                 & {\tts angletypes}     & \multicolumn{3}{l}{(see \tabref{topfile2}, directive {\tts angles})}          & \\
1197                 & {\tts dihedraltypes}$^{(*)}$ & \multicolumn{3}{l}{(see \tabref{topfile2}, directive {\tts dihedrals})}& \\
1198                 & {\tts constrainttypes}& \multicolumn{3}{l}{(see \tabref{topfile2}, directive {\tts constraints})}     & \\
1199 LJ              & {\tts nonbond_params} & 2 & 1 & $V^{(cr)}$; $W^{(cr)}$ & \\
1200 Buckingham      & {\tts nonbond_params} & 2 & 2 & $a$ (\kJmol); $b$ (nm$^{-1})$;  & \\
1201  & & & & $c_6$ (\kJmolnm{6}) & \\
1202 \dline
1203 \multicolumn{6}{c}{~} \\
1204 \multicolumn{6}{c}{\bf \large Molecule definition(s)} \\
1205 \dline
1206 {\em mandatory} & {\tts moleculetype}   & & &   molecule name; $n_{ex}^{(nrexcl)}$  &   \\
1207 \hline
1208 {\em mandatory} & {\tts atoms}          & 1 &   & atom type; residue number;    & type  \\
1209                 &                       &   &   & residue name; atom name;      &       \\
1210                 &                       &   &   & charge group number; $q$ (e); $m$ (u)         & $q,m$ \\
1211 \hline
1212 \multicolumn{6}{|c|}{} \\
1213 \multicolumn{6}{|c|}{intra-molecular interaction and geometry definitions as described
1214 in \tabref{topfile2}} \\
1215 \multicolumn{6}{|c|}{} \\
1216 \dline
1217 \multicolumn{6}{c}{~} \\
1218 \multicolumn{6}{c}{\bf \large System} \\
1219 \dline
1220 {\em mandatory} & {\tts system}         & & &   system name     &       \\
1221 \hline
1222 {\em mandatory} & {\tts molecules}      & & &   \multicolumn{2}{l|}{molecule name; number of molecules} \\
1223 \dline
1224 \multicolumn{6}{c}{~} \\
1225 \multicolumn{6}{c}{\bf \large Inter-molecular interactions} \\
1226 \dline
1227 {\em optional}  & \multicolumn{4}{l}{\tts intermolecular_interactions}  & \\
1228 \hline
1229 \multicolumn{6}{|c|}{one or more bonded interactions as described in \tabref{topfile2}, with two or more atoms,} \\
1230 \multicolumn{6}{|c|}{no interactions that generate exclusions, no constraints, use global atom numbers} \\
1231 \dline
1232 \multicolumn{6}{c}{~} \\
1233 \multicolumn{6}{l}{`\# at' is the required number of atom type indices for this directive} \\
1234 \multicolumn{6}{l}{`f. tp' is the value used to select this function type} \\
1235 \multicolumn{6}{l}{`F. E.' indicates which of the parameters can be interpolated in free energy calculations} \\
1236 \multicolumn{6}{l}{~$^{(cr)}$ the combination rule determines the type of LJ parameters, see~\ssecref{nbpar}}\\
1237 \multicolumn{6}{l}{~$^{(*)}$ for {\tts dihedraltypes} one can specify 4 atoms or the inner (outer for improper) 2 atoms}\\
1238 \multicolumn{6}{l}{~$^{(nrexcl)}$ exclude neighbors $n_{ex}$ bonds away for non-bonded interactions}\\
1239 \multicolumn{6}{l}{For free energy calculations, type, $q$ and $m$  or no parameters should be added}\\
1240 \multicolumn{6}{l}{for topology `B' ($\lambda = 1$) on the same line, after the normal parameters.}
1241 \end{tabular}
1242 }
1243 \caption{The topology ({\tts *.top}) file.}
1244 \label{tab:topfile1}
1245 \end{table}
1246
1247 \newcommand{\fnm}[1]{\footnotemark[#1]}
1248 \renewcommand{\thefootnote}{\fnsymbol{footnote}}
1249 %\renewcommand{\tts}{\tt \small}
1250 \newcommand{\ttss}{\tt \scriptsize}
1251
1252 \begin{landscape}
1253 \begin{longtable}{|p{1.8in}|lcc>{\raggedright}p{2.5in}cc|}
1254 \caption{Details of {\tt [~moleculetype~]} directives}\\
1255 \dline
1256 Name of interaction              & Topology file directive          & num.  & func. & Order of parameters and their units                   & use in     & Cross- \\
1257                                  &                                  & atoms\fnm{1} & type\fnm{2} &                                          & F.E.?\fnm{3} & references \\
1258 \dline
1259 \endhead
1260 \dline
1261 \endfoot
1262 % The footnotetext fields only work inside the body, and not from a
1263 % column with ``p'' formatting'!
1264 bond                               & {\tts bonds}\fnm{4},\fnm{5}    & 2     & 1     & $b_0$ (nm); $k_b$ (\kJmolnm{-2})                      & all        & \ssecref{harmonicbond}
1265 \label{tab:topfile2} \footnotetext[1]{The required number of atom indices for this directive}\footnotetext[2]{The index to use to select this function type}\footnotetext[3]{Indicates which of the parameters can be interpolated in free energy calculations}\footnotetext[4]{This interaction type will be used by {{\tts grompp}} for generating exclusions}\footnotetext[5]{This interaction type can be converted to constraints by {{\tts grompp}}}\footnotetext[7]{The combination rule determines the type of LJ parameters, see~\ssecref{nbpar}}\footnotetext[6]{No connection, and so no exclusions, are generated for this interaction}
1266 \\
1267 G96 bond                           & {\tts bonds}\fnm{4},\fnm{5}    & 2     & 2     & $b_0$ (nm); $k_b$ (\kJmolnm{-4})                      & all        & \ssecref{G96bond} \\
1268 Morse                              & {\tts bonds}\fnm{4},\fnm{5}    & 2     & 3     & $b_0$ (nm); $D$ (\kJmol); $\beta$ (nm$^{-1}$)         & all        & \ssecref{Morsebond} \\
1269 cubic bond                         & {\tts bonds}\fnm{4},\fnm{5}    & 2     & 4     & $b_0$ (nm); $C_{i=2,3}$ (\kJmolnm{-i})                &            & \ssecref{cubicbond} \\
1270 connection                         & {\tts bonds}\fnm{4}            & 2     & 5     &                                                       &            & \tsecref{excl} \\
1271 harmonic potential                 & {\tts bonds}                   & 2     & 6     & $b_0$ (nm); $k_b$ (\kJmolnm{-2})                      & all        & \ssecref{harmonicbond},\tsecref{excl} \\
1272 FENE bond                          & {\tts bonds}\fnm{4}            & 2     & 7     & $b_m$ (nm); $k_b$ (\kJmolnm{-2})                      &            & \ssecref{FENEbond} \\
1273 tabulated bond                     & {\tts bonds}\fnm{4}            & 2     & 8     & table number ($\geq 0$); $k$ (\kJmol)                 & $k$        & \ssecref{tabulatedinteraction} \\
1274 tabulated bond\fnm{6}              & {\tts bonds}                   & 2     & 9     & table number ($\geq 0$); $k$ (\kJmol)                 & $k$        & \ssecref{tabulatedinteraction},\tsecref{excl} \\
1275 restraint potential                & {\tts bonds}                   & 2     & 10    & low, up$_1$, up$_2$ (nm); $k_{dr}$ (\kJmolnm{-2})     & all        & \ssecref{harmonicrestraint} \\
1276 extra LJ or Coulomb                & {\tts pairs}                   & 2     & 1     & $V$\fnm{7}; $W$\fnm{7}                                & all        & \ssecref{pairinteractions} \\
1277 extra LJ or Coulomb                & {\tts pairs}                   & 2     & 2     & fudge QQ (); $q_i$, $q_j$ (e), $V$\fnm{7}; $W$\fnm{7} &            & \ssecref{pairinteractions} \\
1278 extra LJ or Coulomb                & {\tts pairs_nb}                & 2     & 1     & $q_i$, $q_j$ (e); $V$\fnm{7}; $W$\fnm{7}              &            & \ssecref{pairinteractions} \\
1279 angle                              & {\tts angles}\fnm{5}           & 3     & 1     & $\theta_0$ (deg); $k_\theta$ (\kJmolrad{-2})          & all        & \ssecref{harmonicangle} \\
1280 G96 angle                          & {\tts angles}\fnm{5}           & 3     & 2     & $\theta_0$ (deg); $k_\theta$ (\kJmol)                 & all        & \ssecref{G96angle} \\
1281 cross bond-bond                    & {\tts angles}                  & 3     & 3     & $r_{1e}$, $r_{2e}$ (nm); $k_{rr'}$ (\kJmolnm{-2})     &            & \ssecref{bondbondcross} \\
1282 cross bond-angle                   & {\tts angles}                  & 3     & 4     & $r_{1e}$, $r_{2e}$ $r_{3e}$ (nm); $k_{r\theta}$ (\kJmolnm{-2}) &   & \ssecref{bondanglecross} \\
1283 Urey-Bradley                       & {\tts angles}\fnm{5}           & 3     & 5     & $\theta_0$ (deg); $k_\theta$ (\kJmolrad{-2}); $r_{13}$ (nm); $k_{UB}$ (\kJmolnm{-2}) & all & \ssecref{Urey-Bradley} \\
1284 quartic angle                      & {\tts angles}\fnm{5}           & 3     & 6     & $\theta_0$ (deg); $C_{i=0,1,2,3,4}$ (\kJmolrad{-i})   &            & \ssecref{quarticangle} \\
1285 tabulated angle                    & {\tts angles}                  & 3     & 8     & table number ($\geq 0$); $k$ (\kJmol)                 & $k$        & \ssecref{tabulatedinteraction} \\
1286 restricted bending potential       & {\tts angles}                  & 3     & 10    & $\theta_0$ (deg); $k_\theta$ (\kJmol)                 &            & \ssecref{ReB} \\
1287 proper dihedral                    & {\tts dihedrals}               & 4     & 1     & $\phi_s$ (deg); $k_\phi$ (\kJmol); multiplicity       & $\phi,k$   & \ssecref{properdihedral} \\
1288 improper dihedral                  & {\tts dihedrals}               & 4     & 2     & $\xi_0$ (deg); $k_\xi$ (\kJmolrad{-2})                & all        & \ssecref{harmonicimproperdihedral} \\
1289 Ryckaert-Bellemans dihedral        & {\tts dihedrals}               & 4     & 3     & $C_0$, $C_1$, $C_2$, $C_3$, $C_4$, $C_5$ (\kJmol)     & all        & \ssecref{RBdihedral} \\
1290 periodic improper dihedral         & {\tts dihedrals}               & 4     & 4     & $\phi_s$ (deg); $k_\phi$ (\kJmol); multiplicity       & $\phi,k$   & \ssecref{periodicimproperdihedral} \\
1291 Fourier dihedral                   & {\tts dihedrals}               & 4     & 5     & $C_1$, $C_2$, $C_3$, $C_4$ (\kJmol)                   & all        & \ssecref{Fourierdihedral} \\
1292 tabulated dihedral                 & {\tts dihedrals}               & 4     & 8     & table number ($\geq 0$); $k$ (\kJmol)                 & $k$        & \ssecref{tabulatedinteraction} \\
1293 proper dihedral (multiple)         & {\tts dihedrals}               & 4     & 9     & $\phi_s$ (deg); $k_\phi$ (\kJmol); multiplicity       & $\phi,k$   & \ssecref{properdihedral} \\
1294 restricted dihedral                & {\tts dihedrals}               & 4     & 10    & $\phi_0$ (deg); $k_\phi$ (\kJmol)                    &            & \ssecref{ReT} \\
1295 combined bending-torsion potential & {\tts dihedrals}               & 4     & 11    & $a_0$, $a_1$, $a_2$, $a_3$, $a_4$ (\kJmol)       &            & \ssecref{CBT} \\
1296 exclusions                         & {\tts exclusions}              & 1     &       & one or more atom indices                              &            & \tsecref{excl} \\
1297 constraint                         & {\tts constraints}\fnm{4}      & 2     & 1     & $b_0$ (nm)                                            & all        & \sssecref{constraints},\tsecref{constraints} \\
1298 constraint\fnm{6}                  & {\tts constraints}             & 2     & 2     & $b_0$ (nm)                                            & all        & \sssecref{constraints},\tsecref{constraints},\tsecref{excl} \\
1299 SETTLE                             & {\tts settles}                 & 1     & 1     & $d_{\mbox{\sc oh}}$, $d_{\mbox{\sc hh}}$ (nm)         &            & \ssecref{SETTLE},\tsecref{constraints} \\
1300 2-body virtual site                & {\tts virtual_sites2}          & 3     & 1     & $a$ ()                                                &            & \ssecref{vsite2} \\
1301 3-body virtual site                & {\tts virtual_sites3}          & 4     & 1     & $a$, $b$ ()                                           &            & \ssecref{vsite3} \\
1302 3-body virtual site (fd)           & {\tts virtual_sites3}          & 4     & 2     & $a$ (); $d$ (nm)                                      &            & \ssecref{vsite3fd} \\
1303 3-body virtual site (fad)          & {\tts virtual_sites3}          & 4     & 3     & $\theta$ (deg); $d$ (nm)                              &            & \ssecref{vsite3fad} \\
1304 3-body virtual site (out)          & {\tts virtual_sites3}          & 4     & 4     & $a$, $b$ (); $c$ (nm$^{-1}$)                          &            & \ssecref{vsite3out} \\
1305 4-body virtual site (fdn)          & {\tts virtual_sites4}          & 5     & 2     & $a$, $b$ (); $c$ (nm)                                 &            & \ssecref{vsite4fdn} \\
1306 N-body virtual site (COG)          & {\tts virtual_sitesn}          & 1     & 1     & one or more constructing atom indices                 &            & \ssecref{vsiteN} \\
1307 N-body virtual site (COM)          & {\tts virtual_sitesn}          & 1     & 2     & one or more constructing atom indices                 &            & \ssecref{vsiteN} \\
1308 N-body virtual site (COW)          & {\tts virtual_sitesn}          & 1     & 3     & one or more pairs consisting of constructing atom index and weight & & \ssecref{vsiteN} \\
1309 position restraint                 & {\ttss position_restraints}    & 1     & 1     & $k_{x}$, $k_{y}$, $k_{z}$ (\kJmolnm{-2}) & all                     & \ssecref{positionrestraint} \\
1310 flat-bottomed position restraint   & {\ttss position_restraints}    & 1     & 2     & $g$, $r$ (nm), $k$ (\kJmolnm{-2})                     &            & \ssecref{fbpositionrestraint} \\
1311 %restraint potential                & {\tts bonds}                   & 2     & 10    & low, up$_1$, up$_2$ (nm); $k_{dr}$ (\kJmolnm{-2})     &            & \ssecref{} \\
1312 distance restraint                 & {\ttss distance_restraints}    & 2     & 1     & type; label; low, up$_1$, up$_2$ (nm); weight ()      &            & \ssecref{distancerestraint} \\
1313 dihedral restraint                 & {\ttss dihedral_restraints}    & 4     & 1     & $\phi_0$ (deg); $\Delta\phi$ (deg); $k_\mathrm{dihr}$ (\kJmolrad{-2}) &    all   & \ssecref{dihedralrestraint} \\
1314 orientation restraint              & {\ttss orientation_restraints} & 2     & 1     & exp.; label; $\alpha$; $c$ (U nm$^\alpha$); obs. (U); weight (U$^{-1}$) & & \ssecref{orientationrestraint} \\
1315 angle restraint                    & {\ttss angle_restraints}       & 4     & 1     & $\theta_0$ (deg); $k_c$ (\kJmol); multiplicity        & $\theta,k$ & \ssecref{anglerestraint} \\
1316 angle restraint (z)                & {\ttss angle_restraints_z}     & 2     & 1     & $\theta_0$ (deg); $k_c$ (\kJmol); multiplicity        & $\theta,k$ & \ssecref{anglerestraint} \\
1317 \end{longtable}
1318 \end{landscape}
1319
1320 \renewcommand{\thefootnote}{\arabic{footnote}}
1321
1322 %\renewcommand\floatpagefraction{.5}
1323
1324
1325 Description of the file layout:
1326 \begin{itemize}
1327 \item Semicolon (;) and newline characters surround comments
1328 \item On a line ending with $\backslash$ the newline character is ignored.
1329 \item Directives are surrounded by {\tt [} and {\tt ]}
1330 \item The topology hierarchy (which must be followed) consists of three levels:
1331 \begin{itemize}
1332 \item the parameter level, which defines certain force-field specifications 
1333       (see~\tabref{topfile1})
1334 \item the molecule level, which should contain one or more molecule
1335       definitions (see~\tabref{topfile2})
1336 \item the system level, containing only system-specific information 
1337       ({\tt [~system~]} and {\tt [~molecules~]})
1338 \end{itemize}
1339 \item Items should be separated by spaces or tabs, not commas
1340 \item Atoms in molecules should be numbered consecutively starting at 1
1341 \item Atoms in the same charge group must be listed consecutively
1342 \item The file is parsed only once, which implies that no forward
1343       references can be treated: items must be defined before they
1344       can be used
1345 \item Exclusions can be generated from the bonds or
1346       overridden manually
1347 \item The bonded force types can be generated from the atom types or
1348       overridden per bond
1349 \item It is possible to apply multiple bonded interactions of the same type
1350       on the same atoms
1351 \item Descriptive comment lines and empty lines are highly recommended
1352 \item Starting with {\gromacs} version 3.1.3, all directives at the
1353       parameter level can be used multiple times and there are no
1354       restrictions on the order, except that an atom type needs to be
1355       defined before it can be used in other parameter definitions
1356 \item If parameters for a certain interaction are defined multiple times
1357       for the same combination of atom types the last definition is used;
1358       starting with {\gromacs} version 3.1.3 {\tt grompp} generates a
1359       warning for parameter redefinitions with different values
1360 \item Using one of the {\tt [~atoms~]}, {\tt [~bonds~]},
1361       {\tt [~pairs~]}, {\tt [~angles~]}, etc. without having used
1362       {\tt [~moleculetype~]}
1363       before is meaningless and generates a warning
1364 \item Using {\tt [~molecules~]} without having used
1365       {\tt [~system~]} before is meaningless and generates a warning.
1366 \item After {\tt [~system~]} the only allowed directive is {\tt [~molecules~]}
1367 \item Using an unknown string in {\tt [ ]} causes all the data until
1368       the next directive to be ignored and generates a warning
1369 \end{itemize}
1370
1371 Here is an example of a topology file, {\tt urea.top}:
1372
1373 {\small
1374 \begin{verbatim}
1375 ;
1376 ;       Example topology file
1377 ;
1378 ; The force-field files to be included
1379 #include "amber99.ff/forcefield.itp"
1380
1381 [ moleculetype ]
1382 ; name  nrexcl
1383 Urea         3
1384
1385 [ atoms ]
1386    1  C  1  URE      C      1     0.880229  12.01000   ; amber C  type
1387    2  O  1  URE      O      2    -0.613359  16.00000   ; amber O  type
1388    3  N  1  URE     N1      3    -0.923545  14.01000   ; amber N  type
1389    4  H  1  URE    H11      4     0.395055   1.00800   ; amber H  type
1390    5  H  1  URE    H12      5     0.395055   1.00800   ; amber H  type
1391    6  N  1  URE     N2      6    -0.923545  14.01000   ; amber N  type
1392    7  H  1  URE    H21      7     0.395055   1.00800   ; amber H  type
1393    8  H  1  URE    H22      8     0.395055   1.00800   ; amber H  type
1394
1395 [ bonds ]
1396     1   2
1397     1   3       
1398     1   6
1399     3   4
1400     3   5
1401     6   7
1402     6   8
1403
1404 [ dihedrals ] 
1405 ;   ai    aj    ak    al funct  definition
1406      2     1     3     4   9     
1407      2     1     3     5   9     
1408      2     1     6     7   9     
1409      2     1     6     8   9     
1410      3     1     6     7   9     
1411      3     1     6     8   9     
1412      6     1     3     4   9     
1413      6     1     3     5   9     
1414
1415 [ dihedrals ] 
1416      3     6     1     2   4     
1417      1     4     3     5   4     
1418      1     7     6     8   4
1419
1420 [ position_restraints ]
1421 ; you wouldn't normally use this for a molecule like Urea,
1422 ; but we include it here for didactic purposes
1423 ; ai   funct    fc
1424    1     1     1000    1000    1000 ; Restrain to a point
1425    2     1     1000       0    1000 ; Restrain to a line (Y-axis)
1426    3     1     1000       0       0 ; Restrain to a plane (Y-Z-plane)
1427
1428 [ dihedral_restraints ]
1429 ; ai   aj    ak    al  type  phi  dphi  fc
1430     3    6     1    2     1  180     0  10
1431     1    4     3    5     1  180     0  10
1432
1433 ; Include TIP3P water topology
1434 #include "amber99/tip3p.itp"
1435
1436 [ system ]
1437 Urea in Water
1438
1439 [ molecules ]
1440 ;molecule name   nr.
1441 Urea             1
1442 SOL              1000
1443 \end{verbatim}}
1444
1445 Here follows the explanatory text.
1446
1447 {\bf {\tt \#include "amber99.ff/forcefield.itp"} :} this includes the
1448 information for the force field you are using, including
1449 bonded and non-bonded parameters. This example uses the AMBER99 force
1450 field, but your simulation may use a different force field.
1451 {\tt grompp} will automatically go and find this file and copy-and-paste
1452 its content. That content can be seen in
1453 \linebreak {\tt share/top/amber99.ff/forcefield.itp}, and it is
1454
1455 {\small
1456 \begin{verbatim}
1457 #define _FF_AMBER
1458 #define _FF_AMBER99
1459
1460 [ defaults ]
1461 ; nbfunc        comb-rule       gen-pairs       fudgeLJ fudgeQQ
1462 1               2               yes             0.5     0.8333
1463
1464 #include "ffnonbonded.itp"
1465 #include "ffbonded.itp"
1466 #include "gbsa.itp"
1467 \end{verbatim}}
1468
1469 The two {\tt \#define} statements set up the conditions so that
1470 future parts of the topology can know that the AMBER 99 force
1471 field is in use.
1472
1473 {\bf {\tt [~defaults~]} :}
1474 \begin{itemize}
1475 \item {\tt nbfunc} is the non-bonded function type. Use 1 (Lennard-Jones) or 2 (Buckingham)
1476 \item {\tt comb-rule} is the number of the \normindex{combination rule} (see \ssecref{nbpar}).
1477 \item {\tt gen-pairs} is for pair generation. The default is `no', {\ie} 
1478 get 1-4 parameters from the pairtypes list. When parameters
1479 are not present in the list, stop with a fatal error.
1480 Setting `yes' generates 1-4 parameters that are not present in the pair list
1481 from normal Lennard-Jones parameters using {\tt fudgeLJ}
1482 \item {\tt fudgeLJ} is the factor by which to multiply Lennard-Jones 1-4 interactions, default 1
1483 \item {\tt fudgeQQ} is the factor by which to multiply electrostatic 1-4 interactions, default 1
1484 \item $N$ is the power for the repulsion term in a 6-$N$ potential (with 
1485 nonbonded-type Lennard-Jones only), starting with {\gromacs} version 4.5,
1486 {\tt mdrun} also reads and applies $N$, for values not equal to 12 tabulated
1487 interaction functions are used
1488 (in older version you would have to use user tabulated interactions).
1489 \end{itemize}
1490 {\bf Note} that {\tt gen-pairs}, {\tt fudgeLJ}, {\tt fudgeQQ}, and $N$ are optional.
1491 {\tt fudgeLJ} is only used when generate pairs is set to `yes', and
1492 {\tt fudgeQQ} is always used. However, if you
1493 want to specify $N$ you need to give a value for the other parameters as well.
1494
1495 Then some other {\tt \#include} statements add in the large amount of data needed
1496 to describe the rest of the force field. We will skip these and return to {\tt urea.top}.
1497 There we will see
1498
1499 % move these figures so they end up on facing pages 
1500 % (first figure on even page)
1501 %\input{topolfig}
1502
1503 {\bf {\tt [~moleculetype~]} :} defines the name of your molecule in
1504 this {\tt *.top} and nrexcl = 3 stands for excluding non-bonded
1505 interactions between atoms that are no further than 3 bonds away.
1506
1507 {\bf {\tt [~atoms~]} :} defines the molecule, where {\tt nr} and
1508 {\tt type} are fixed, the rest is user defined. So {\tt atom} can be named
1509 as you like, {\tt cgnr} made larger or smaller (if possible, the total
1510 charge of a charge group should be zero), and charges can be changed
1511 here too.
1512
1513 {\bf {\tt [~bonds~]} :} no comment.
1514
1515 {\bf {\tt [~pairs~]} :} LJ and Coulomb 1-4 interactions
1516
1517 {\bf {\tt [~angles~]} :} no comment
1518
1519 {\bf {\tt [~dihedrals~]} :} in this case there are 9 proper dihedrals
1520 (funct = 1), 3 improper (funct = 4) and no Ryckaert-Bellemans type
1521 dihedrals. If you want to include Ryckaert-Bellemans type dihedrals
1522 in a topology, do the following (in case of {\eg} decane):
1523 \begin{verbatim}
1524 [ dihedrals ]
1525 ;  ai    aj    ak    al funct       c0       c1       c2
1526     1    2     3     4     3 
1527     2    3     4     5     3
1528 \end{verbatim}
1529 In the original implementation of the potential for
1530 alkanes~\cite{Ryckaert78} no 1-4 interactions were used, which means
1531 that in order to implement that particular force field you need to remove the 1-4
1532 interactions from the {\tt [~pairs~]} section of your topology. In
1533 most modern force fields, like OPLS/AA or Amber the rules are
1534 different, and the Ryckaert-Bellemans potential is used as a cosine
1535 series in combination with 1-4 interactions.
1536
1537 {\bf {\tt [~position_restraints~]} :} harmonically restrain the selected particles
1538 to reference positions (\ssecref{positionrestraint}). 
1539 The reference positions are read from a 
1540 separate coordinate file by {\tt \normindex{grompp}}.
1541
1542
1543 {\bf {\tt [~dihedral_restraints~]} :} restrain selected dihedrals to a reference value.
1544 The implementation of dihedral restraints is described in section \ssecref{dihedralrestraint} of the manual.
1545 The parameters specified in the [dihedral_restraints] directive are as follows:
1546 \begin{itemize}
1547 \item {\tt type} has only one possible value which is 1
1548 \item {\tt phi} is the value of $\phi_0$ in \eqnref{dphi} and \eqnref{dihre} of the manual.
1549 \item {\tt dphi} is the value of $\Delta\phi$ in \eqnref{dihre} of the manual.
1550 \item {\tt fc} is the force constant $k_{dihr}$ in \eqnref{dihre} of the manual.
1551 \end{itemize}
1552
1553 {\bf {\tt \#include "tip3p.itp"} :} includes a topology file that was already
1554 constructed (see section~\ssecref{molitp}).
1555
1556 {\bf {\tt [~system~]} :} title of your system, user-defined
1557
1558 {\bf {\tt [~molecules~]} :} this defines the total number of (sub)molecules
1559 in your system that are defined in this {\tt *.top}. In this
1560 example file, it stands for 1 urea molecule dissolved in 1000 water
1561 molecules. The molecule type SOL is defined in the {\tt tip3p.itp} file.
1562 Each name here must correspond to a name given with {\tt [~moleculetype~]}
1563 earlier in the topology. The order of the blocks of molecule types and
1564 the numbers of such molecules must match the coordinate file that
1565 accompanies the topology when supplied to {\tt \normindex{grompp}}.
1566 The blocks of molecules do not need to be contiguous, but some
1567 tools (e.g. {\tt \normindex{genion}}) may act only on the first or
1568 last such block of a particular molecule type. Also, these blocks
1569 have nothing to do with the definition of \normindex{groups}
1570 (see \secref{groupconcept} and \secref{usinggroups}).
1571
1572 \subsection{Molecule.itp file}
1573 \label{subsec:molitp}
1574 If you construct a topology file you will use frequently (like the water
1575 molecule, {\tt tip3p.itp}, which is already constructed for you) it is
1576 good to make a {\tt molecule.itp} file. This only lists the
1577 information of one particular molecule and allows you to re-use the
1578 {\tt [ moleculetype ]} in multiple systems without re-invoking
1579 {\tt pdb2gmx} or manually copying and pasting. An example
1580 {\tt urea.itp} follows: 
1581
1582 {\small
1583 \begin{verbatim}
1584 [ moleculetype ]
1585 ; molname       nrexcl
1586 URE             3
1587
1588 [ atoms ]
1589    1  C  1  URE      C      1     0.880229  12.01000   ; amber C  type
1590 ...
1591    8  H  1  URE    H22      8     0.395055   1.00800   ; amber H  type
1592
1593 [ bonds ]
1594     1   2
1595 ...
1596     6   8
1597 [ dihedrals ] 
1598 ;   ai    aj    ak    al funct  definition
1599      2     1     3     4   9     
1600 ...
1601      6     1     3     5   9     
1602 [ dihedrals ] 
1603      3     6     1     2   4     
1604      1     4     3     5   4     
1605      1     7     6     8   4
1606 \end{verbatim}}
1607
1608 Using {\tt *.itp} files results in a very short {\tt *.top} file:
1609
1610 {\small
1611 \begin{verbatim}
1612 ;
1613 ;       Example topology file
1614 ;
1615 ; The force field files to be included
1616 #include "amber99.ff/forcefield.itp"
1617
1618 #include "urea.itp"
1619
1620 ; Include TIP3P water topology
1621 #include "amber99/tip3p.itp"
1622
1623 [ system ]
1624 Urea in Water
1625
1626 [ molecules ]
1627 ;molecule name   nr.
1628 Urea             1
1629 SOL              1000
1630 \end{verbatim}}
1631
1632 \subsection{Ifdef statements}
1633 \label{subsec:ifdef}
1634 A very powerful feature in {\gromacs} is the use of {\tt \#ifdef}
1635 statements in your {\tt *.top} file. By making use of this statement,
1636 and associated {\tt \#define} statements like were seen in
1637 \linebreak {\tt amber99.ff/forcefield.itp} earlier,
1638 different parameters for one molecule can be used in the same
1639 {\tt *.top} file. An example is given for TFE, where there is an option to
1640 use different charges on the atoms: charges derived by De Loof
1641 {\etal}~\cite{Loof92} or by Van Buuren and
1642 Berendsen~\cite{Buuren93a}. In fact, you can use much of the functionality of the
1643 C preprocessor, {\tt cpp}, because {\tt grompp} contains similar pre-processing
1644 functions to scan the file.  The
1645 way to make use of the {\tt \#ifdef} option is as follows:
1646 \begin{itemize}
1647 \item either use the option {\tt define = -DDeLoof} in the
1648       {\tt *.mdp} file (containing {\tt grompp} input
1649       parameters), or use the line {\tt \#define DeLoof}
1650       early in your {\tt *.top} or {\tt *.itp} file; and
1651 \item put the {\tt \#ifdef} statements in your {\tt *.top}, as
1652       shown below: 
1653 \end{itemize}
1654
1655 {\small
1656 \begin{verbatim}
1657 ...
1658
1659
1660
1661 [ atoms ]
1662 ; nr     type     resnr    residu     atom      cgnr      charge        mass
1663 #ifdef DeLoof
1664 ; Use Charges from DeLoof
1665    1        C        1        TFE        C         1        0.74        
1666    2        F        1        TFE        F         1       -0.25        
1667    3        F        1        TFE        F         1       -0.25        
1668    4        F        1        TFE        F         1       -0.25        
1669    5      CH2        1        TFE      CH2         1        0.25        
1670    6       OA        1        TFE       OA         1       -0.65        
1671    7       HO        1        TFE       HO         1        0.41        
1672 #else
1673 ; Use Charges from VanBuuren
1674    1        C        1        TFE        C         1        0.59        
1675    2        F        1        TFE        F         1       -0.2         
1676    3        F        1        TFE        F         1       -0.2         
1677    4        F        1        TFE        F         1       -0.2         
1678    5      CH2        1        TFE      CH2         1        0.26        
1679    6       OA        1        TFE       OA         1       -0.55        
1680    7       HO        1        TFE       HO         1        0.3         
1681 #endif
1682
1683 [ bonds ]
1684 ;  ai    aj funct           c0           c1
1685     6     7     1 1.000000e-01 3.138000e+05 
1686     1     2     1 1.360000e-01 4.184000e+05 
1687     1     3     1 1.360000e-01 4.184000e+05 
1688     1     4     1 1.360000e-01 4.184000e+05 
1689     1     5     1 1.530000e-01 3.347000e+05 
1690     5     6     1 1.430000e-01 3.347000e+05 
1691 ...
1692 \end{verbatim}}
1693
1694 This mechanism is used by {\tt pdb2gmx} to implement optional position
1695 restraints (\ssecref{positionrestraint}) by {\tt \#include}-ing an {\tt .itp} file whose contents
1696 will be meaningful only if a particular {\tt \#define} is set (and spelled
1697 correctly!)
1698
1699 \subsection{Topologies for free energy calculations}
1700 \index{free energy topologies}
1701 Free energy differences between two systems, A and B, can be calculated as
1702 described in \secref{fecalc}.
1703 Systems A and B are described by topologies
1704 consisting of the same number of molecules with the same number of
1705 atoms. Masses and non-bonded interactions can be perturbed by adding B
1706 parameters under the {\tt [~atoms~]} directive. Bonded interactions can be 
1707 perturbed by adding B parameters to the bonded types or the bonded
1708 interactions. The parameters that can be perturbed are listed in  
1709 Tables \ref{tab:topfile1} and \ref{tab:topfile2}.
1710 The $\lambda$-dependence of the interactions is described
1711 in section \secref{feia}.
1712 The bonded parameters that are used (on the line of the bonded
1713 interaction definition, or the ones looked up on atom types
1714 in the bonded type lists) is explained in \tabref{topfe}.
1715 In most cases, things should work intuitively.
1716 When the A and B atom types in a bonded interaction
1717 are not all identical and parameters are not present for the B-state,
1718 either on the line or in the bonded types,
1719 {\tt grompp} uses the A-state parameters and issues a warning.
1720 For free energy calculations, all or no parameters for topology B
1721 ($\lambda = 1$) should be added on the same line, after the normal
1722 parameters, in the same order as the normal parameters.
1723 From {\gromacs} 4.6 onward, if $\lambda$ is treated as a vector, then
1724 the {\tt bonded-lambdas} component controls all bonded terms that are
1725 not explicitly labeled as restraints.  Restrain terms are controlled
1726 by the {\tt restraint-lambdas} component.
1727
1728 \begin{table}
1729 \centerline{
1730 \begin{tabular}{|c|cc|cc|cc|c|}
1731 \dline
1732 B-state atom types & \multicolumn{2}{c|}{parameters} & \multicolumn{4}{c|}{parameters in bonded types} & \\
1733 all identical to      & \multicolumn{2}{c|}{on line} & \multicolumn{2}{c|}{A atom types} & \multicolumn{2}{c|}{B atom types} & message \\
1734 A-state atom types & A & B & A & B & A & B & \\
1735 \dline
1736     & +AB & $-$ &  x  &  x  &     &     & \\
1737     & +A  & +B  &  x  &  x  &     &     & \\
1738 yes & $-$ & $-$ & $-$ & $-$ &     &     & error \\
1739     & $-$ & $-$ & +AB & $-$ &     &     & \\
1740     & $-$ & $-$ & +A  & +B  &     &     & \\
1741 \hline
1742     & +AB & $-$ &  x  &  x  &  x  &  x  & warning \\
1743     & +A  & +B  &  x  &  x  &  x  &  x  & \\
1744     & $-$ & $-$ & $-$ & $-$ &  x  &  x  & error \\
1745 no  & $-$ & $-$ & +AB & $-$ & $-$ & $-$ & warning \\
1746     & $-$ & $-$ & +A  & +B  & $-$ & $-$ & warning \\
1747     & $-$ & $-$ & +A  &  x  & +B  & $-$ & \\
1748     & $-$ & $-$ & +A  &  x  &  +  & +B  & \\
1749 \dline
1750 \end{tabular}
1751 }
1752 \caption{The bonded parameters that are used for free energy topologies,
1753 on the line of the bonded interaction definition or looked up
1754 in the bond types section based on atom types. A and B indicate the
1755 parameters used for state A and B respectively, + and $-$ indicate
1756 the (non-)presence of parameters in the topology, x indicates that
1757 the presence has no influence.}
1758 \label{tab:topfe}
1759 \end{table}
1760
1761 Below is an example of a topology which changes from 200 propanols to
1762 200 pentanes using the \gromosv{96} force field.\\
1763
1764 {\small
1765 \begin{verbatim}
1766  
1767 ; Include force field parameters
1768 #include "gromos43a1.ff/forcefield.itp"
1769
1770 [ moleculetype ]
1771 ; Name            nrexcl
1772 PropPent          3
1773
1774 [ atoms ]
1775 ; nr type resnr residue atom cgnr  charge    mass  typeB chargeB  massB
1776   1    H    1     PROP    PH    1   0.398    1.008  CH3     0.0  15.035
1777   2   OA    1     PROP    PO    1  -0.548  15.9994  CH2     0.0  14.027
1778   3  CH2    1     PROP   PC1    1   0.150   14.027  CH2     0.0  14.027
1779   4  CH2    1     PROP   PC2    2   0.000   14.027
1780   5  CH3    1     PROP   PC3    2   0.000   15.035
1781
1782 [ bonds ]
1783 ;  ai    aj funct    par_A  par_B 
1784     1     2     2    gb_1   gb_26
1785     2     3     2    gb_17  gb_26
1786     3     4     2    gb_26  gb_26
1787     4     5     2    gb_26
1788
1789 [ pairs ]
1790 ;  ai    aj funct
1791     1     4     1
1792     2     5     1
1793
1794 [ angles ]
1795 ;  ai    aj    ak funct    par_A   par_B
1796     1     2     3     2    ga_11   ga_14
1797     2     3     4     2    ga_14   ga_14
1798     3     4     5     2    ga_14   ga_14
1799
1800 [ dihedrals ]
1801 ;  ai    aj    ak    al funct    par_A   par_B
1802     1     2     3     4     1    gd_12   gd_17
1803     2     3     4     5     1    gd_17   gd_17
1804
1805 [ system ]
1806 ; Name
1807 Propanol to Pentane
1808
1809 [ molecules ]
1810 ; Compound        #mols
1811 PropPent          200
1812 \end{verbatim}}
1813
1814 Atoms that are not perturbed, {\tt PC2} and {\tt PC3}, do not need B-state parameter
1815 specifications, since the B parameters will be copied from the A parameters.
1816 Bonded interactions between atoms that are not perturbed do not need B
1817 parameter specifications, as is the case for the last bond in the example topology.
1818 Topologies using the OPLS/AA force field need no bonded parameters at all,
1819 since both the A and B parameters are determined by the atom types.
1820 Non-bonded interactions involving one or two perturbed atoms use the 
1821 free-energy perturbation functional forms.
1822 Non-bonded interactions between two non-perturbed atoms use the normal
1823 functional forms.
1824 This means that when, for instance, only the charge of a particle is
1825 perturbed, its Lennard-Jones interactions will also be affected when
1826 lambda is not equal to zero or one.
1827
1828 {\bf Note} that this topology uses the \gromosv{96} force field, in which the bonded
1829 interactions are not determined by the atom types. The bonded interaction
1830 strings are converted by the C-preprocessor. The force-field parameter
1831 files contain lines like:
1832
1833 {\small
1834 \begin{verbatim}
1835 #define gb_26       0.1530  7.1500e+06
1836
1837 #define gd_17     0.000       5.86          3
1838 \end{verbatim}}
1839
1840 \subsection{Constraint forces\index{constraint force}}
1841 \label{subsec:constraintforce}
1842 The constraint force between two atoms in one molecule can be calculated
1843 with the free energy perturbation code by adding a constraint between the
1844 two atoms, with a different length in the A and B topology. When the B length
1845 is 1 nm longer than the A length and lambda is kept constant at zero,
1846 the derivative of the Hamiltonian with respect to lambda is the constraint
1847 force. For constraints between molecules, the pull code can be used,
1848 see \secref{pull}.
1849 Below is an example for calculating the constraint force at 0.7 nm
1850 between two methanes in water, by combining the two methanes into one ``molecule.''
1851 {\bf Note} that the definition of a ``molecule'' in {\gromacs} does not necessarily
1852 correspond to the chemical definition of a molecule.  In {\gromacs}, a ``molecule''
1853 can be defined as any group of atoms that one wishes to consider simultaneously.
1854 The added constraint is of function type 2, which means that it is not
1855 used for generating exclusions (see~\secref{excl}). 
1856 Note that the constraint free energy term is included in the derivative term, and is
1857 specifically included in the {\tt bonded-lambdas} component. However, the free
1858 energy for changing constraints is {\em not} included in the potential energy
1859 differences used for BAR and MBAR, as this requires reevaluating the energy at
1860 each of the constraint components.  This functionality is planned for later versions.\\
1861
1862 {\small
1863 \begin{verbatim}
1864 ; Include force-field parameters
1865 #include "gromos43a1.ff/forcefield.itp"
1866
1867 [ moleculetype ]
1868 ; Name            nrexcl
1869 Methanes               1
1870
1871 [ atoms ]
1872 ; nr   type   resnr  residu   atom    cgnr     charge    mass
1873    1    CH4     1     CH4      C1       1          0    16.043
1874    2    CH4     1     CH4      C2       2          0    16.043
1875 [ constraints ]
1876 ;  ai    aj funct   length_A  length_B
1877     1     2     2        0.7       1.7
1878
1879 #include "gromos43a1.ff/spc.itp"
1880
1881 [ system ]
1882 ; Name
1883 Methanes in Water
1884
1885 [ molecules ]
1886 ; Compound        #mols
1887 Methanes              1
1888 SOL                2002
1889 \end{verbatim}}
1890
1891 \subsection{Coordinate file}
1892 \label{subsec:grofile}
1893 Files with the {\tt .gro} file extension contain a molecular structure in 
1894 \gromosv{87} format. A sample piece is included below:
1895
1896 {\small
1897 \begin{verbatim}
1898 MD of 2 waters, reformat step, PA aug-91
1899     6
1900     1WATER  OW1    1   0.126   1.624   1.679  0.1227 -0.0580  0.0434
1901     1WATER  HW2    2   0.190   1.661   1.747  0.8085  0.3191 -0.7791
1902     1WATER  HW3    3   0.177   1.568   1.613 -0.9045 -2.6469  1.3180
1903     2WATER  OW1    4   1.275   0.053   0.622  0.2519  0.3140 -0.1734
1904     2WATER  HW2    5   1.337   0.002   0.680 -1.0641 -1.1349  0.0257
1905     2WATER  HW3    6   1.326   0.120   0.568  1.9427 -0.8216 -0.0244
1906    1.82060   1.82060   1.82060
1907 \end{verbatim}}
1908
1909 This format is fixed, {\ie} all columns are in a fixed position. If you
1910 want to read such a file in your own program without using the
1911 {\gromacs} libraries you can use the following formats:
1912
1913 {\bf C-format:} {\tt "\%5i\%5s\%5s\%5i\%8.3f\%8.3f\%8.3f\%8.4f\%8.4f\%8.4f"}
1914
1915 Or to be more precise, with title {\em etc.} it looks like this:
1916
1917 \begin{verbatim}
1918   "%s\n", Title
1919   "%5d\n", natoms
1920   for (i=0; (i<natoms); i++) {
1921     "%5d%-5s%5s%5d%8.3f%8.3f%8.3f%8.4f%8.4f%8.4f\n",
1922       residuenr,residuename,atomname,atomnr,x,y,z,vx,vy,vz
1923   }
1924   "%10.5f%10.5f%10.5f%10.5f%10.5f%10.5f%10.5f%10.5f%10.5f\n",
1925     box[X][X],box[Y][Y],box[Z][Z],
1926     box[X][Y],box[X][Z],box[Y][X],box[Y][Z],box[Z][X],box[Z][Y]
1927 \end{verbatim}
1928
1929 {\bf Fortran format:} {\tt (i5,2a5,i5,3f8.3,3f8.4)}
1930
1931 So {\tt confin.gro} is the {\gromacs} coordinate file and is almost
1932 the same as the \gromosv{87} file (for {\gromos} users: when used with
1933 {\tt ntx=7}).  The only difference is the box for which {\gromacs} uses a
1934 tensor, not a vector.
1935
1936
1937
1938 \section{Force field organization \index{force field organization}}
1939 \label{sec:fforganization}
1940
1941 \subsection{Force-field files}
1942 \label{subsec:fffiles}
1943 Many force fields are available by default.
1944 Force fields are detected by the presence of {\tt <name>.ff} directories
1945 in the {\tt \$GMXLIB/share/gromacs/top} sub-directory and/or the working directory.
1946 The information regarding the location of the force field files is printed
1947 by {\tt pdb2gmx} so you can easily keep track of which version of a force field
1948 is being called, in case you have made modifications in one location or another.
1949 The force fields included with {\gromacs} are:
1950
1951 {\small
1952 \begin{itemize}
1953  \item AMBER03 protein, nucleic AMBER94 (Duan et al., J. Comp. Chem. 24, 1999-2012, 2003)
1954  \item AMBER94 force field (Cornell et al., JACS 117, 5179-5197, 1995)
1955  \item AMBER96 protein, nucleic AMBER94 (Kollman et al., Acc. Chem. Res. 29, 461-469, 1996)
1956  \item AMBER99 protein, nucleic AMBER94 (Wang et al., J. Comp. Chem. 21, 1049-1074, 2000)
1957  \item AMBER99SB protein, nucleic AMBER94 (Hornak et al., Proteins 65, 712-725, 2006)
1958  \item AMBER99SB-ILDN protein, nucleic AMBER94 (Lindorff-Larsen et al., Proteins 78, 1950-58, 2010)
1959  \item AMBERGS force field (Garcia \& Sanbonmatsu, PNAS 99, 2782-2787, 2002)
1960  \item CHARMM27 all-atom force field (CHARM22 plus CMAP for proteins)
1961  \item GROMOS96 43a1 force field
1962  \item GROMOS96 43a2 force field (improved alkane dihedrals)
1963  \item GROMOS96 45a3 force field (Schuler JCC 2001 22 1205)
1964  \item GROMOS96 53a5 force field (JCC 2004 vol 25 pag 1656)
1965  \item GROMOS96 53a6 force field (JCC 2004 vol 25 pag 1656)
1966  \item GROMOS96 54a7 force field (Eur. Biophys. J. (2011), 40,, 843-856, DOI: 10.1007/s00249-011-0700-9)
1967  \item OPLS-AA/L all-atom force field (2001 aminoacid dihedrals)
1968 \end{itemize}} 
1969  
1970 A force field is included at the beginning of a topology file with an
1971 {\tt \#include} statement followed by {\tt <name>.ff/forcefield.itp}.
1972 This statement includes the force-field file,
1973 which, in turn, may include other force-field files. All the force fields
1974 are organized in the same way. An example of the
1975 {\tt amber99.ff/forcefield.itp} was shown in \ssecref{topfile}.
1976
1977 For each force field, there several files which are only used by {\tt pdb2gmx}.
1978 These are: residue databases ({\tt .rtp}, see~\ssecref{rtp})
1979 the hydrogen database ({\tt .hdb}, see~\ssecref{hdb}), two termini databases
1980 ({\tt .n.tdb} and {\tt .c.tdb}, see~\ssecref{tdb}) and
1981 the atom type database ({\tt .atp}, see~\ssecref{atomtype}), which contains only the masses.  Other optional
1982 files are described in~\secref{pdb2gmxfiles}.
1983
1984
1985 \subsection{Changing force-field parameters\index{force field}}
1986 If one wants to change the parameters of few bonded interactions in
1987 a molecule, this is most easily accomplished by typing the parameters
1988 behind the definition of the bonded interaction directly in the {\tt *.top} file 
1989 under the {\tt [~moleculetype~]} section (see \ssecref{topfile} for the format
1990 and units).
1991 If one wants to change the parameters for all instances of a certain
1992 interaction one can change them in the force-field file or add a
1993 new {\tt [~???types~]} section after including the force field.
1994 When parameters for a certain interaction are defined multiple times,
1995 the last definition is used. As of {\gromacs} version 3.1.3, a warning is
1996 generated when parameters are redefined with a different value.
1997 Changing the Lennard-Jones parameters of an atom type is not
1998 recommended, because in the {\gromos} force fields
1999 the Lennard-Jones parameters for several combinations of atom types
2000 are not generated according to the standard combination rules.
2001 Such combinations (and possibly others that do follow the
2002 combination rules) are defined in the {\tt [~nonbond_params~]}
2003 section, and changing the Lennard-Jones parameters of an atom type
2004 has no effect on these combinations.
2005
2006 \subsection{Adding atom types\swapindexquiet{adding}{atom types}}
2007 As of {\gromacs} version 3.1.3, atom types can be added in an extra
2008 {\tt [~atomtypes~]} section after the the inclusion of the normal
2009 force field. After the definition of the new atom type(s), additional
2010 non-bonded and pair parameters can be defined.
2011 In pre-3.1.3 versions of {\gromacs}, the new atom types needed to be
2012 added in the {\tt [~atomtypes~]} section of the force-field files,
2013 because all non-bonded parameters above the last {\tt [~atomtypes~]}
2014 section would be overwritten using the standard combination rules.
2015
2016 % LocalWords:  parameterized fffiles ptype polarizable gromacs atp ype arameter
2017 % LocalWords:  lll carboxyl OA hydroxyl NL porphyrin OPLS CP HCR OWT fd funct
2018 % LocalWords:  grompp statprop atomtype rtp esidue opology pdb gmx kJ mol gro
2019 % LocalWords:  grofile dihedrals bon itp func kb th cth cq cp mult Ryckaert aj
2020 % LocalWords:  Bellemans ak alkanes alkane llrllrllr LJ der nb topfile llllll
2021 % LocalWords:  llll nonbond params ij pairtypes fecalc moleculetype indices mdp
2022 % LocalWords:  constraintforce SPC molname nrexcl nr ren HW doh dhh aminoacids
2023 % LocalWords:  dat basename rna dna arn hdb sn rtpo gmxfiles molitp ndx ARG CYS
2024 % LocalWords:  defaultgroups impropers chargegroup bondedtypes hydrogens ARGN
2025 % LocalWords:  preprocessor protonation specbond protonated arginine aspartic
2026 % LocalWords:  ASPH GLU glutamic GLUH HISD histidine HISE HISH LYSN LYS IUPAC
2027 % LocalWords:  wildcards xlateat asparagine HD HH cis deprotonated oxygens COOH
2028 % LocalWords:  llllc tp cr QQ atomtypes bondtypes angletypes dihedraltypes FENE
2029 % LocalWords:  constrainttypes intra nbpar morse dr Coul rr UB dih constr hh ai
2030 % LocalWords:  vsite sitesn construc restr ffgmx resnr residu cgnr al fc spc gb
2031 % LocalWords:  FudgeLJ FudgeQQ nonbonded mdrun decane posre Ifdef ifdef TFE cpp
2032 % LocalWords:  Loof Buuren Berendsen DDeloof DeLoof VanBuuren endif feia topfe
2033 % LocalWords:  propanols pentanes ffG PropPent typeB chargeB massB ga gd mols
2034 % LocalWords:  Propanol Pentane methanes aug natoms residuenr residuename vx vy
2035 % LocalWords:  atomname atomnr vz Fortran confin ntx GROMOS nbfunc GROningen ff
2036 % LocalWords:  fudgeLJ fudgeQQ ffgmxnb ffgmxbon tdb ffbonded ffnonbonded nonbond
2037 % LocalWords:  MAchine BIOSON Groningen Spoel Drunen Comp Phys Comm trr AA fdn
2038 % LocalWords:  aliphatic CHARMM polarisability quadrupole tt normvsbds Waals jj
2039 % LocalWords:  pairinteractions num Buckingham rcl trans Intramolecular Lennard
2040 % LocalWords:  excl gen solute unscaled moltype intramol dgimplement Qiu HCT rt
2041 % LocalWords:  Onufriev OBC LINCS doc xxx residuetypes polyatomic co rotatable
2042 % LocalWords:  heme cysteine lysine CH NH LP amine nitrenyl ethynyl vsd MCH MNH
2043 % LocalWords:  chainsep resA atomA nbondsA resB atomB nbondsB newresA newresB
2044 % LocalWords:  rad deg lcc cc nm intramolecular forcefield PME Ewald
2045 % LocalWords:  solvation et groupconcept PHE TYR TRP equilibrated pre
2046 % LocalWords:  macromolecule disulfide harmonicbond Morsebond vsiteN
2047 % LocalWords:  cubicbond FENEbond tabulatedinteraction harmonicangle
2048 % LocalWords:  harmonicrestraint bondbondcross bondanglecross genion
2049 % LocalWords:  quarticangle properdihedral harmonicimproperdihedral
2050 % LocalWords:  RBdihedral periodicimproperdihedral Fourierdihedral
2051 % LocalWords:  positionrestraint distancerestraint dihedralrestraint
2052 % LocalWords:  orientationrestraint anglerestraint usinggroups ing
2053 % LocalWords:  DDeLoof MBAR Duan JACS Kollman Acc Hornak ILDN AMBERGS
2054 % LocalWords:  Lindorff Sanbonmatsu PNAS CMAP Schuler JCC pag
2055 % LocalWords:  aminoacid