Added pull potential flat-bottom-high
[gromacs.git] / docs / manual / special.tex
blobec0fac0ce7786aa920fc895bad7584a0811610ed
2 % This file is part of the GROMACS molecular simulation package.
4 % Copyright (c) 2013,2014,2015,2016, 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.
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.
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.
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.
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.
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.
35 \chapter{Special Topics}
36 \label{ch:special}
39 \section{Free energy implementation}
40 \label{sec:dgimplement}
41 For free energy calculations, there are two things that must be
42 specified; the end states, and the pathway connecting the end states.
43 The end states can be specified in two ways. The most straightforward
44 is through the specification of end states in the topology file. Most
45 potential forms support both an $A$ state and a $B$ state. Whenever both
46 states are specified, then the $A$ state corresponds to the initial free
47 energy state, and the $B$ state corresponds to the final state.
49 In some cases, the end state can also be defined in some cases without
50 altering the topology, solely through the {\tt .mdp} file, through the use
51 of the {\tt couple-moltype},{\tt couple-lambda0}, {\tt couple-lambda1}, and
52 {\tt couple-intramol} mdp keywords. Any molecule type selected in
53 {\tt couple-moltype} will automatically have a $B$ state implicitly
54 constructed (and the $A$ state redefined) according to the {\tt couple-lambda}
55 keywords. {\tt couple-lambda0} and {\tt couple-lambda1} define the non-bonded
56 parameters that are present in the $A$ state ({\tt couple-lambda0})
57 and the $B$ state ({\tt couple-lambda1}). The choices are 'q','vdw', and
58 'vdw-q'; these indicate the Coulombic, van der Waals, or both parameters
59 that are turned on in the respective state.
61 Once the end states are defined, then the path between the end states
62 has to be defined. This path is defined solely in the .mdp file.
63 Starting in 4.6, $\lambda$ is a vector of components, with Coulombic,
64 van der Waals, bonded, restraint, and mass components all able to be
65 adjusted independently. This makes it possible to turn off the
66 Coulombic term linearly, and then the van der Waals using soft core,
67 all in the same simulation. This is especially useful for replica
68 exchange or expanded ensemble simulations, where it is important to
69 sample all the way from interacting to non-interacting states in the
70 same simulation to improve sampling.
72 {\tt fep-lambdas} is the default array of $\lambda$ values ranging
73 from 0 to 1. All of the other lambda arrays use the values in this
74 array if they are not specified. The previous behavior, where the
75 pathway is controlled by a single $\lambda$ variable, can be preserved
76 by using only {\tt fep-lambdas} to define the pathway.
78 For example, if you wanted to first to change the Coulombic terms,
79 then the van der Waals terms, changing bonded at the same time rate as
80 the van der Waals, but changing the restraints throughout the first
81 two-thirds of the simulation, then you could use this $\lambda$ vector:
83 \begin{verbatim}
84 coul-lambdas = 0.0 0.2 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0
85 vdw-lambdas = 0.0 0.0 0.0 0.0 0.4 0.5 0.6 0.7 0.8 1.0
86 bonded-lambdas = 0.0 0.0 0.0 0.0 0.4 0.5 0.6 0.7 0.8 1.0
87 restraint-lambdas = 0.0 0.0 0.1 0.2 0.3 0.5 0.7 1.0 1.0 1.0
88 \end{verbatim}
90 This is also equivalent to:
92 \begin{verbatim}
93 fep-lambdas = 0.0 0.0 0.0 0.0 0.4 0.5 0.6 0.7 0.8 1.0
94 coul-lambdas = 0.0 0.2 0.5 1.0 1.0 1.0 1.0 1.0 1.0 1.0
95 restraint-lambdas = 0.0 0.0 0.1 0.2 0.3 0.5 0.7 1.0 1.0 1.0
96 \end{verbatim}
97 The {\tt fep-lambda array}, in this case, is being used as the default to
98 fill in the bonded and van der Waals $\lambda$ arrays. Usually, it's best to fill
99 in all arrays explicitly, just to make sure things are properly
100 assigned.
102 If you want to turn on only restraints going from $A$ to $B$, then it would be:
103 \begin{verbatim}
104 restraint-lambdas = 0.0 0.1 0.2 0.4 0.6 1.0
105 \end{verbatim}
106 and all of the other components of the $\lambda$ vector would be left in the $A$ state.
108 To compute free energies with a vector $\lambda$ using
109 thermodynamic integration, then the TI equation becomes vector equation:
110 \beq
111 \Delta F = \int \langle \nabla H \rangle \cdot d\vec{\lambda}
112 \eeq
113 or for finite differences:
114 \beq
115 \Delta F \approx \int \sum \langle \nabla H \rangle \cdot \Delta\lambda
116 \eeq
118 The external {\tt pymbar} script downloaded from https://SimTK.org/home/pymbar can
119 compute this integral automatically from the {\gromacs} dhdl.xvg output.
121 \section{Potential of mean force}
123 A potential of mean force (PMF) is a potential that is obtained
124 by integrating the mean force from an ensemble of configurations.
125 In {\gromacs}, there are several different methods to calculate the mean force.
126 Each method has its limitations, which are listed below.
127 \begin{itemize}
128 \item{\bf pull code:} between the centers of mass of molecules or groups of molecules.
129 \item{\bf free-energy code with harmonic bonds or constraints:} between single atoms.
130 \item{\bf free-energy code with position restraints:} changing the conformation of a relatively immobile group of atoms.
131 \item{\bf pull code in limited cases:} between groups of atoms that are
132 part of a larger molecule for which the bonds are constrained with
133 SHAKE or LINCS. If the pull group if relatively large,
134 the pull code can be used.
135 \end{itemize}
136 The pull and free-energy code a described in more detail
137 in the following two sections.
139 \subsubsection{Entropic effects}
140 When a distance between two atoms or the centers of mass of two groups
141 is constrained or restrained, there will be a purely entropic contribution
142 to the PMF due to the rotation of the two groups~\cite{RMNeumann1980a}.
143 For a system of two non-interacting masses the potential of mean force is:
144 \beq
145 V_{pmf}(r) = -(n_c - 1) k_B T \log(r)
146 \eeq
147 where $n_c$ is the number of dimensions in which the constraint works
148 (i.e. $n_c=3$ for a normal constraint and $n_c=1$ when only
149 the $z$-direction is constrained).
150 Whether one needs to correct for this contribution depends on what
151 the PMF should represent. When one wants to pull a substrate
152 into a protein, this entropic term indeed contributes to the work to
153 get the substrate into the protein. But when calculating a PMF
154 between two solutes in a solvent, for the purpose of simulating
155 without solvent, the entropic contribution should be removed.
156 {\bf Note} that this term can be significant; when at 300K the distance is halved,
157 the contribution is 3.5 kJ~mol$^{-1}$.
159 \section{Non-equilibrium pulling}
160 When the distance between two groups is changed continuously,
161 work is applied to the system, which means that the system is no longer
162 in equilibrium. Although in the limit of very slow pulling
163 the system is again in equilibrium, for many systems this limit
164 is not reachable within reasonable computational time.
165 However, one can use the Jarzynski relation~\cite{Jarzynski1997a}
166 to obtain the equilibrium free-energy difference $\Delta G$
167 between two distances from many non-equilibrium simulations:
168 \begin{equation}
169 \Delta G_{AB} = -k_BT \log \left\langle e^{-\beta W_{AB}} \right\rangle_A
170 \label{eq:Jarz}
171 \end{equation}
172 where $W_{AB}$ is the work performed to force the system along one path
173 from state A to B, the angular bracket denotes averaging over
174 a canonical ensemble of the initial state A and $\beta=1/k_B T$.
177 \section{The pull code}
178 \index{center-of-mass pulling}
179 \label{sec:pull}
180 The pull code applies forces or constraints between the centers
181 of mass of one or more pairs of groups of atoms.
182 Each pull reaction coordinate is called a ``coordinate'' and it operates
183 on usually two, but sometimes more, pull groups. A pull group can be part of one or more pull
184 coordinates. Furthermore, a coordinate can also operate on a single group
185 and an absolute reference position in space.
186 The distance between a pair of groups can be determined
187 in 1, 2 or 3 dimensions, or can be along a user-defined vector.
188 The reference distance can be constant or can change linearly with time.
189 Normally all atoms are weighted by their mass, but an additional
190 weighting factor can also be used.
191 \begin{figure}
192 \centerline{\includegraphics[width=6cm,angle=270]{plots/pull}}
193 \caption{Schematic picture of pulling a lipid out of a lipid bilayer
194 with umbrella pulling. $V_{rup}$ is the velocity at which the spring is
195 retracted, $Z_{link}$ is the atom to which the spring is attached and
196 $Z_{spring}$ is the location of the spring.}
197 \label{fig:pull}
198 \end{figure}
200 Several different pull types, i.e. ways to apply the pull force, are supported,
201 and in all cases the reference distance can be constant
202 or linearly changing with time.
203 \begin{enumerate}
204 \item{\textbf{Umbrella pulling}\swapindexquiet{umbrella}{pulling}}
205 A harmonic potential is applied between
206 the centers of mass of two groups.
207 Thus, the force is proportional to the displacement.
208 \item{\textbf{Constraint pulling\swapindexquiet{constraint}{pulling}}}
209 The distance between the centers of mass of two groups is constrained.
210 The constraint force can be written to a file.
211 This method uses the SHAKE algorithm but only needs 1 iteration to be
212 exact if only two groups are constrained.
213 \item{\textbf{Constant force pulling}}
214 A constant force is applied between the centers of mass of two groups.
215 Thus, the potential is linear.
216 In this case there is no reference distance of pull rate.
217 \item{\textbf{Flat bottom pulling}}
218 Like umbrella pulling, but the potential and force are zero for
219 coordinate values below ({\tt pull-coord?-type = flat-bottom}) or above
220 ({\tt pull-coord?-type = flat-bottom-high}) a reference value.
221 This is useful for restraining e.g. the distance
222 between two molecules to a certain region.
223 \end{enumerate}
225 \subsubsection{Definition of the center of mass}
227 In {\gromacs}, there are three ways to define the center of mass of a group.
228 The standard way is a ``plain'' center of mass, possibly with additional
229 weighting factors. With periodic boundary conditions it is no longer
230 possible to uniquely define the center of mass of a group of atoms.
231 Therefore, a reference atom is used. For determining the center of mass,
232 for all other atoms in the group, the closest periodic image to the reference
233 atom is used. This uniquely defines the center of mass.
234 By default, the middle (determined by the order in the topology) atom
235 is used as a reference atom, but the user can also select any other atom
236 if it would be closer to center of the group.
238 For a layered system, for instance a lipid bilayer, it may be of interest
239 to calculate the PMF of a lipid as function of its distance
240 from the whole bilayer. The whole bilayer can be taken as reference
241 group in that case, but it might also be of interest to define the
242 reaction coordinate for the PMF more locally. The {\tt .mdp} option
243 {\tt pull-coord?-geometry = cylinder} does not
244 use all the atoms of the reference group, but instead dynamically only those
245 within a cylinder with radius {\tt pull-cylinder-r} around the pull vector going
246 through the pull group. This only
247 works for distances defined in one dimension, and the cylinder is
248 oriented with its long axis along this one dimension. To avoid jumps in
249 the pull force, contributions of atoms are weighted as a function of distance
250 (in addition to the mass weighting):
251 \bea
252 w(r < r_\mathrm{cyl}) & = &
253 1-2 \left(\frac{r}{r_\mathrm{cyl}}\right)^2 + \left(\frac{r}{r_\mathrm{cyl}}\right)^4 \\
254 w(r \geq r_\mathrm{cyl}) & = & 0
255 \eea
256 Note that the radial dependence on the weight causes a radial force on
257 both cylinder group and the other pull group. This is an undesirable,
258 but unavoidable effect. To minimize this effect, the cylinder radius should
259 be chosen sufficiently large. The effective mass is 0.47 times that of
260 a cylinder with uniform weights and equal to the mass of uniform cylinder
261 of 0.79 times the radius.
263 \begin{figure}
264 \centerline{\includegraphics[width=6cm]{plots/pullref}}
265 \caption{Comparison of a plain center of mass reference group versus a cylinder
266 reference group applied to interface systems. C is the reference group.
267 The circles represent the center of mass of two groups plus the reference group,
268 $d_c$ is the reference distance.}
269 \label{fig:pullref}
270 \end{figure}
272 For a group of molecules in a periodic system, a plain reference group
273 might not be well-defined. An example is a water slab that is connected
274 periodically in $x$ and $y$, but has two liquid-vapor interfaces along $z$.
275 In such a setup, water molecules can evaporate from the liquid and they
276 will move through the vapor, through the periodic boundary, to the other
277 interface. Such a system is inherently periodic and there is no proper way
278 of defining a ``plain'' center of mass along $z$. A proper solution is to using
279 a cosine shaped weighting profile for all atoms in the reference group.
280 The profile is a cosine with a single period in the unit cell. Its phase
281 is optimized to give the maximum sum of weights, including mass weighting.
282 This provides a unique and continuous reference position that is nearly
283 identical to the plain center of mass position in case all atoms are all
284 within a half of the unit-cell length. See ref \cite{Engin2010a} for details.
286 When relative weights $w_i$ are used during the calculations, either
287 by supplying weights in the input or due to cylinder geometry
288 or due to cosine weighting,
289 the weights need to be scaled to conserve momentum:
290 \beq
291 w'_i = w_i
292 \left. \sum_{j=1}^N w_j \, m_j \right/ \sum_{j=1}^N w_j^2 \, m_j
293 \eeq
294 where $m_j$ is the mass of atom $j$ of the group.
295 The mass of the group, required for calculating the constraint force, is:
296 \beq
297 M = \sum_{i=1}^N w'_i \, m_i
298 \eeq
299 The definition of the weighted center of mass is:
300 \beq
301 \ve{r}_{com} = \left. \sum_{i=1}^N w'_i \, m_i \, \ve{r}_i \right/ M
302 \eeq
303 From the centers of mass the AFM, constraint, or umbrella force $\ve{F}_{\!com}$
304 on each group can be calculated.
305 The force on the center of mass of a group is redistributed to the atoms
306 as follows:
307 \beq
308 \ve{F}_{\!i} = \frac{w'_i \, m_i}{M} \, \ve{F}_{\!com}
309 \eeq
311 \subsubsection{Definition of the pull direction}
313 The most common setup is to pull along the direction of the vector containing
314 the two pull groups, this is selected with
315 {\tt pull-coord?-geometry = distance}. You might want to pull along a certain
316 vector instead, which is selected with {\tt pull-coord?-geometry = direction}.
317 But this can cause unwanted torque forces in the system, unless you pull against a reference group with (nearly) fixed orientation, e.g. a membrane protein embedded in a membrane along x/y while pulling along z. If your reference group does not have a fixed orientation, you should probably use
318 {\tt pull-coord?-geometry = direction-relative}, see \figref{pulldirrel}.
319 Since the potential now depends on the coordinates of two additional groups defining the orientation, the torque forces will work on these two groups.
321 \begin{figure}
322 \centerline{\includegraphics[width=5cm]{plots/pulldirrel}}
323 \caption{The pull setup for geometry {\tt direction-relative}. The ``normal'' pull groups are 1 and 2. Groups 3 and 4 define the pull direction and thus the direction of the normal pull forces (red). This leads to reaction forces (blue) on groups 3 and 4, which are perpendicular to the pull direction. Their magnitude is given by the ``normal'' pull force times the ratio of $d_p$ and the distance between groups 3 and 4.}
324 \label{fig:pulldirrel}
325 \end{figure}
328 \subsubsection{Limitations}
329 There is one theoretical limitation:
330 strictly speaking, constraint forces can only be calculated between
331 groups that are not connected by constraints to the rest of the system.
332 If a group contains part of a molecule of which the bond lengths
333 are constrained, the pull constraint and LINCS or SHAKE bond constraint
334 algorithms should be iterated simultaneously. This is not done in {\gromacs}.
335 This means that for simulations with {\tt constraints = all-bonds}
336 in the {\tt .mdp} file pulling is, strictly speaking,
337 limited to whole molecules or groups of molecules.
338 In some cases this limitation can be avoided by using the free energy code,
339 see \secref{fepmf}.
340 In practice, the errors caused by not iterating the two constraint
341 algorithms can be negligible when the pull group consists of a large
342 amount of atoms and/or the pull force is small.
343 In such cases, the constraint correction displacement of the pull group
344 is small compared to the bond lengths.
348 \section{\normindex{Enforced Rotation}}
349 \index{rotational pulling|see{enforced rotation}}
350 \index{pulling, rotational|see{enforced rotation}}
351 \label{sec:rotation}
353 \mathchardef\mhyphen="2D
354 \newcommand{\rotiso }{^\mathrm{iso}}
355 \newcommand{\rotisopf }{^\mathrm{iso\mhyphen pf}}
356 \newcommand{\rotpm }{^\mathrm{pm}}
357 \newcommand{\rotpmpf }{^\mathrm{pm\mhyphen pf}}
358 \newcommand{\rotrm }{^\mathrm{rm}}
359 \newcommand{\rotrmpf }{^\mathrm{rm\mhyphen pf}}
360 \newcommand{\rotrmtwo }{^\mathrm{rm2}}
361 \newcommand{\rotrmtwopf }{^\mathrm{rm2\mhyphen pf}}
362 \newcommand{\rotflex }{^\mathrm{flex}}
363 \newcommand{\rotflext }{^\mathrm{flex\mhyphen t}}
364 \newcommand{\rotflextwo }{^\mathrm{flex2}}
365 \newcommand{\rotflextwot}{^\mathrm{flex2\mhyphen t}}
367 This module can be used to enforce the rotation of a group of atoms, as {\eg}
368 a protein subunit. There are a variety of rotation potentials, among them
369 complex ones that allow flexible adaptations of both the rotated subunit as
370 well as the local rotation axis during the simulation. An example application
371 can be found in ref. \cite{Kutzner2011}.
373 \begin{figure}
374 \centerline{\includegraphics[width=13cm]{plots/rotation.pdf}}
375 \caption[Fixed and flexible axis rotation]{Comparison of fixed and flexible axis
376 rotation.
377 {\sf A:} Rotating the sketched shape inside the white tubular cavity can create
378 artifacts when a fixed rotation axis (dashed) is used. More realistically, the
379 shape would revolve like a flexible pipe-cleaner (dotted) inside the bearing (gray).
380 {\sf B:} Fixed rotation around an axis \ve{v} with a pivot point
381 specified by the vector \ve{u}.
382 {\sf C:} Subdividing the rotating fragment into slabs with separate rotation
383 axes ($\uparrow$) and pivot points ($\bullet$) for each slab allows for
384 flexibility. The distance between two slabs with indices $n$ and $n+1$ is $\Delta x$.}
385 \label{fig:rotation}
386 \end{figure}
388 \begin{figure}
389 \centerline{\includegraphics[width=13cm]{plots/equipotential.pdf}}
390 \caption{Selection of different rotation potentials and definition of notation.
391 All four potentials $V$ (color coded) are shown for a single atom at position
392 $\ve{x}_j(t)$.
393 {\sf A:} Isotropic potential $V\rotiso$,
394 {\sf B:} radial motion potential $V\rotrm$ and flexible potential
395 $V\rotflex$,
396 {\sf C--D:} radial motion\,2 potential $V\rotrmtwo$ and
397 flexible\,2 potential $V\rotflextwo$ for $\epsilon' = 0$\,nm$^2$ {\sf (C)}
398 and $\epsilon' = 0.01$\,nm$^2$ {\sf (D)}. The rotation axis is perpendicular to
399 the plane and marked by $\otimes$. The light gray contours indicate Boltzmann factors
400 $e^{-V/(k_B T)}$ in the $\ve{x}_j$-plane for $T=300$\,K and
401 $k=200$\,kJ/(mol$\cdot$nm$^2$). The green arrow shows the direction of the
402 force $\ve{F}_{\!j}$ acting on atom $j$; the blue dashed line indicates the
403 motion of the reference position.}
404 \label{fig:equipotential}
405 \end{figure}
407 \subsection{Fixed Axis Rotation}
408 \subsubsection{Stationary Axis with an Isotropic Potential}
409 In the fixed axis approach (see \figref{rotation}B), torque on a group of $N$
410 atoms with positions $\ve{x}_i$ (denoted ``rotation group'') is applied by
411 rotating a reference set of atomic positions -- usually their initial positions
412 $\ve{y}_i^0$ -- at a constant angular velocity $\omega$ around an axis
413 defined by a direction vector $\hat{\ve{v}}$ and a pivot point \ve{u}.
414 To that aim, each atom with position $\ve{x}_i$ is attracted by a
415 ``virtual spring'' potential to its moving reference position
416 $\ve{y}_i = \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{u})$,
417 where $\mathbf{\Omega}(t)$ is a matrix that describes the rotation around the
418 axis. In the simplest case, the ``springs'' are described by a harmonic
419 potential,
420 \beq
421 V\rotiso = \frac{k}{2} \sum_{i=1}^{N} w_i \left[ \mathbf{\Omega}(t)
422 (\ve{y}_i^0 - \ve{u}) - (\ve{x}_i - \ve{u}) \right]^2 ,
423 \label{eqn:potiso}
424 \eeq
425 with optional mass-weighted prefactors $w_i = N \, m_i/M$ with total mass
426 $M = \sum_{i=1}^N m_i$.
427 The rotation matrix $\mathbf{\Omega}(t)$ is
428 \newcommand{\omcost}{\,\xi\,} % abbreviation
429 \begin{displaymath}
430 \mathbf{\Omega}(t) =
431 \left(
432 \begin{array}{ccc}
433 \cos\omega t + v_x^2\omcost & v_x v_y\omcost - v_z\sin\omega t & v_x v_z\omcost + v_y\sin\omega t\\
434 v_x v_y\omcost + v_z\sin\omega t & \cos\omega t + v_y^2\omcost & v_y v_z\omcost - v_x\sin\omega t\\
435 v_x v_z\omcost - v_y\sin\omega t & v_y v_z\omcost + v_x\sin\omega t & \cos\omega t + v_z^2\omcost \\
436 \end{array}
437 \right) ,
438 \end{displaymath}
439 where $v_x$, $v_y$, and $v_z$ are the components of the normalized rotation vector
440 $\hat{\ve{v}}$, and $\omcost := 1-\cos(\omega t)$. As illustrated in
441 \figref{equipotential}A for a single atom $j$, the
442 rotation matrix $\mathbf{\Omega}(t)$ operates on the initial reference positions
443 $\ve{y}_j^0 = \ve{x}_j(t_0)$ of atom $j$ at $t=t_0$. At a later
444 time $t$, the reference position has rotated away from its initial place
445 (along the blue dashed line), resulting in the force
446 \beq
447 \ve{F}_{\!j}\rotiso
448 = -\nabla_{\!j} \, V\rotiso
449 = k \, w_j \left[
450 \mathbf{\Omega}(t) (\ve{y}_j^0 - \ve{u}) - (\ve{x}_j - \ve{u} ) \right] ,
451 \label{eqn:force_fixed}
452 \eeq
453 which is directed towards the reference position.
456 \subsubsection{Pivot-Free Isotropic Potential}
457 Instead of a fixed pivot vector \ve{u} this potential uses the center of
458 mass $\ve{x}_c$ of the rotation group as pivot for the rotation axis,
459 \beq
460 \ve{x}_c = \frac{1}{M} \sum_{i=1}^N m_i \ve{x}_i
461 \label{eqn:com}
462 \mbox{\hspace{4ex}and\hspace{4ex}}
463 \ve{y}_c^0 = \frac{1}{M} \sum_{i=1}^N m_i \ve{y}_i^0 \ ,
464 \eeq
465 which yields the ``pivot-free'' isotropic potential
466 \beq
467 V\rotisopf = \frac{k}{2} \sum_{i=1}^{N} w_i \left[ \mathbf{\Omega}(t)
468 (\ve{y}_i^0 - \ve{y}_c^0) - (\ve{x}_i - \ve{x}_c) \right]^2 ,
469 \label{eqn:potisopf}
470 \eeq
471 with forces
472 \beq
473 \mathbf{F}_{\!j}\rotisopf = k \, w_j
474 \left[
475 \mathbf{\Omega}(t) ( \ve{y}_j^0 - \ve{y}_c^0)
476 - ( \ve{x}_j - \ve{x}_c )
477 \right] .
478 \label{eqn:force_isopf}
479 \eeq
480 Without mass-weighting, the pivot $\ve{x}_c$ is the geometrical center of
481 the group.
482 \label{sec:fixed}
484 \subsubsection{Parallel Motion Potential Variant}
485 The forces generated by the isotropic potentials
486 (\eqnsref{potiso}{potisopf}) also contain components parallel
487 to the rotation axis and thereby restrain motions along the axis of either the
488 whole rotation group (in case of $V\rotiso$) or within
489 the rotation group (in case of $V\rotisopf$). For cases where
490 unrestrained motion along the axis is preferred, we have implemented a
491 ``parallel motion'' variant by eliminating all components parallel to the
492 rotation axis for the potential. This is achieved by projecting the distance
493 vectors between reference and actual positions
494 \beq
495 \ve{r}_i = \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{u}) - (\ve{x}_i - \ve{u})
496 \eeq
497 onto the plane perpendicular to the rotation vector,
498 \beq
499 \label{eqn:project}
500 \ve{r}_i^\perp := \ve{r}_i - (\ve{r}_i \cdot \hat{\ve{v}})\hat{\ve{v}} \ ,
501 \eeq
502 yielding
503 \bea
504 \nonumber
505 V\rotpm &=& \frac{k}{2} \sum_{i=1}^{N} w_i ( \ve{r}_i^\perp )^2 \\
506 &=& \frac{k}{2} \sum_{i=1}^{N} w_i
507 \left\lbrace
508 \mathbf{\Omega}(t)
509 (\ve{y}_i^0 - \ve{u}) - (\ve{x}_i - \ve{u}) \right. \nonumber \\
510 && \left. - \left\lbrace
511 \left[ \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{u}) - (\ve{x}_i - \ve{u}) \right] \cdot\hat{\ve{v}}
512 \right\rbrace\hat{\ve{v}} \right\rbrace^2 ,
513 \label{eqn:potpm}
514 \eea
515 and similarly
516 \beq
517 \ve{F}_{\!j}\rotpm = k \, w_j \, \ve{r}_j^\perp .
518 \label{eqn:force_pm}
519 \eeq
521 \subsubsection{Pivot-Free Parallel Motion Potential}
522 Replacing in \eqnref{potpm} the fixed pivot \ve{u} by the center
523 of mass $\ve{x_c}$ yields the pivot-free variant of the parallel motion
524 potential. With
525 \beq
526 \ve{s}_i = \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{y}_c^0) - (\ve{x}_i - \ve{x}_c)
527 \eeq
528 the respective potential and forces are
529 \bea
530 V\rotpmpf &=& \frac{k}{2} \sum_{i=1}^{N} w_i ( \ve{s}_i^\perp )^2 \ , \\
531 \label{eqn:potpmpf}
532 \ve{F}_{\!j}\rotpmpf &=& k \, w_j \, \ve{s}_j^\perp .
533 \label{eqn:force_pmpf}
534 \eea
536 \subsubsection{Radial Motion Potential}
537 In the above variants, the minimum of the rotation potential is either a single
538 point at the reference position $\ve{y}_i$ (for the isotropic potentials) or a
539 single line through $\ve{y}_i$ parallel to the rotation axis (for the
540 parallel motion potentials). As a result, radial forces restrict radial motions
541 of the atoms. The two subsequent types of rotation potentials, $V\rotrm$
542 and $V\rotrmtwo$, drastically reduce or even eliminate this effect. The first
543 variant, $V\rotrm$ (\figref{equipotential}B), eliminates all force
544 components parallel to the vector connecting the reference atom and the
545 rotation axis,
546 \beq
547 V\rotrm = \frac{k}{2} \sum_{i=1}^{N} w_i \left[
548 \ve{p}_i
549 \cdot(\ve{x}_i - \ve{u}) \right]^2 ,
550 \label{eqn:potrm}
551 \eeq
552 with
553 \beq
554 \ve{p}_i :=
555 \frac{\hat{\ve{v}}\times \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{u})} {\| \hat{\ve{v}}\times \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{u})\|} \ .
556 \eeq
557 This variant depends only on the distance $\ve{p}_i \cdot (\ve{x}_i -
558 \ve{u})$ of atom $i$ from the plane spanned by $\hat{\ve{v}}$ and
559 $\mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{u})$. The resulting force is
560 \beq
561 \mathbf{F}_{\!j}\rotrm =
562 -k \, w_j \left[ \ve{p}_j\cdot(\ve{x}_j - \ve{u}) \right] \,\ve{p}_j \, .
563 \label{eqn:potrm_force}
564 \eeq
566 \subsubsection{Pivot-Free Radial Motion Potential}
567 Proceeding similar to the pivot-free isotropic potential yields a pivot-free
568 version of the above potential. With
569 \beq
570 \ve{q}_i :=
571 \frac{\hat{\ve{v}}\times \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{y}_c^0)} {\| \hat{\ve{v}}\times \mathbf{\Omega}(t) (\ve{y}_i^0 - \ve{y}_c^0)\|} \, ,
572 \eeq
573 the potential and force for the pivot-free variant of the radial motion potential read
574 \bea
575 V\rotrmpf & = & \frac{k}{2} \sum_{i=1}^{N} w_i \left[
576 \ve{q}_i
577 \cdot(\ve{x}_i - \ve{x}_c)
578 \right]^2 \, , \\
579 \label{eqn:potrmpf}
580 \mathbf{F}_{\!j}\rotrmpf & = &
581 -k \, w_j \left[ \ve{q}_j\cdot(\ve{x}_j - \ve{x}_c) \right] \,\ve{q}_j
582 + k \frac{m_j}{M} \sum_{i=1}^{N} w_i \left[
583 \ve{q}_i\cdot(\ve{x}_i - \ve{x}_c) \right]\,\ve{q}_i \, .
584 \label{eqn:potrmpf_force}
585 \eea
587 \subsubsection{Radial Motion 2 Alternative Potential}
588 As seen in \figref{equipotential}B, the force resulting from
589 $V\rotrm$ still contains a small, second-order radial component. In most
590 cases, this perturbation is tolerable; if not, the following
591 alternative, $V\rotrmtwo$, fully eliminates the radial contribution to the
592 force, as depicted in \figref{equipotential}C,
593 \beq
594 V\rotrmtwo =
595 \frac{k}{2} \sum_{i=1}^{N} w_i\,
596 \frac{\left[ (\hat{\ve{v}} \times ( \ve{x}_i - \ve{u} ))
597 \cdot \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{u}) \right]^2}
598 {\| \hat{\ve{v}} \times (\ve{x}_i - \ve{u}) \|^2 +
599 \epsilon'} \, ,
600 \label{eqn:potrm2}
601 \eeq
602 where a small parameter $\epsilon'$ has been introduced to avoid singularities.
603 For $\epsilon'=0$\,nm$^2$, the equipotential planes are spanned by $\ve{x}_i -
604 \ve{u}$ and $\hat{\ve{v}}$, yielding a force
605 perpendicular to $\ve{x}_i - \ve{u}$, thus not contracting or
606 expanding structural parts that moved away from or toward the rotation axis.
608 Choosing a small positive $\epsilon'$ ({\eg},
609 $\epsilon'=0.01$\,nm$^2$, \figref{equipotential}D) in the denominator of
610 \eqnref{potrm2} yields a well-defined potential and continuous forces also
611 close to the rotation axis, which is not the case for $\epsilon'=0$\,nm$^2$
612 (\figref{equipotential}C). With
613 \bea
614 \ve{r}_i & := & \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{u})\\
615 \ve{s}_i & := & \frac{\hat{\ve{v}} \times (\ve{x}_i -
616 \ve{u} ) }{ \| \hat{\ve{v}} \times (\ve{x}_i - \ve{u})
617 \| } \equiv \; \Psi_{i} \;\; {\hat{\ve{v}} \times
618 (\ve{x}_i-\ve{u} ) }\\
619 \Psi_i^{*} & := & \frac{1}{ \| \hat{\ve{v}} \times
620 (\ve{x}_i-\ve{u}) \|^2 + \epsilon'}
621 \eea
622 the force on atom $j$ reads
623 \beq
624 \ve{F}_{\!j}\rotrmtwo =
625 - k\;
626 \left\lbrace w_j\;
627 (\ve{s}_j\cdot\ve{r}_{\!j})\;
628 \left[ \frac{\Psi_{\!j}^* }{\Psi_{\!j} } \ve{r}_{\!j}
629 - \frac{\Psi_{\!j}^{*2}}{\Psi_{\!j}^3}
630 (\ve{s}_j\cdot\ve{r}_{\!j})\ve{s}_j \right]
631 \right\rbrace \times \hat{\ve{v}} .
632 \label{eqn:potrm2_force}
633 \eeq
635 \subsubsection{Pivot-Free Radial Motion 2 Potential}
636 The pivot-free variant of the above potential is
637 \beq
638 V\rotrmtwopf =
639 \frac{k}{2} \sum_{i=1}^{N} w_i\,
640 \frac{\left[ (\hat{\ve{v}} \times ( \ve{x}_i - \ve{x}_c ))
641 \cdot \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{y}_c) \right]^2}
642 {\| \hat{\ve{v}} \times (\ve{x}_i - \ve{x}_c) \|^2 +
643 \epsilon'} \, .
644 \label{eqn:potrm2pf}
645 \eeq
646 With
647 \bea
648 \ve{r}_i & := & \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{y}_c)\\
649 \ve{s}_i & := & \frac{\hat{\ve{v}} \times (\ve{x}_i -
650 \ve{x}_c ) }{ \| \hat{\ve{v}} \times (\ve{x}_i - \ve{x}_c)
651 \| } \equiv \; \Psi_{i} \;\; {\hat{\ve{v}} \times
652 (\ve{x}_i-\ve{x}_c ) }\\ \Psi_i^{*} & := & \frac{1}{ \| \hat{\ve{v}} \times
653 (\ve{x}_i-\ve{x}_c) \|^2 + \epsilon'}
654 \eea
655 the force on atom $j$ reads
656 \bea
657 \nonumber
658 \ve{F}_{\!j}\rotrmtwopf & = &
659 - k\;
660 \left\lbrace w_j\;
661 (\ve{s}_j\cdot\ve{r}_{\!j})\;
662 \left[ \frac{\Psi_{\!j}^* }{\Psi_{\!j} } \ve{r}_{\!j}
663 - \frac{\Psi_{\!j}^{*2}}{\Psi_{\!j}^3}
664 (\ve{s}_j\cdot\ve{r}_{\!j})\ve{s}_j \right]
665 \right\rbrace \times \hat{\ve{v}}\\
667 + k\;\frac{m_j}{M} \left\lbrace \sum_{i=1}^{N}
668 w_i\;(\ve{s}_i\cdot\ve{r}_i) \;
669 \left[ \frac{\Psi_i^* }{\Psi_i } \ve{r}_i
670 - \frac{\Psi_i^{*2}}{\Psi_i^3} (\ve{s}_i\cdot\ve{r}_i )\;
671 \ve{s}_i \right] \right\rbrace \times \hat{\ve{v}} \, .
672 \label{eqn:potrm2pf_force}
673 \eea
675 \subsection{Flexible Axis Rotation}
676 As sketched in \figref{rotation}A--B, the rigid body behavior of
677 the fixed axis rotation scheme is a drawback for many applications. In
678 particular, deformations of the rotation group are suppressed when the
679 equilibrium atom positions directly depend on the reference positions.
680 To avoid this limitation, \eqnsref{potrmpf}{potrm2pf}
681 will now be generalized towards a ``flexible axis'' as sketched in
682 \figref{rotation}C. This will be achieved by subdividing the
683 rotation group into a set of equidistant slabs perpendicular to
684 the rotation vector, and by applying a separate rotation potential to each
685 of these slabs. \figref{rotation}C shows the midplanes of the slabs
686 as dotted straight lines and the centers as thick black dots.
688 To avoid discontinuities in the potential and in the forces, we define
689 ``soft slabs'' by weighing the contributions of each
690 slab $n$ to the total potential function $V\rotflex$ by a Gaussian
691 function
692 \beq
693 \label{eqn:gaussian}
694 g_n(\ve{x}_i) = \Gamma \ \mbox{exp} \left(
695 -\frac{\beta_n^2(\ve{x}_i)}{2\sigma^2} \right) ,
696 \eeq
697 centered at the midplane of the $n$th slab. Here $\sigma$ is the width
698 of the Gaussian function, $\Delta x$ the distance between adjacent slabs, and
699 \beq
700 \beta_n(\ve{x}_i) := \ve{x}_i \cdot \hat{\ve{v}} - n \, \Delta x \, .
701 \eeq
703 \begin{figure}
704 \centerline{\includegraphics[width=6.5cm]{plots/gaussians.pdf}}
705 \caption{Gaussian functions $g_n$ centered at $n \, \Delta x$ for a slab
706 distance $\Delta x = 1.5$ nm and $n \geq -2$. Gaussian function $g_0$ is
707 highlighted in bold; the dashed line depicts the sum of the shown Gaussian
708 functions.}
709 \label{fig:gaussians}
710 \end{figure}
712 A most convenient choice is $\sigma = 0.7 \Delta x$ and
713 \begin{displaymath}
714 1/\Gamma = \sum_{n \in Z}
715 \mbox{exp}
716 \left(-\frac{(n - \frac{1}{4})^2}{2\cdot 0.7^2}\right)
717 \approx 1.75464 \, ,
718 \end{displaymath}
719 which yields a nearly constant sum, essentially independent of $\ve{x}_i$
720 (dashed line in \figref{gaussians}), {\ie},
721 \beq
722 \sum_{n \in Z} g_n(\ve{x}_i) = 1 + \epsilon(\ve{x}_i) \, ,
723 \label{eqn:normal}
724 \eeq
725 with $ | \epsilon(\ve{x}_i) | < 1.3\cdot 10^{-4}$. This choice also
726 implies that the individual contributions to the force from the slabs add up to
727 unity such that no further normalization is required.
729 To each slab center $\ve{x}_c^n$, all atoms contribute by their
730 Gaussian-weighted (optionally also mass-weighted) position vectors
731 $g_n(\ve{x}_i) \, \ve{x}_i$. The
732 instantaneous slab centers $\ve{x}_c^n$ are calculated from the
733 current positions $\ve{x}_i$,
734 \beq
735 \label{eqn:defx0}
736 \ve{x}_c^n =
737 \frac{\sum_{i=1}^N g_n(\ve{x}_i) \, m_i \, \ve{x}_i}
738 {\sum_{i=1}^N g_n(\ve{x}_i) \, m_i} \, ,\\
739 \eeq
740 while the reference centers $\ve{y}_c^n$ are calculated from the reference
741 positions $\ve{y}_i^0$,
742 \beq
743 \label{eqn:defy0}
744 \ve{y}_c^n =
745 \frac{\sum_{i=1}^N g_n(\ve{y}_i^0) \, m_i \, \ve{y}_i^0}
746 {\sum_{i=1}^N g_n(\ve{y}_i^0) \, m_i} \, .
747 \eeq
748 Due to the rapid decay of $g_n$, each slab
749 will essentially involve contributions from atoms located within $\approx
750 3\Delta x$ from the slab center only.
752 \subsubsection{Flexible Axis Potential}
753 We consider two flexible axis variants. For the first variant,
754 the slab segmentation procedure with Gaussian weighting is applied to the radial
755 motion potential (\eqnref{potrmpf}\,/\,\figref{equipotential}B),
756 yielding as the contribution of slab $n$
757 \begin{displaymath}
758 V^n =
759 \frac{k}{2} \sum_{i=1}^{N} w_i \, g_n(\ve{x}_i)
760 \left[
761 \ve{q}_i^n
762 \cdot
763 (\ve{x}_i - \ve{x}_c^n)
764 \right]^2 ,
765 \label{eqn:flexpot}
766 \end{displaymath}
767 and a total potential function
768 \beq
769 V\rotflex = \sum_n V^n \, .
770 \label{eqn:potflex}
771 \eeq
772 Note that the global center of mass $\ve{x}_c$ used in
773 \eqnref{potrmpf} is now replaced by $\ve{x}_c^n$, the center of mass of
774 the slab. With
775 \bea
776 \ve{q}_i^n & := & \frac{\hat{\ve{v}} \times
777 \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{y}_c^n) }{ \| \hat{\ve{v}}
778 \times \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{y}_c^n) \| } \\
779 b_i^n & := & \ve{q}_i^n \cdot (\ve{x}_i - \ve{x}_c^n) \, ,
780 \eea
781 the resulting force on atom $j$ reads
782 \bea
783 \nonumber\hspace{-15mm}
784 \ve{F}_{\!j}\rotflex &=&
785 - \, k \, w_j \sum_n g_n(\ve{x}_j) \, b_j^n \left\lbrace \ve{q}_j^n -
786 b_j^n \frac{\beta_n(\ve{x}_j)}{2\sigma^2} \hat{\ve{v}} \right\rbrace \\ & &
787 + \, k \, m_j \sum_n \frac{g_n(\ve{x}_j)}{\sum_h g_n(\ve{x}_h)}
788 \sum_{i=1}^{N} w_i \, g_n(\ve{x}_i) \, b_i^n \left\lbrace
789 \ve{q}_i^n -\frac{\beta_n(\ve{x}_j)}{\sigma^2}
790 \left[ \ve{q}_i^n \cdot (\ve{x}_j - \ve{x}_c^n )\right]
791 \hat{\ve{v}} \right\rbrace .
792 \label{eqn:potflex_force}
793 \eea
795 Note that for $V\rotflex$, as defined, the slabs are fixed in space and so
796 are the reference centers $\ve{y}_c^n$. If during the simulation the
797 rotation group moves too far in $\ve{v}$ direction, it may enter a
798 region where -- due to the lack of nearby reference positions -- no reference
799 slab centers are defined, rendering the potential evaluation impossible.
800 We therefore have included a slightly modified version of this potential that
801 avoids this problem by attaching the midplane of slab $n=0$ to the center of mass
802 of the rotation group, yielding slabs that move with the rotation group.
803 This is achieved by subtracting the center of mass $\ve{x}_c$ of the
804 group from the positions,
805 \beq
806 \tilde{\ve{x}}_i = \ve{x}_i - \ve{x}_c \, , \mbox{\ \ \ and \ \ }
807 \tilde{\ve{y}}_i^0 = \ve{y}_i^0 - \ve{y}_c^0 \, ,
808 \label{eqn:trafo}
809 \eeq
810 such that
811 \bea
812 V\rotflext
813 & = & \frac{k}{2} \sum_n \sum_{i=1}^{N} w_i \, g_n(\tilde{\ve{x}}_i)
814 \left[ \frac{\hat{\ve{v}} \times \mathbf{\Omega}(t)(\tilde{\ve{y}}_i^0
815 - \tilde{\ve{y}}_c^n) }{ \| \hat{\ve{v}} \times
816 \mathbf{\Omega}(t)(\tilde{\ve{y}}_i^0 -
817 \tilde{\ve{y}}_c^n) \| }
818 \cdot
819 (\tilde{\ve{x}}_i - \tilde{\ve{x}}_c^n)
820 \right]^2 .
821 \label{eqn:potflext}
822 \eea
823 To simplify the force derivation, and for efficiency reasons, we here assume
824 $\ve{x}_c$ to be constant, and thus $\partial \ve{x}_c / \partial x =
825 \partial \ve{x}_c / \partial y = \partial \ve{x}_c / \partial z = 0$. The
826 resulting force error is small (of order $O(1/N)$ or $O(m_j/M)$ if
827 mass-weighting is applied) and can therefore be tolerated. With this assumption,
828 the forces $\ve{F}\rotflext$ have the same form as
829 \eqnref{potflex_force}.
831 \subsubsection{Flexible Axis 2 Alternative Potential}
832 In this second variant, slab segmentation is applied to $V\rotrmtwo$
833 (\eqnref{potrm2pf}), resulting in a flexible axis potential without radial
834 force contributions (\figref{equipotential}C),
835 \beq
836 V\rotflextwo =
837 \frac{k}{2} \sum_{i=1}^{N} \sum_n w_i\,g_n(\ve{x}_i)
838 \frac{\left[ (\hat{\ve{v}} \times ( \ve{x}_i - \ve{x}_c^n ))
839 \cdot \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{y}_c^n) \right]^2}
840 {\| \hat{\ve{v}} \times (\ve{x}_i - \ve{x}_c^n) \|^2 +
841 \epsilon'} \, .
842 \label{eqn:potflex2}
843 \eeq
844 With
845 \bea
846 \ve{r}_i^n & := & \mathbf{\Omega}(t)(\ve{y}_i^0 - \ve{y}_c^n)\\
847 \ve{s}_i^n & := & \frac{\hat{\ve{v}} \times (\ve{x}_i -
848 \ve{x}_c^n ) }{ \| \hat{\ve{v}} \times (\ve{x}_i - \ve{x}_c^n)
849 \| } \equiv \; \psi_{i} \;\; {\hat{\ve{v}} \times (\ve{x}_i-\ve{x}_c^n ) }\\
850 \psi_i^{*} & := & \frac{1}{ \| \hat{\ve{v}} \times (\ve{x}_i-\ve{x}_c^n) \|^2 + \epsilon'}\\
851 W_j^n & := & \frac{g_n(\ve{x}_j)\,m_j}{\sum_h g_n(\ve{x}_h)\,m_h}\\
852 \ve{S}^n & := &
853 \sum_{i=1}^{N} w_i\;g_n(\ve{x}_i)
854 \; (\ve{s}_i^n\cdot\ve{r}_i^n)
855 \left[ \frac{\psi_i^* }{\psi_i } \ve{r}_i^n
856 - \frac{\psi_i^{*2}}{\psi_i^3} (\ve{s}_i^n\cdot\ve{r}_i^n )\;
857 \ve{s}_i^n \right] \label{eqn:Sn}
858 \eea
859 the force on atom $j$ reads
860 \bea
861 \nonumber
862 \ve{F}_{\!j}\rotflextwo & = &
863 - k\;
864 \left\lbrace \sum_n w_j\;g_n(\ve{x}_j)\;
865 (\ve{s}_j^n\cdot\ve{r}_{\!j}^n)\;
866 \left[ \frac{\psi_j^* }{\psi_j } \ve{r}_{\!j}^n
867 - \frac{\psi_j^{*2}}{\psi_j^3} (\ve{s}_j^n\cdot\ve{r}_{\!j}^n)\;
868 \ve{s}_{\!j}^n \right] \right\rbrace \times \hat{\ve{v}} \\
869 \nonumber
871 + k \left\lbrace \sum_n W_{\!j}^n \, \ve{S}^n \right\rbrace \times
872 \hat{\ve{v}}
873 - k \left\lbrace \sum_n W_{\!j}^n \; \frac{\beta_n(\ve{x}_j)}{\sigma^2} \frac{1}{\psi_j}\;\;
874 \ve{s}_j^n \cdot
875 \ve{S}^n \right\rbrace \hat{\ve{v}}\\
876 & &
877 + \frac{k}{2} \left\lbrace \sum_n w_j\;g_n(\ve{x}_j)
878 \frac{\beta_n(\ve{x}_j)}{\sigma^2}
879 \frac{\psi_j^*}{\psi_j^2}( \ve{s}_j^n \cdot \ve{r}_{\!j}^n )^2 \right\rbrace
880 \hat{\ve{v}} .
881 \label{eqn:potflex2_force}
882 \eea
884 Applying transformation (\ref{eqn:trafo}) yields a ``translation-tolerant''
885 version of the flexible\,2 potential, $V\rotflextwot$. Again,
886 assuming that $\partial \ve{x}_c / \partial x$, $\partial \ve{x}_c /
887 \partial y$, $\partial \ve{x}_c / \partial z$ are small, the
888 resulting equations for $V\rotflextwot$ and $\ve{F}\rotflextwot$ are
889 similar to those of $V\rotflextwo$ and $\ve{F}\rotflextwo$.
891 \subsection{Usage}
892 To apply enforced rotation, the particles $i$ that are to
893 be subjected to one of the rotation potentials are defined via index groups
894 {\tt rot-group0}, {\tt rot-group1}, etc., in the {\tt .mdp} input file.
895 The reference positions $\ve{y}_i^0$ are
896 read from a special {\tt .trr} file provided to {\tt grompp}. If no such file is found,
897 $\ve{x}_i(t=0)$ are used as reference positions and written to {\tt .trr} such
898 that they can be used for subsequent setups. All parameters of the potentials
899 such as $k$, $\epsilon'$, etc. (\tabref{vars}) are provided as {\tt .mdp}
900 parameters; {\tt rot-type} selects the type of the potential.
901 The option {\tt rot-massw} allows to choose whether or not to use
902 mass-weighted averaging.
903 For the flexible potentials, a cutoff value $g_n^\mathrm{min}$
904 (typically $g_n^\mathrm{min}=0.001$) makes shure that only
905 significant contributions to $V$ and \ve{F} are evaluated, {\ie} terms with
906 $g_n(\ve{x}) < g_n^\mathrm{min}$ are omitted.
907 \tabref{quantities} summarizes observables that are written
908 to additional output files and which are described below.
911 \begin{table}[tbp]
912 \caption{Parameters used by the various rotation potentials.
913 {\sf x}'s indicate which parameter is actually used for a given potential.}
914 %\small
916 \newcommand{\kunit}{$\frac{\mathrm{kJ}}{\mathrm{mol} \cdot \mathrm{nm}^2}$}
917 \newcommand{\smtt}[1]{{\hspace{-0.5ex}\small #1\hspace{-0.5ex}}}
918 \label{tab:vars}
919 \begin{center}
920 \begin{tabular}{l>{$}l<{$}rccccccc}
921 \hline
922 parameter & & & $k$ & $\hat{\ve{v}}$ & $\ve{u}$ & $\omega$ & $\epsilon'$ & $\Delta x$ & $g_n^\mathrm{min}$ \\
923 \multicolumn{3}{l}{{\tt .mdp} input variable name} & \smtt{k} & \smtt{vec} & \smtt{pivot} & \smtt{rate} & \smtt{eps} & \smtt{slab-dist} & \smtt{min-gauss} \\
924 unit & & & \kunit & - & nm & $^\circ$/ps & nm$^2$ & nm & \,-\, \\[1mm]
925 \hline \multicolumn{2}{l}{fixed axis potentials:} & eqn.\\
926 isotropic & V\rotiso & (\ref{eqn:potiso}) & {\sf x} & {\sf x} & {\sf x} & {\sf x} & - & - & - \\
927 --- pivot-free & V\rotisopf & (\ref{eqn:potisopf}) & {\sf x} & {\sf x} & - & {\sf x} & - & - & - \\
928 parallel motion & V\rotpm & (\ref{eqn:potpm}) & {\sf x} & {\sf x} & {\sf x} & {\sf x} & - & - & - \\
929 --- pivot-free & V\rotpmpf & (\ref{eqn:potpmpf}) & {\sf x} & {\sf x} & - & {\sf x} & - & - & - \\
930 radial motion & V\rotrm & (\ref{eqn:potrm}) & {\sf x} & {\sf x} & {\sf x} & {\sf x} & - & - & - \\
931 --- pivot-free & V\rotrmpf & (\ref{eqn:potrmpf}) & {\sf x} & {\sf x} & - & {\sf x} & - & - & - \\
932 radial motion\,2 & V\rotrmtwo & (\ref{eqn:potrm2}) & {\sf x} & {\sf x} & {\sf x} & {\sf x} & {\sf x} & - & - \\
933 --- pivot-free & V\rotrmtwopf & (\ref{eqn:potrm2pf}) & {\sf x} & {\sf x} & - & {\sf x} & {\sf x} & - & - \\ \hline
934 \multicolumn{2}{l}{flexible axis potentials:} & eqn.\\
935 flexible & V\rotflex & (\ref{eqn:potflex}) & {\sf x} & {\sf x} & - & {\sf x} & - & {\sf x} & {\sf x} \\
936 --- transl. tol. & V\rotflext & (\ref{eqn:potflext}) & {\sf x} & {\sf x} & - & {\sf x} & - & {\sf x} & {\sf x} \\
937 flexible\,2 & V\rotflextwo & (\ref{eqn:potflex2}) & {\sf x} & {\sf x} & - & {\sf x} & {\sf x} & {\sf x} & {\sf x} \\
938 --- transl. tol. & V\rotflextwot & - & {\sf x} & {\sf x} & - & {\sf x} & {\sf x} & {\sf x} & {\sf x} \\
939 \hline
940 \end{tabular}
941 \end{center}
942 \end{table}
944 \begin{table}
945 \caption{Quantities recorded in output files during enforced rotation simulations.
946 All slab-wise data is written every {\tt nstsout} steps, other rotation data every {\tt nstrout} steps.}
947 \label{tab:quantities}
948 \begin{center}
949 \begin{tabular}{llllcc}
950 \hline
951 quantity & unit & equation & output file & fixed & flexible\\ \hline
952 $V(t)$ & kJ/mol & see \ref{tab:vars} & {\tt rotation} & {\sf x} & {\sf x} \\
953 $\theta_\mathrm{ref}(t)$ & degrees & $\theta_\mathrm{ref}(t)=\omega t$ & {\tt rotation} & {\sf x} & {\sf x} \\
954 $\theta_\mathrm{av}(t)$ & degrees & (\ref{eqn:avangle}) & {\tt rotation} & {\sf x} & - \\
955 $\theta_\mathrm{fit}(t)$, $\theta_\mathrm{fit}(t,n)$ & degrees & (\ref{eqn:rmsdfit}) & {\tt rotangles} & - & {\sf x} \\
956 $\ve{y}_0(n)$, $\ve{x}_0(t,n)$ & nm & (\ref{eqn:defx0}, \ref{eqn:defy0})& {\tt rotslabs} & - & {\sf x} \\
957 $\tau(t)$ & kJ/mol & (\ref{eqn:torque}) & {\tt rotation} & {\sf x} & - \\
958 $\tau(t,n)$ & kJ/mol & (\ref{eqn:torque}) & {\tt rottorque} & - & {\sf x} \\ \hline
959 \end{tabular}
960 \end{center}
961 \end{table}
964 \subsubsection*{Angle of Rotation Groups: Fixed Axis}
965 For fixed axis rotation, the average angle $\theta_\mathrm{av}(t)$ of the
966 group relative to the reference group is determined via the distance-weighted
967 angular deviation of all rotation group atoms from their reference positions,
968 \beq
969 \theta_\mathrm{av} = \left. \sum_{i=1}^{N} r_i \ \theta_i \right/ \sum_{i=1}^N r_i \ .
970 \label{eqn:avangle}
971 \eeq
972 Here, $r_i$ is the distance of the reference position to the rotation axis, and
973 the difference angles $\theta_i$ are determined from the atomic positions,
974 projected onto a plane perpendicular to the rotation axis through pivot point
975 $\ve{u}$ (see \eqnref{project} for the definition of $\perp$),
976 \beq
977 \cos \theta_i =
978 \frac{(\ve{y}_i-\ve{u})^\perp \cdot (\ve{x}_i-\ve{u})^\perp}
979 { \| (\ve{y}_i-\ve{u})^\perp \cdot (\ve{x}_i-\ve{u})^\perp
980 \| } \ .
981 \eeq
983 The sign of $\theta_\mathrm{av}$ is chosen such that
984 $\theta_\mathrm{av} > 0$ if the actual structure rotates ahead of the reference.
986 \subsubsection*{Angle of Rotation Groups: Flexible Axis}
987 For flexible axis rotation, two outputs are provided, the angle of the
988 entire rotation group, and separate angles for the segments in the slabs.
989 The angle of the entire rotation group is determined by an RMSD fit
990 of $\ve{x}_i$
991 to the reference positions $\ve{y}_i^0$ at $t=0$, yielding $\theta_\mathrm{fit}$
992 as the angle by which the reference has to be rotated around $\hat{\ve{v}}$
993 for the optimal fit,
994 \beq
995 \mathrm{RMSD} \big( \ve{x}_i,\ \mathbf{\Omega}(\theta_\mathrm{fit})
996 \ve{y}_i^0 \big) \stackrel{!}{=} \mathrm{min} \, .
997 \label{eqn:rmsdfit}
998 \eeq
999 To determine the local angle for each slab $n$, both reference and actual
1000 positions are weighted with the Gaussian function of slab $n$, and
1001 $\theta_\mathrm{fit}(t,n)$ is calculated as in \eqnref{rmsdfit}) from the
1002 Gaussian-weighted positions.
1004 For all angles, the {\tt .mdp} input option {\tt rot-fit-method} controls
1005 whether a normal RMSD fit is performed or whether for the fit each
1006 position $\ve{x}_i$ is put at the same distance to the rotation axis as its
1007 reference counterpart $\ve{y}_i^0$. In the latter case, the RMSD
1008 measures only angular differences, not radial ones.
1011 \subsubsection*{Angle Determination by Searching the Energy Minimum}
1012 Alternatively, for {\tt rot-fit-method = potential}, the angle of the rotation
1013 group is determined as the angle for which the rotation potential energy is minimal.
1014 Therefore, the used rotation potential is additionally evaluated for a set of angles
1015 around the current reference angle. In this case, the {\tt rotangles.log} output file
1016 contains the values of the rotation potential at the chosen set of angles, while
1017 {\tt rotation.xvg} lists the angle with minimal potential energy.
1020 \subsubsection*{Torque}
1021 \label{torque}
1022 The torque $\ve{\tau}(t)$ exerted by the rotation potential is calculated for fixed
1023 axis rotation via
1024 \beq
1025 \ve{\tau}(t) = \sum_{i=1}^{N} \ve{r}_i(t) \times \ve{f}_{\!i}^\perp(t) ,
1026 \label{eqn:torque}
1027 \eeq
1028 where $\ve{r}_i(t)$ is the distance vector from the rotation axis to
1029 $\ve{x}_i(t)$ and $\ve{f}_{\!i}^\perp(t)$ is the force component
1030 perpendicular to $\ve{r}_i(t)$ and $\hat{\ve{v}}$. For flexible axis
1031 rotation, torques $\ve{\tau}_{\!n}$ are calculated for each slab using the
1032 local rotation axis of the slab and the Gaussian-weighted positions.
1035 \section{\normindex{Computational Electrophysiology}}
1036 \label{sec:compel}
1038 The Computational Electrophysiology (CompEL) protocol \cite{Kutzner2011b} allows the simulation of
1039 ion flux through membrane channels, driven by transmembrane potentials or ion
1040 concentration gradients. Just as in real cells, CompEL establishes transmembrane
1041 potentials by sustaining a small imbalance of charges $\Delta q$ across the membrane,
1042 which gives rise to a potential difference $\Delta U$ according to the membrane capacitance:
1043 \beq
1044 \Delta U = \Delta q / C_{membrane}
1045 \eeq
1046 The transmembrane electric field and concentration gradients are controlled by
1047 {\tt .mdp} options, which allow the user to set reference counts for the ions on either side
1048 of the membrane. If a difference between the actual and the reference numbers persists
1049 over a certain time span, specified by the user, a number of ion/water pairs are
1050 exchanged between the compartments until the reference numbers are restored.
1051 Alongside the calculation of channel conductance and ion selectivity, CompEL simulations also
1052 enable determination of the channel reversal potential, an important
1053 characteristic obtained in electrophysiology experiments.
1055 In a CompEL setup, the simulation system is divided into two compartments {\bf A} and {\bf B}
1056 with independent ion concentrations. This is best achieved by using double bilayer systems with
1057 a copy (or copies) of the channel/pore of interest in each bilayer (\figref{compelsetup} A, B).
1058 If the channel axes point in the same direction, channel flux is observed
1059 simultaneously at positive and negative potentials in this way, which is for instance
1060 important for studying channel rectification.
1062 \begin{figure}
1063 \centerline{\includegraphics[width=13.5cm]{plots/compelsetup.pdf}}
1064 \caption{Typical double-membrane setup for CompEL simulations (A, B).
1065 Ion\,/\,water molecule exchanges will be performed as needed
1066 between the two light blue volumes around the dotted black lines (A).
1067 Plot (C) shows the potential difference $\Delta U$ resulting
1068 from the selected charge imbalance $\Delta q_{ref}$ between the compartments.}
1069 \label{fig:compelsetup}
1070 \end{figure}
1072 The potential difference $\Delta U$ across the membrane is easily calculated with the
1073 {\tt gmx potential} utility. By this, the potential drop along $z$ or the
1074 pore axis is exactly known in each time interval of the simulation (\figref{compelsetup} C).
1075 Type and number of ions $n_i$ of charge $q_i$, traversing the channel in the simulation,
1076 are written to the {\tt swapions.xvg} output file, from which the average channel
1077 conductance $G$ in each interval $\Delta t$ is determined by:
1078 \beq
1079 G = \frac{\sum_{i} n_{i}q_{i}}{\Delta t \, \Delta U} \, .
1080 \eeq
1081 The ion selectivity is calculated as the number flux ratio of different species.
1082 Best results are obtained by averaging these values over several overlapping time intervals.
1084 The calculation of reversal potentials is best achieved using a small set of simulations in which a given
1085 transmembrane concentration gradient is complemented with small ion imbalances of varying magnitude. For
1086 example, if one compartment contains 1\,M salt and the other 0.1\,M, and given charge neutrality otherwise,
1087 a set of simulations with $\Delta q = 0\,e$, $\Delta q = 2\,e$, $\Delta q = 4\,e$ could
1088 be used. Fitting a straight line through the current-voltage relationship of all obtained
1089 $I$-$U$ pairs near zero current will then yield $U_{rev}$.
1091 \subsection{Usage}
1092 The following {\tt .mdp} options control the CompEL protocol:
1093 {\small
1094 \begin{verbatim}
1095 swapcoords = Z ; Swap positions: no, X, Y, Z
1096 swap-frequency = 100 ; Swap attempt frequency
1097 \end{verbatim}}
1098 Choose {\tt Z} if your membrane is in the $xy$-plane (\figref{compelsetup}).
1099 Ions will be exchanged between compartments depending on their $z$-positions alone.
1100 {\tt swap-frequency} determines how often a swap attempt will be made.
1101 This step requires that the positions of the split groups, the ions, and possibly the solvent molecules are
1102 communicated between the parallel processes, so if chosen too small it can decrease the simulation
1103 performance. The {\tt Position swapping} entry in the cycle and time accounting
1104 table at the end of the {\tt md.log} file summarizes the amount of runtime spent
1105 in the swap module.
1106 {\small
1107 \begin{verbatim}
1108 split-group0 = channel0 ; Defines compartment boundary
1109 split-group1 = channel1 ; Defines other compartment boundary
1110 massw-split0 = no ; use mass-weighted center?
1111 massw-split1 = no
1112 \end{verbatim}}
1113 {\tt split-group0} and {\tt split-group1} are two index groups that define the boundaries
1114 between the two compartments, which are usually the centers of the channels.
1115 If {\tt massw-split0} or {\tt massw-split1} are set to {\tt yes}, the center of mass
1116 of each index group is used as boundary, here in $z$-direction. Otherwise, the
1117 geometrical centers will be used ($\times$ in \figref{compelsetup} A). If, such as here, a membrane
1118 channel is selected as split group, the center of the channel will define the dividing
1119 plane between the compartments (dashed horizontal lines). All index groups
1120 must be defined in the index file.
1122 If, to restore the requested ion counts, an ion from one compartment has to be exchanged
1123 with a water molecule from the other compartment, then those molecules are swapped
1124 which have the largest distance to the compartment-defining boundaries
1125 (dashed horizontal lines). Depending on the ion concentration,
1126 this effectively results in exchanges of molecules between the light blue volumes.
1127 If a channel is very asymmetric in $z$-direction and would extend into one of the
1128 swap volumes, one can offset the swap exchange plane with the {\tt bulk-offset}
1129 parameter. A value of 0.0 means no offset $b$, values $-1.0 < b < 0$ move the
1130 swap exchange plane closer to the lower, values $0 < b < 1.0$ closer to the upper
1131 membrane. \figref{compelsetup} A (left) depicts that for the {\bf A} compartment.
1133 {\small
1134 \begin{verbatim}
1135 solvent-group = SOL ; Group containing the solvent molecules
1136 iontypes = 3 ; Number of different ion types to control
1137 iontype0-name = NA ; Group name of the ion type
1138 iontype0-in-A = 51 ; Reference count of ions of type 0 in A
1139 iontype0-in-B = 35 ; Reference count of ions of type 0 in B
1140 iontype1-name = K
1141 iontype1-in-A = 10
1142 iontype1-in-B = 38
1143 iontype2-name = CL
1144 iontype2-in-A = -1
1145 iontype2-in-B = -1
1146 \end{verbatim}}
1147 The group name of solvent molecules acting as exchange partners for the ions
1148 has to be set with {\tt solvent-group}. The number of different ionic species under
1149 control of the CompEL protocol is given by the {\tt iontypes} parameter,
1150 while {\tt iontype0-name} gives the name of the index group containing the
1151 atoms of this ionic species. The reference number of ions of this type
1152 can be set with the {\tt iontype0-in-A} and {\tt iontype0-in-B} options
1153 for compartments {\bf A} and {\bf B}, respectively. Obviously,
1154 the sum of {\tt iontype0-in-A} and {\tt iontype0-in-B} needs to equal the number
1155 of ions in the group defined by {\tt iontype0-name}.
1156 A reference number of {\tt -1} means: use the number of ions as found at the beginning
1157 of the simulation as the reference value.
1159 {\small
1160 \begin{verbatim}
1161 coupl-steps = 10 ; Average over these many swap steps
1162 threshold = 1 ; Do not swap if < threshold
1163 \end{verbatim}}
1164 If {\tt coupl-steps} is set to 1, then the momentary ion distribution determines
1165 whether ions are exchanged. {\tt coupl-steps} \textgreater\ 1 will use the time-average
1166 of ion distributions over the selected number of attempt steps instead. This can be useful, for example,
1167 when ions diffuse near compartment boundaries, which would lead to numerous unproductive
1168 ion exchanges. A {\tt threshold} of 1 means that a swap is performed if the average ion
1169 count in a compartment differs by at least 1 from the requested values. Higher thresholds
1170 will lead to toleration of larger differences. Ions are exchanged until the requested
1171 number $\pm$ the threshold is reached.
1173 {\small
1174 \begin{verbatim}
1175 cyl0-r = 5.0 ; Split cylinder 0 radius (nm)
1176 cyl0-up = 0.75 ; Split cylinder 0 upper extension (nm)
1177 cyl0-down = 0.75 ; Split cylinder 0 lower extension (nm)
1178 cyl1-r = 5.0 ; same for other channel
1179 cyl1-up = 0.75
1180 cyl1-down = 0.75
1181 \end{verbatim}}
1182 The cylinder options are used to define virtual geometric cylinders around the
1183 channel's pore to track how many ions of which type have passed each channel.
1184 Ions will be counted as having traveled through a channel
1185 according to the definition of the channel's cylinder radius, upper and lower extension,
1186 relative to the location of the respective split group. This will not affect the actual
1187 flux or exchange, but will provide you with the ion permeation numbers across
1188 each of the channels. Note that an ion can only be counted as passing through a particular
1189 channel if it is detected \emph{within} the defined split cylinder in a swap step.
1190 If {\tt swap-frequency} is chosen too high, a particular ion may be detected in compartment {\bf A}
1191 in one swap step, and in compartment {\bf B} in the following swap step, so it will be unclear
1192 through which of the channels it has passed.
1194 A double-layered system for CompEL simulations can be easily prepared by
1195 duplicating an existing membrane/channel MD system in the direction of the membrane
1196 normal (typically $z$) with {\tt gmx editconf -translate 0 0 <l_z>}, where {\tt l_z}
1197 is the box length in that direction. If you have already defined index groups for
1198 the channel for the single-layered system, {\tt gmx make_ndx -n index.ndx -twin} will
1199 provide you with the groups for the double-layered system.
1201 To suppress large fluctuations of the membranes along the swap direction,
1202 it may be useful to apply a harmonic potential (acting only in the swap dimension)
1203 between each of the two channel and/or bilayer centers using umbrella pulling
1204 (see section~\ref{sec:pull}).
1206 \subsection*{Multimeric channels}
1207 If a split group consists of more than one molecule, the correct PBC image of all molecules
1208 with respect to each other has to be chosen such that the channel center can be correctly
1209 determined. \gromacs\ assumes that the starting structure in the {\tt .tpr}
1210 file has the correct PBC representation. Set the following environment variable
1211 to check whether that is the case:
1212 \begin{itemize}
1213 \item {\tt GMX_COMPELDUMP}: output the starting structure after it has been made whole to
1214 {\tt .pdb} file.
1215 \end{itemize}
1218 \section{Calculating a PMF using the free-energy code}
1219 \label{sec:fepmf}
1220 \index{potentials of mean force}
1221 \index{free energy calculations}
1222 The free-energy coupling-parameter approach (see~\secref{fecalc})
1223 provides several ways to calculate potentials of mean force.
1224 A potential of mean force between two atoms can be calculated
1225 by connecting them with a harmonic potential or a constraint.
1226 For this purpose there are special potentials that avoid the generation of
1227 extra exclusions, see~\secref{excl}.
1228 When the position of the minimum or the constraint length is 1 nm more
1229 in state B than in state A, the restraint or constraint force is given
1230 by $\partial H/\partial \lambda$.
1231 The distance between the atoms can be changed as a function of $\lambda$
1232 and time by setting {\tt delta-lambda} in the {\tt .mdp} file.
1233 The results should be identical (although not numerically
1234 due to the different implementations) to the results of the pull code
1235 with umbrella sampling and constraint pulling.
1236 Unlike the pull code, the free energy code can also handle atoms that
1237 are connected by constraints.
1239 Potentials of mean force can also be calculated using position restraints.
1240 With position restraints, atoms can be linked to a position in space
1241 with a harmonic potential (see \ssecref{positionrestraint}).
1242 These positions can be made a function of the coupling parameter $\lambda$.
1243 The positions for the A and the B states are supplied to {\tt grompp} with
1244 the {\tt -r} and {\tt -rb} options, respectively.
1245 One could use this approach to do \normindex{targeted MD};
1246 note that we do not encourage the use of targeted MD for proteins.
1247 A protein can be forced from one conformation to another by using
1248 these conformations as position restraint coordinates for state A and B.
1249 One can then slowly change $\lambda$ from 0 to 1.
1250 The main drawback of this approach is that the conformational freedom
1251 of the protein is severely limited by the position restraints,
1252 independent of the change from state A to B.
1253 Also, the protein is forced from state A to B in an almost straight line,
1254 whereas the real pathway might be very different.
1255 An example of a more fruitful application is a solid system or a liquid
1256 confined between walls where one wants to measure the force required
1257 to change the separation between the boundaries or walls.
1258 Because the boundaries (or walls) already need to be fixed,
1259 the position restraints do not limit the system in its sampling.
1261 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1262 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1263 %%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1264 \newcommand{\amine}{\sf -NH$_2$}
1265 \newcommand{\amines}{\sf -NH-}
1266 \newcommand{\aminep}{\sf -NH$_3^+$}
1267 \section{Removing fastest \swapindex{degrees of}{freedom}}
1268 \label{sec:rmfast}
1269 The maximum time step in MD simulations is limited by the smallest
1270 oscillation period that can be found in the simulated
1271 system. Bond-stretching vibrations are in their quantum-mechanical
1272 ground state and are therefore better represented by a constraint
1273 instead of a harmonic potential.
1275 For the remaining degrees of freedom, the shortest oscillation period
1276 (as measured from a simulation) is 13~fs for bond-angle vibrations
1277 involving hydrogen atoms. Taking as a guideline that with a Verlet
1278 (leap-frog) integration scheme a minimum of 5 numerical integration
1279 steps should be performed per period of a harmonic oscillation in
1280 order to integrate it with reasonable accuracy, the maximum time step
1281 will be about 3~fs. Disregarding these very fast oscillations of
1282 period 13~fs, the next shortest periods are around 20~fs, which will
1283 allow a maximum time step of about 4~fs.
1285 Removing the bond-angle degrees of freedom from hydrogen atoms can
1286 best be done by defining them as \normindex{virtual interaction sites}
1287 instead of normal atoms. Whereas a normal atom is connected to the molecule
1288 with bonds, angles and dihedrals, a virtual site's position is calculated
1289 from the position of three nearby heavy atoms in a predefined manner
1290 (see also \secref{virtual_sites}). For the hydrogens in water and in
1291 hydroxyl, sulfhydryl, or amine groups, no degrees of freedom can be
1292 removed, because rotational freedom should be preserved. The only
1293 other option available to slow down these motions is to increase the
1294 mass of the hydrogen atoms at the expense of the mass of the connected
1295 heavy atom. This will increase the moment of inertia of the water
1296 molecules and the hydroxyl, sulfhydryl, or amine groups, without
1297 affecting the equilibrium properties of the system and without
1298 affecting the dynamical properties too much. These constructions will
1299 shortly be described in \secref{vsitehydro} and have previously
1300 been described in full detail~\cite{feenstra99}.
1302 Using both virtual sites and \swapindex{modified}{mass}es, the next
1303 bottleneck is likely to be formed by the improper dihedrals (which are
1304 used to preserve planarity or chirality of molecular groups) and the
1305 peptide dihedrals. The peptide dihedral cannot be changed without
1306 affecting the physical behavior of the protein. The improper dihedrals
1307 that preserve planarity mostly deal with aromatic residues. Bonds,
1308 angles, and dihedrals in these residues can also be replaced with
1309 somewhat elaborate virtual site constructions.
1311 All modifications described in this section can be performed using the
1312 {\gromacs} topology building tool {\tt \normindex{pdb2gmx}}. Separate
1313 options exist to increase hydrogen masses, virtualize all hydrogen atoms,
1314 or also virtualize all aromatic residues. {\bf Note} that when all hydrogen
1315 atoms are virtualized, those inside the aromatic residues will be
1316 virtualized as well, {\ie} hydrogens in the aromatic residues are treated
1317 differently depending on the treatment of the aromatic residues.
1319 Parameters for the virtual site constructions for the hydrogen atoms are
1320 inferred from the force-field parameters ({\em vis}. bond lengths and
1321 angles) directly by {\tt \normindex{grompp}} while processing the
1322 topology file. The constructions for the aromatic residues are based
1323 on the bond lengths and angles for the geometry as described in the
1324 force fields, but these parameters are hard-coded into {\tt
1325 \normindex{pdb2gmx}} due to the complex nature of the construction
1326 needed for a whole aromatic group.
1328 \subsection{Hydrogen bond-angle vibrations}
1329 \label{sec:vsitehydro}
1330 \subsubsection{Construction of virtual sites} %%%%%%%%%%%%%%%%%%%%%%%%%
1331 \begin{figure}
1332 \centerline{\includegraphics[width=11cm]{plots/dumtypes}}
1333 \caption[Virtual site constructions for hydrogen atoms.]{The different
1334 types of virtual site constructions used for hydrogen atoms. The atoms
1335 used in the construction of the virtual site(s) are depicted as black
1336 circles, virtual sites as gray ones. Hydrogens are smaller than heavy
1337 atoms. {\sf A}: fixed bond angle, note that here the hydrogen is not a
1338 virtual site; {\sf B}: in the plane of three atoms, with fixed distance;
1339 {\sf C}: in the plane of three atoms, with fixed angle and distance;
1340 {\sf D}: construction for amine groups ({\amine} or {\aminep}), see
1341 text for details.}
1342 \label{fig:vsitehydro}
1343 \end{figure}
1345 The goal of defining hydrogen atoms as virtual sites is to remove all
1346 high-frequency degrees of freedom from them. In some cases, not all
1347 degrees of freedom of a hydrogen atom should be removed, {\eg} in the
1348 case of hydroxyl or amine groups the rotational freedom of the
1349 hydrogen atom(s) should be preserved. Care should be taken that no
1350 unwanted correlations are introduced by the construction of virtual
1351 sites, {\eg} bond-angle vibration between the constructing atoms could
1352 translate into hydrogen bond-length vibration. Additionally, since
1353 virtual sites are by definition massless, in order to preserve total
1354 system mass, the mass of each hydrogen atom that is treated as virtual
1355 site should be added to the bonded heavy atom.
1357 Taking into account these considerations, the hydrogen atoms in a
1358 protein naturally fall into several categories, each requiring a
1359 different approach (see also \figref{vsitehydro}).
1361 \begin{itemize}
1363 \item{\em hydroxyl ({\sf -OH}) or sulfhydryl ({\sf -SH})
1364 hydrogen:\/} The only internal degree of freedom in a hydroxyl group
1365 that can be constrained is the bending of the {\sf C-O-H} angle. This
1366 angle is fixed by defining an additional bond of appropriate length,
1367 see \figref{vsitehydro}A. Doing so removes the high-frequency angle bending,
1368 but leaves the dihedral rotational freedom. The same goes for a
1369 sulfhydryl group. {\bf Note} that in these cases the hydrogen is not treated
1370 as a virtual site.
1372 \item{\em single amine or amide ({\amines}) and aromatic hydrogens
1373 ({\sf -CH-}):\/} The position of these hydrogens cannot be constructed
1374 from a linear combination of bond vectors, because of the flexibility
1375 of the angle between the heavy atoms. Instead, the hydrogen atom is
1376 positioned at a fixed distance from the bonded heavy atom on a line
1377 going through the bonded heavy atom and a point on the line through
1378 both second bonded atoms, see \figref{vsitehydro}B.
1380 \item{\em planar amine ({\amine}) hydrogens:\/} The method used for
1381 the single amide hydrogen is not well suited for planar amine groups,
1382 because no suitable two heavy atoms can be found to define the
1383 direction of the hydrogen atoms. Instead, the hydrogen is constructed
1384 at a fixed distance from the nitrogen atom, with a fixed angle to the
1385 carbon atom, in the plane defined by one of the other heavy atoms, see
1386 \figref{vsitehydro}C.
1388 \item{\em amine group (umbrella {\amine} or {\aminep}) hydrogens:\/}
1389 Amine hydrogens with rotational freedom cannot be constructed as virtual
1390 sites from the heavy atoms they are connected to, since this would
1391 result in loss of the rotational freedom of the amine group. To
1392 preserve the rotational freedom while removing the hydrogen bond-angle
1393 degrees of freedom, two ``dummy masses'' are constructed with the same
1394 total mass, moment of inertia (for rotation around the {\sf C-N} bond)
1395 and center of mass as the amine group. These dummy masses have no
1396 interaction with any other atom, except for the fact that they are
1397 connected to the carbon and to each other, resulting in a rigid
1398 triangle. From these three particles, the positions of the nitrogen and
1399 hydrogen atoms are constructed as linear combinations of the two
1400 carbon-mass vectors and their outer product, resulting in an amine
1401 group with rotational freedom intact, but without other internal
1402 degrees of freedom. See \figref{vsitehydro}D.
1404 \end{itemize}
1406 \begin{figure}
1407 \centerline{\includegraphics[width=15cm]{plots/dumaro}}
1408 \caption[Virtual site constructions for aromatic residues.]{The
1409 different types of virtual site constructions used for aromatic
1410 residues. The atoms used in the construction of the virtual site(s) are
1411 depicted as black circles, virtual sites as gray ones. Hydrogens are
1412 smaller than heavy atoms. {\sf A}: phenylalanine; {\sf B}: tyrosine
1413 (note that the hydroxyl hydrogen is {\em not} a virtual site); {\sf C}:
1414 tryptophan; {\sf D}: histidine.}
1415 \label{fig:vistearo}
1416 \end{figure}
1418 \subsection{Out-of-plane vibrations in aromatic groups}
1419 \label{sec:vsitearo}
1420 The planar arrangements in the side chains of the aromatic residues
1421 lends itself perfectly to a virtual-site construction, giving a
1422 perfectly planar group without the inherently unstable constraints
1423 that are necessary to keep normal atoms in a plane. The basic approach
1424 is to define three atoms or dummy masses with constraints between them
1425 to fix the geometry and create the rest of the atoms as simple virtual
1426 sites type (see \secref{virtual_sites}) from these three. Each of
1427 the aromatic residues require a different approach:
1429 \begin{itemize}
1431 \item{\em Phenylalanine:\/} {\sf C}$_\gamma$, {\sf C}$_{{\epsilon}1}$,
1432 and {\sf C}$_{{\epsilon}2}$ are kept as normal atoms, but with each a
1433 mass of one third the total mass of the phenyl group. See
1434 \figref{vsitehydro}A.
1436 \item{\em Tyrosine:\/} The ring is treated identically to the
1437 phenylalanine ring. Additionally, constraints are defined between {\sf
1438 C}$_{{\epsilon}1}$, {\sf C}$_{{\epsilon}2}$, and {\sf O}$_{\eta}$.
1439 The original improper dihedral angles will keep both triangles (one
1440 for the ring and one with {\sf O}$_{\eta}$) in a plane, but due to the
1441 larger moments of inertia this construction will be much more
1442 stable. The bond-angle in the hydroxyl group will be constrained by a
1443 constraint between {\sf C}$_\gamma$ and {\sf H}$_{\eta}$. {\bf Note} that
1444 the hydrogen is not treated as a virtual site. See
1445 \figref{vsitehydro}B.
1447 \item{\em Tryptophan:\/} {\sf C}$_\beta$ is kept as a normal atom
1448 and two dummy masses are created at the center of mass of each of the
1449 rings, each with a mass equal to the total mass of the respective ring
1450 ({\sf C}$_{{\delta}2}$ and {\sf C}$_{{\epsilon}2}$ are each
1451 counted half for each ring). This keeps the overall center of mass and
1452 the moment of inertia almost (but not quite) equal to what it was. See
1453 \figref{vsitehydro}C.
1455 \item{\em Histidine:\/} {\sf C}$_\gamma$, {\sf C}$_{{\epsilon}1}$
1456 and {\sf N}$_{{\epsilon}2}$ are kept as normal atoms, but with masses
1457 redistributed such that the center of mass of the ring is
1458 preserved. See \figref{vsitehydro}D.
1460 \end{itemize}
1462 \section{Viscosity calculation\index{viscosity}}
1464 The shear viscosity is a property of liquids that can be determined easily
1465 by experiment. It is useful for parameterizing a force field
1466 because it is a kinetic property, while most other properties
1467 which are used for parameterization are thermodynamic.
1468 The viscosity is also an important property, since it influences
1469 the rates of conformational changes of molecules solvated in the liquid.
1471 The viscosity can be calculated from an equilibrium simulation using
1472 an Einstein relation:
1473 \beq
1474 \eta = \frac{1}{2}\frac{V}{k_B T} \lim_{t \rightarrow \infty}
1475 \frac{\mbox{d}}{\mbox{d} t} \left\langle
1476 \left( \int_{t_0}^{{t_0}+t} P_{xz}(t') \mbox{d} t' \right)^2
1477 \right\rangle_{t_0}
1478 \eeq
1479 This can be done with {\tt g_energy}.
1480 This method converges very slowly~\cite{Hess2002a}, and as such
1481 a nanosecond simulation might not
1482 be long enough for an accurate determination of the viscosity.
1483 The result is very dependent on the treatment of the electrostatics.
1484 Using a (short) cut-off results in large noise on the off-diagonal
1485 pressure elements, which can increase the calculated viscosity by an order
1486 of magnitude.
1488 {\gromacs} also has a non-equilibrium method for determining
1489 the viscosity~\cite{Hess2002a}.
1490 This makes use of the fact that energy, which is fed into system by
1491 external forces, is dissipated through viscous friction. The generated heat
1492 is removed by coupling to a heat bath. For a Newtonian liquid adding a
1493 small force will result in a velocity gradient according to the following
1494 equation:
1495 \beq
1496 a_x(z) + \frac{\eta}{\rho} \frac{\partial^2 v_x(z)}{\partial z^2} = 0
1497 \eeq
1498 Here we have applied an acceleration $a_x(z)$ in the $x$-direction, which
1499 is a function of the $z$-coordinate.
1500 In {\gromacs} the acceleration profile is:
1501 \beq
1502 a_x(z) = A \cos\left(\frac{2\pi z}{l_z}\right)
1503 \eeq
1504 where $l_z$ is the height of the box. The generated velocity profile is:
1505 \beq
1506 v_x(z) = V \cos\left(\frac{2\pi z}{l_z}\right)
1507 \eeq
1508 \beq
1509 V = A \frac{\rho}{\eta}\left(\frac{l_z}{2\pi}\right)^2
1510 \eeq
1511 The viscosity can be calculated from $A$ and $V$:
1512 \beq
1513 \label{visc}
1514 \eta = \frac{A}{V}\rho \left(\frac{l_z}{2\pi}\right)^2
1515 \eeq
1517 In the simulation $V$ is defined as:
1518 \beq
1519 V = \frac{\displaystyle \sum_{i=1}^N m_i v_{i,x} 2 \cos\left(\frac{2\pi z}{l_z}\right)}
1520 {\displaystyle \sum_{i=1}^N m_i}
1521 \eeq
1522 The generated velocity profile is not coupled to the heat bath. Moreover,
1523 the velocity profile is excluded from the kinetic energy.
1524 One would like $V$ to be as large as possible to get good statistics.
1525 However, the shear rate should not be so high that the system gets too far
1526 from equilibrium. The maximum shear rate occurs where the cosine is zero,
1527 the rate being:
1528 \beq
1529 \mbox{sh}_{\max} = \max_z \left| \frac{\partial v_x(z)}{\partial z} \right|
1530 = A \frac{\rho}{\eta} \frac{l_z}{2\pi}
1531 \eeq
1532 For a simulation with: $\eta=10^{-3}$ [kg\,m$^{-1}$\,s$^{-1}$],
1533 $\rho=10^3$\,[kg\,m$^{-3}$] and $l_z=2\pi$\,[nm],
1534 $\mbox{sh}_{\max}=1$\,[ps\,nm$^{-1}$] $A$.
1535 This shear rate should be smaller than one over the longest
1536 correlation time in the system. For most liquids, this will be the rotation
1537 correlation time, which is around 10 ps. In this case, $A$ should
1538 be smaller than 0.1\,[nm\,ps$^{-2}$].
1539 When the shear rate is too high, the observed viscosity will be too low.
1540 Because $V$ is proportional to the square of the box height,
1541 the optimal box is elongated in the $z$-direction.
1542 In general, a simulation length of 100 ps is enough to obtain an
1543 accurate value for the viscosity.
1545 The heat generated by the viscous friction is removed by coupling to a heat
1546 bath. Because this coupling is not instantaneous the real temperature of the
1547 liquid will be slightly lower than the observed temperature.
1548 Berendsen derived this temperature shift~\cite{Berendsen91}, which can
1549 be written in terms of the shear rate as:
1550 \beq
1551 T_s = \frac{\eta\,\tau}{2 \rho\,C_v} \mbox{sh}_{\max}^2
1552 \eeq
1553 where $\tau$ is the coupling time for the Berendsen thermostat and
1554 $C_v$ is the heat capacity. Using the values of the example above,
1555 $\tau=10^{-13}$ [s] and $C_v=2 \cdot 10^3$\,[J kg$^{-1}$\,K$^{-1}$], we
1556 get: $T_s=25$\,[K\,ps$^{-2}$]\,sh$_{\max}^2$. When we want the shear
1557 rate to be smaller than $1/10$\,[ps$^{-1}$], $T_s$ is smaller than
1558 0.25\,[K], which is negligible.
1560 {\bf Note} that the system has to build up the velocity profile when starting
1561 from an equilibrium state. This build-up time is of the order of the
1562 correlation time of the liquid.
1564 Two quantities are written to the energy file, along with their averages
1565 and fluctuations: $V$ and $1/\eta$, as obtained from (\ref{visc}).
1567 \section{Tabulated interaction functions\index{tabulated interaction functions}}
1568 \subsection{Cubic splines for potentials}
1569 \label{subsec:cubicspline}
1570 In some of the inner loops of {\gromacs}, look-up tables are used
1571 for computation of potential and forces.
1572 The tables are interpolated using a cubic
1573 spline algorithm.
1574 There are separate tables for electrostatic, dispersion, and repulsion
1575 interactions,
1576 but for the sake of caching performance these have been combined
1577 into a single array.
1578 The cubic spline interpolation for $x_i \leq x < x_{i+1}$ looks like this:
1579 \beq
1580 V_s(x) = A_0 + A_1 \,\epsilon + A_2 \,\epsilon^2 + A_3 \,\epsilon^3
1581 \label{eqn:spline}
1582 \eeq
1583 where the table spacing $h$ and fraction $\epsilon$ are given by:
1584 \bea
1585 h &=& x_{i+1} - x_i \\
1586 \epsilon&=& (x - x_i)/h
1587 \eea
1588 so that $0 \le \epsilon < 1$.
1589 From this, we can calculate the derivative in order to determine the forces:
1590 \beq
1591 -V_s'(x) ~=~
1592 -\frac{{\rm d}V_s(x)}{{\rm d}\epsilon}\frac{{\rm d}\epsilon}{{\rm d}x} ~=~
1593 -(A_1 + 2 A_2 \,\epsilon + 3 A_3 \,\epsilon^2)/h
1594 \eeq
1595 The four coefficients are determined from the four conditions
1596 that $V_s$ and $-V_s'$ at both ends of each interval should match
1597 the exact potential $V$ and force $-V'$.
1598 This results in the following errors for each interval:
1599 \bea
1600 |V_s - V |_{max} &=& V'''' \frac{h^4}{384} + O(h^5) \\
1601 |V_s' - V' |_{max} &=& V'''' \frac{h^3}{72\sqrt{3}} + O(h^4) \\
1602 |V_s''- V''|_{max} &=& V'''' \frac{h^2}{12} + O(h^3)
1603 \eea
1604 V and V' are continuous, while V'' is the first discontinuous
1605 derivative.
1606 The number of points per nanometer is 500 and 2000
1607 for mixed- and double-precision versions of {\gromacs}, respectively.
1608 This means that the errors in the potential and force will usually
1609 be smaller than the mixed precision accuracy.
1611 {\gromacs} stores $A_0$, $A_1$, $A_2$ and $A_3$.
1612 The force routines get a table with these four parameters and
1613 a scaling factor $s$ that is equal to the number of points per nm.
1614 ({\bf Note} that $h$ is $s^{-1}$).
1615 The algorithm goes a little something like this:
1616 \begin{enumerate}
1617 \item Calculate distance vector (\ve{r}$_{ij}$) and distance r$_{ij}$
1618 \item Multiply r$_{ij}$ by $s$ and truncate to an integer value $n_0$
1619 to get a table index
1620 \item Calculate fractional component ($\epsilon$ = $s$r$_{ij} - n_0$)
1621 and $\epsilon^2$
1622 \item Do the interpolation to calculate the potential $V$ and the scalar force $f$
1623 \item Calculate the vector force \ve{F} by multiplying $f$ with \ve{r}$_{ij}$
1624 \end{enumerate}
1626 {\bf Note} that table look-up is significantly {\em
1627 slower} than computation of the most simple Lennard-Jones and Coulomb
1628 interaction. However, it is much faster than the shifted Coulomb
1629 function used in conjunction with the PPPM method. Finally, it is much
1630 easier to modify a table for the potential (and get a graphical
1631 representation of it) than to modify the inner loops of the MD
1632 program.
1634 \subsection{User-specified potential functions}
1635 \label{subsec:userpot}
1636 You can also use your own potential functions\index{potential function} without
1637 editing the {\gromacs} code. The potential function should be according to the
1638 following equation
1639 \beq
1640 V(r_{ij}) ~=~ \frac{q_i q_j}{4 \pi\epsilon_0} f(r_{ij}) + C_6 \,g(r_{ij}) + C_{12} \,h(r_{ij})
1641 \eeq
1642 where $f$, $g$, and $h$ are user defined functions. {\bf Note} that if $g(r)$ represents a
1643 normal dispersion interaction, $g(r)$ should be $<$ 0. C$_6$, C$_{12}$
1644 and the charges are read from the topology. Also note that combination
1645 rules are only supported for Lennard-Jones and Buckingham, and that
1646 your tables should match the parameters in the binary topology.
1648 When you add the following lines in your {\tt .mdp} file:
1650 {\small
1651 \begin{verbatim}
1652 rlist = 1.0
1653 coulombtype = User
1654 rcoulomb = 1.0
1655 vdwtype = User
1656 rvdw = 1.0
1657 \end{verbatim}}
1659 {\tt mdrun} will read a single non-bonded table file,
1660 or multiple when {\tt energygrp-table} is set (see below).
1661 The name of the file(s) can be set with the {\tt mdrun} option {\tt -table}.
1662 The table file should contain seven columns of table look-up data in the
1663 order: $x$, $f(x)$, $-f'(x)$, $g(x)$, $-g'(x)$, $h(x)$, $-h'(x)$.
1664 The $x$ should run from 0 to $r_c+1$ (the value of {\tt table_extension} can be
1665 changed in the {\tt .mdp} file).
1666 You can choose the spacing you like; for the standard tables {\gromacs}
1667 uses a spacing of 0.002 and 0.0005 nm when you run in mixed
1668 and double precision, respectively. In this
1669 context, $r_c$ denotes the maximum of the two cut-offs {\tt rvdw} and
1670 {\tt rcoulomb} (see above). These variables need not be the same (and
1671 need not be 1.0 either). Some functions used for potentials contain a
1672 singularity at $x = 0$, but since atoms are normally not closer to each
1673 other than 0.1 nm, the function value at $x = 0$ is not important.
1674 Finally, it is also
1675 possible to combine a standard Coulomb with a modified LJ potential
1676 (or vice versa). One then specifies {\eg} {\tt coulombtype = Cut-off} or
1677 {\tt coulombtype = PME}, combined with {\tt vdwtype = User}. The table file must
1678 always contain the 7 columns however, and meaningful data (i.e. not
1679 zeroes) must be entered in all columns. A number of pre-built table
1680 files can be found in the {\tt GMXLIB} directory for 6-8, 6-9, 6-10, 6-11, and 6-12
1681 Lennard-Jones potentials combined with a normal Coulomb.
1683 If you want to have different functional forms between different
1684 groups of atoms, this can be set through energy groups.
1685 Different tables can be used for non-bonded interactions between
1686 different energy groups pairs through the {\tt .mdp} option {\tt energygrp-table}
1687 (see details in the User Guide).
1688 Atoms that should interact with a different potential should
1689 be put into different energy groups.
1690 Between group pairs which are not listed in {\tt energygrp-table},
1691 the normal user tables will be used. This makes it easy to use
1692 a different functional form between a few types of atoms.
1694 \section{Mixed Quantum-Classical simulation techniques}
1696 In a molecular mechanics (MM) force field, the influence of electrons
1697 is expressed by empirical parameters that are assigned on the basis of
1698 experimental data, or on the basis of results from high-level quantum
1699 chemistry calculations. These are valid for the ground state of a
1700 given covalent structure, and the MM approximation is usually
1701 sufficiently accurate for ground-state processes in which the overall
1702 connectivity between the atoms in the system remains
1703 unchanged. However, for processes in which the connectivity does
1704 change, such as chemical reactions, or processes that involve multiple
1705 electronic states, such as photochemical conversions, electrons can no
1706 longer be ignored, and a quantum mechanical description is required
1707 for at least those parts of the system in which the reaction takes
1708 place.
1710 One approach to the simulation of chemical reactions in solution, or
1711 in enzymes, is to use a combination of quantum mechanics (QM) and
1712 molecular mechanics (MM). The reacting parts of the system are treated
1713 quantum mechanically, with the remainder being modeled using the
1714 force field. The current version of {\gromacs} provides interfaces to
1715 several popular Quantum Chemistry packages (MOPAC~\cite{mopac},
1716 GAMESS-UK~\cite{gamess-uk}, Gaussian~\cite{g03} and CPMD~\cite{Car85a}).
1718 {\gromacs} interactions between the two subsystems are
1719 either handled as described by Field {\em et al.}~\cite{Field90a} or
1720 within the ONIOM approach by Morokuma and coworkers~\cite{Maseras96a,
1721 Svensson96a}.
1723 \subsection{Overview}
1725 Two approaches for describing the interactions between the QM and MM
1726 subsystems are supported in this version:
1728 \begin{enumerate}
1729 \item{\textbf{Electronic Embedding}} The electrostatic interactions
1730 between the electrons of the QM region and the MM atoms and between
1731 the QM nuclei and the MM atoms are included in the Hamiltonian for
1732 the QM subsystem: \beq H^{QM/MM} =
1733 H^{QM}_e-\sum_i^n\sum_J^M\frac{e^2Q_J}{4\pi\epsilon_0r_{iJ}}+\sum_A^N\sum_J^M\frac{e^2Z_AQ_J}{e\pi\epsilon_0R_{AJ}},
1734 \eeq where $n$ and $N$ are the number of electrons and nuclei in the
1735 QM region, respectively, and $M$ is the number of charged MM
1736 atoms. The first term on the right hand side is the original
1737 electronic Hamiltonian of an isolated QM system. The first of the
1738 double sums is the total electrostatic interaction between the QM
1739 electrons and the MM atoms. The total electrostatic interaction of the
1740 QM nuclei with the MM atoms is given by the second double sum. Bonded
1741 interactions between QM and MM atoms are described at the MM level by
1742 the appropriate force-field terms. Chemical bonds that connect the two
1743 subsystems are capped by a hydrogen atom to complete the valence of
1744 the QM region. The force on this atom, which is present in the QM
1745 region only, is distributed over the two atoms of the bond. The cap
1746 atom is usually referred to as a link atom.
1748 \item{\textbf{ONIOM}} In the ONIOM approach, the energy and gradients
1749 are first evaluated for the isolated QM subsystem at the desired level
1750 of {\it{ab initio}} theory. Subsequently, the energy and gradients of
1751 the total system, including the QM region, are computed using the
1752 molecular mechanics force field and added to the energy and gradients
1753 calculated for the isolated QM subsystem. Finally, in order to correct
1754 for counting the interactions inside the QM region twice, a molecular
1755 mechanics calculation is performed on the isolated QM subsystem and
1756 the energy and gradients are subtracted. This leads to the following
1757 expression for the total QM/MM energy (and gradients likewise): \beq
1758 E_{tot} = E_{I}^{QM}
1759 +E_{I+II}^{MM}-E_{I}^{MM}, \eeq where the
1760 subscripts I and II refer to the QM and MM subsystems,
1761 respectively. The superscripts indicate at what level of theory the
1762 energies are computed. The ONIOM scheme has the
1763 advantage that it is not restricted to a two-layer QM/MM description,
1764 but can easily handle more than two layers, with each layer described
1765 at a different level of theory.
1766 \end{enumerate}
1768 \subsection{Usage}
1770 To make use of the QM/MM functionality in {\gromacs}, one needs to:
1772 \begin{enumerate}
1773 \item introduce link atoms at the QM/MM boundary, if needed;
1774 \item specify which atoms are to be treated at a QM level;
1775 \item specify the QM level, basis set, type of QM/MM interface and so on.
1776 \end{enumerate}
1778 \subsubsection{Adding link atoms}
1780 At the bond that connects the QM and MM subsystems, a link atoms is
1781 introduced. In {\gromacs} the link atom has special atomtype, called
1782 LA. This atomtype is treated as a hydrogen atom in the QM calculation,
1783 and as a virtual site in the force-field calculation. The link atoms, if
1784 any, are part of the system, but have no interaction with any other
1785 atom, except that the QM force working on it is distributed over the
1786 two atoms of the bond. In the topology, the link atom (LA), therefore,
1787 is defined as a virtual site atom:
1789 {\small
1790 \begin{verbatim}
1791 [ virtual_sites2 ]
1792 LA QMatom MMatom 1 0.65
1793 \end{verbatim}}
1795 See~\secref{vsitetop} for more details on how virtual sites are
1796 treated. The link atom is replaced at every step of the simulation.
1798 In addition, the bond itself is replaced by a constraint:
1800 {\small
1801 \begin{verbatim}
1802 [ constraints ]
1803 QMatom MMatom 2 0.153
1804 \end{verbatim}}
1806 {\bf Note} that, because in our system the QM/MM bond is a carbon-carbon
1807 bond (0.153 nm), we use a constraint length of 0.153 nm, and dummy
1808 position of 0.65. The latter is the ratio between the ideal C-H
1809 bond length and the ideal C-C bond length. With this ratio, the link
1810 atom is always 0.1 nm away from the {\tt QMatom}, consistent with the
1811 carbon-hydrogen bond length. If the QM and MM subsystems are connected
1812 by a different kind of bond, a different constraint and a different
1813 dummy position, appropriate for that bond type, are required.
1815 \subsubsection{Specifying the QM atoms}
1817 Atoms that should be treated at a QM level of theory, including the
1818 link atoms, are added to the index file. In addition, the chemical
1819 bonds between the atoms in the QM region are to be defined as
1820 connect bonds (bond type 5) in the topology file:
1822 {\small
1823 \begin{verbatim}
1824 [ bonds ]
1825 QMatom1 QMatom2 5
1826 QMatom2 QMatom3 5
1827 \end{verbatim}}
1829 \subsubsection{Specifying the QM/MM simulation parameters}
1831 In the {\tt .mdp} file, the following parameters control a QM/MM simulation.
1833 \begin{description}
1835 \item[\tt QMMM = no]\mbox{}\\ If this is set to {\tt yes}, a QM/MM
1836 simulation is requested. Several groups of atoms can be described at
1837 different QM levels separately. These are specified in the QMMM-grps
1838 field separated by spaces. The level of {\it{ab initio}} theory at which the
1839 groups are described is specified by {\tt QMmethod} and {\tt QMbasis}
1840 Fields. Describing the groups at different levels of theory is only
1841 possible with the ONIOM QM/MM scheme, specified by {\tt QMMMscheme}.
1843 \item[\tt QMMM-grps =]\mbox{}\\groups to be described at the QM level
1845 \item[\tt QMMMscheme = normal]\mbox{}\\Options are {\tt normal} and
1846 {\tt ONIOM}. This selects the QM/MM interface. {\tt normal} implies
1847 that the QM subsystem is electronically embedded in the MM
1848 subsystem. There can only be one {\tt QMMM-grps} that is modeled at
1849 the {\tt QMmethod} and {\tt QMbasis} level of {\it{ ab initio}}
1850 theory. The rest of the system is described at the MM level. The QM
1851 and MM subsystems interact as follows: MM point charges are included
1852 in the QM one-electron Hamiltonian and all Lennard-Jones interactions
1853 are described at the MM level. If {\tt ONIOM} is selected, the
1854 interaction between the subsystem is described using the ONIOM method
1855 by Morokuma and co-workers. There can be more than one QMMM-grps each
1856 modeled at a different level of QM theory (QMmethod and QMbasis).
1858 \item[\tt QMmethod = ]\mbox{}\\Method used to compute the energy
1859 and gradients on the QM atoms. Available methods are AM1, PM3, RHF,
1860 UHF, DFT, B3LYP, MP2, CASSCF, MMVB and CPMD. For CASSCF, the number of
1861 electrons and orbitals included in the active space is specified by
1862 {\tt CASelectrons} and {\tt CASorbitals}. For CPMD, the plane-wave
1863 cut-off is specified by the {\tt planewavecutoff} keyword.
1865 \item[\tt QMbasis = ]\mbox{}\\Gaussian basis set used to expand the
1866 electronic wave-function. Only Gaussian basis sets are currently
1867 available, i.e. STO-3G, 3-21G, 3-21G*, 3-21+G*, 6-21G, 6-31G, 6-31G*,
1868 6-31+G*, and 6-311G. For CPMD, which uses plane wave expansion rather
1869 than atom-centered basis functions, the {\tt planewavecutoff} keyword
1870 controls the plane wave expansion.
1872 \item[\tt QMcharge = ]\mbox{}\\The total charge in {\it{e}} of the {\tt
1873 QMMM-grps}. In case there are more than one {\tt QMMM-grps}, the total
1874 charge of each ONIOM layer needs to be specified separately.
1876 \item[\tt QMmult = ]\mbox{}\\The multiplicity of the {\tt
1877 QMMM-grps}. In case there are more than one {\tt QMMM-grps}, the
1878 multiplicity of each ONIOM layer needs to be specified separately.
1880 \item[\tt CASorbitals = ]\mbox{}\\The number of orbitals to be
1881 included in the active space when doing a CASSCF computation.
1883 \item[\tt CASelectrons = ]\mbox{}\\The number of electrons to be
1884 included in the active space when doing a CASSCF computation.
1886 \item[\tt SH = no]\mbox{}\\If this is set to yes, a QM/MM MD
1887 simulation on the excited state-potential energy surface and enforce a
1888 diabatic hop to the ground-state when the system hits the conical
1889 intersection hyperline in the course the simulation. This option only
1890 works in combination with the CASSCF method.
1892 \end{description}
1894 \subsection{Output}
1896 The energies and gradients computed in the QM calculation are added to
1897 those computed by {\gromacs}. In the {\tt .edr} file there is a section
1898 for the total QM energy.
1900 \subsection{Future developments}
1902 Several features are currently under development to increase the
1903 accuracy of the QM/MM interface. One useful feature is the use of
1904 delocalized MM charges in the QM computations. The most important
1905 benefit of using such smeared-out charges is that the Coulombic
1906 potential has a finite value at interatomic distances. In the point
1907 charge representation, the partially-charged MM atoms close to the QM
1908 region tend to ``over-polarize'' the QM system, which leads to artifacts
1909 in the calculation.
1911 What is needed as well is a transition state optimizer.
1913 \section{Using VMD plug-ins for trajectory file I/O}
1914 \index{VMD plug-ins}\index{trajectory file}{\gromacs} tools are able
1915 to use the plug-ins found in an existing installation of
1916 \href{http://www.ks.uiuc.edu/Research/vmd}{VMD} in order to read and
1917 write trajectory files in formats that are not native to
1918 {\gromacs}. You will be able to supply an AMBER DCD-format trajectory
1919 filename directly to {\gromacs} tools, for example.
1921 This requires a VMD installation not older than version 1.8, that your
1922 system provides the dlopen function so that programs can determine at
1923 run time what plug-ins exist, and that you build shared libraries when
1924 building {\gromacs}. CMake will find the vmd executable in your path, and
1925 from it, or the environment variable {\tt VMDDIR} at configuration or
1926 run time, locate the plug-ins. Alternatively, the {\tt VMD_PLUGIN_PATH}
1927 can be used at run time to specify a path where these plug-ins can be
1928 found. Note that these plug-ins are in a binary format, and that format
1929 must match the architecture of the machine attempting to use them.
1932 \section{\normindex{Interactive Molecular Dynamics}}
1933 {\gromacs} supports the interactive molecular dynamics (IMD) protocol as implemented
1934 by \href{http://www.ks.uiuc.edu/Research/vmd}{VMD} to control a running simulation
1935 in NAMD. IMD allows to monitor a running {\gromacs} simulation from a VMD client.
1936 In addition, the user can interact with the simulation by pulling on atoms, residues
1937 or fragments with a mouse or a force-feedback device. Additional information about
1938 the {\gromacs} implementation and an exemplary {\gromacs} IMD system can be found
1939 \href{http://www.mpibpc.mpg.de/grubmueller/interactivemd}{on this homepage}.
1941 \subsection{Simulation input preparation}
1942 The {\gromacs} implementation allows transmission and interaction with a part of the
1943 running simulation only, e.g.\ in cases where no water molecules should be transmitted
1944 or pulled. The group is specified via the {\tt .mdp} option {\tt IMD-group}. When
1945 {\tt IMD-group} is empty, the IMD protocol is disabled and cannot be enabled via the
1946 switches in {\tt mdrun}. To interact with the entire system, {\tt IMD-group} can
1947 be set to {\tt System}. When using {\tt grompp}, a {\tt .gro} file
1948 to be used as VMD input is written out ({\tt -imd} switch of {\tt grompp}).
1950 \subsection{Starting the simulation}
1951 Communication between VMD and {\gromacs} is achieved via TCP sockets and thus enables
1952 controlling an {\tt mdrun} running locally or on a remote cluster. The port for the
1953 connection can be specified with the {\tt -imdport} switch of {\tt mdrun}, 8888 is
1954 the default. If a port number of 0 or smaller is provided, {\gromacs} automatically
1955 assigns a free port to use with IMD.
1957 Every $N$ steps, the {\tt mdrun} client receives the applied forces from VMD and sends the new
1958 positions to the client. VMD permits increasing or decreasing the communication frequency
1959 interactively.
1960 By default, the simulation starts and runs even if no IMD client is connected. This
1961 behavior is changed by the {\tt -imdwait} switch of {\tt mdrun}. After startup and
1962 whenever the client has disconnected, the integration stops until reconnection of
1963 the client.
1964 When the {\tt -imdterm} switch is used, the simulation can be terminated by pressing
1965 the stop button in VMD. This is disabled by default.
1966 Finally, to allow interacting with the simulation (i.e. pulling from VMD) the {\tt -imdpull}
1967 switch has to be used.
1968 Therefore, a simulation can only be monitored but not influenced from the VMD client
1969 when none of {\tt -imdwait}, {\tt -imdterm} or {\tt -imdpull} are set. However, since
1970 the IMD protocol requires no authentication, it is not advisable to run simulations on
1971 a host directly reachable from an insecure environment. Secure shell forwarding of TCP
1972 can be used to connect to running simulations not directly reachable from the interacting host.
1973 Note that the IMD command line switches of {\tt mdrun} are hidden by default and show
1974 up in the help text only with {\tt gmx mdrun -h -hidden}.
1976 \subsection{Connecting from VMD}
1977 In VMD, first the structure corresponding to the IMD group has to be loaded ({\it File
1978 $\rightarrow$ New Molecule}). Then the IMD connection window has to be used
1979 ({\it Extensions $\rightarrow$ Simulation $\rightarrow$ IMD Connect (NAMD)}). In the IMD
1980 connection window, hostname and port have to be specified and followed by pressing
1981 {\it Connect}. {\it Detach Sim} allows disconnecting without terminating the simulation, while
1982 {\it Stop Sim} ends the simulation on the next neighbor searching step (if allowed by
1983 {\tt -imdterm}).
1985 The timestep transfer rate allows adjusting the communication frequency between simulation
1986 and IMD client. Setting the keep rate loads every $N^\mathrm{th}$ frame into VMD instead
1987 of discarding them when a new one is received. The displayed energies are in SI units
1988 in contrast to energies displayed from NAMD simulations.
1991 % LocalWords: PMF pmf kJ mol Jarzynski BT bilayer rup mdp AFM fepmf fecalc rb
1992 % LocalWords: posre grompp fs Verlet dihedrals hydrogens hydroxyl sulfhydryl
1993 % LocalWords: vsitehydro planarity chirality pdb gmx virtualize virtualized xz
1994 % LocalWords: vis massless tryptophan histidine phenyl parameterizing ij PPPM
1995 % LocalWords: parameterization Berendsen rlist coulombtype rcoulomb vdwtype LJ
1996 % LocalWords: rvdw energygrp mdrun pre GMXLIB MOPAC GAMESS CPMD ONIOM
1997 % LocalWords: Morokuma iJ AQ AJ initio atomtype QMatom MMatom QMMM grps et al
1998 % LocalWords: QMmethod QMbasis QMMMscheme RHF DFT LYP CASSCF MMVB CASelectrons
1999 % LocalWords: CASorbitals planewavecutoff STO QMcharge QMmult diabatic edr
2000 % LocalWords: hyperline delocalized Coulombic indices nm ccc th ps
2001 % LocalWords: GTX CPUs GHz md sd bd vv avek tcoupl andersen tc OPLSAA GROMOS
2002 % LocalWords: OBC obc CCMA tol pbc xyz barostat pcoupl acc gpu PLUGIN Cmake GX
2003 % LocalWords: MSVC gcc installpath cmake DGMX DCMAKE functionalities GPGPU GTS
2004 % LocalWords: OpenCL DeviceID gromacs gpus html GeForce Quadro FX Plex CX GF
2005 % LocalWords: VMD DCD GROningen MAchine BIOSON Groningen der Spoel Drunen Comp
2006 % LocalWords: Phys Comm moltype intramol vdw Waals fep multivector pf
2007 % LocalWords: pymbar dhdl xvg LINCS Entropic entropic solutes ref com iso pm
2008 % LocalWords: rm prefactors equipotential potiso potisopf potpm trr
2009 % LocalWords: potrm potrmpf midplanes midplane gaussians potflex vars massw av
2010 % LocalWords: shure observables rccccccc vec eps dist min eqn transl nstsout
2011 % LocalWords: nstrout rotangles rotslabs rottorque RMSD rmsdfit excl NH amine
2012 % LocalWords: positionrestraint es SH phenylalanine solvated sh nanometer QM
2013 % LocalWords: Lennard Buckingham UK Svensson ab vsitetop co UHF MP interatomic
2014 % LocalWords: cg grp coords SPC userpot
2015 % LocalWords: itp sitesn atomtypes ff csg ndx Tesla CHARMM AA gauss
2016 % LocalWords: CMAP nocmap fators Monte performant lib LD DIR llllcc
2017 % LocalWords: CMake homepage DEV overclocking GT dlopen vmd VMDDIR
2018 % LocalWords: versa PME atomperatom graining forcefields hy spc OPLS
2019 % LocalWords: topol multi