Remove COM from SETTLE
[gromacs.git] / docs / reference-manual / algorithms / constraint-algorithms.rst
blobcf5f70f94ecfebf8e9de05a487224d685bf8fdac
1 Constraint algorithms
2 ---------------------
4 Constraints can be imposed in |Gromacs| using LINCS (default) or the
5 traditional SHAKE method.
8 .. _shake:
10 SHAKE
11 ~~~~~
13 The SHAKE \ :ref:`46 <refRyckaert77>` algorithm changes a
14 set of unconstrained coordinates :math:`\mathbf{r}^{'}` to
15 a set of coordinates :math:`\mathbf{r}''` that fulfill a
16 list of distance constraints, using a set :math:`\mathbf{r}` reference, as
18 .. math:: {\rm SHAKE}(\mathbf{r}^{'} \rightarrow \mathbf{r}'';\, \mathbf{r})
19           :label: eqnshakebase
21 This action is consistent with solving a set of Lagrange multipliers in
22 the constrained equations of motion. SHAKE needs a *relative tolerance*;
23 it will continue until all constraints are satisfied within that
24 relative tolerance. An error message is given if SHAKE cannot reset the
25 coordinates because the deviation is too large, or if a given number of
26 iterations is surpassed.
28 Assume the equations of motion must fulfill :math:`K` holonomic
29 constraints, expressed as
31 .. math:: \sigma_k(\mathbf{r}_1 \ldots \mathbf{r}_N) = 0; \;\; k=1 \ldots K.
32           :label: eqnshakemotconstr
34 For example, :math:`(\mathbf{r}_1 - \mathbf{r}_2)^2 - b^2 = 0`.
35 Then the forces are defined as
37 .. math:: - \frac{\partial}{\partial \mathbf{r}_i} \left( V + \sum_{k=1}^K \lambda_k
38           \sigma_k \right),
39           :label: eqnshakeforce
41 where :math:`\lambda_k` are Lagrange multipliers which must be solved
42 to fulfill the constraint equations. The second part of this sum
43 determines the *constraint forces* :math:`\mathbf{G}_i`, defined by
45 .. math:: \mathbf{G}_i = -\sum_{k=1}^K \lambda_k \frac{\partial \sigma_k}{\partial
46           \mathbf{r}_i}
47           :label: eqnshakeconstrforces
49 The displacement due to the constraint forces in the leap-frog or
50 Verlet algorithm is equal to :math:`(\mathbf{G}_i/m_i)({{\Delta t}})^2`. Solving the
51 Lagrange multipliers (and hence the displacements) requires the solution
52 of a set of coupled equations of the second degree. These are solved
53 iteratively by SHAKE. :ref:`settle` 
55 .. _settle:
57 SETTLE
58 ~~~~~~
60 For the special case of rigid
61 water molecules, that often make up more than 80% of the simulation
62 system we have implemented the SETTLE algorithm \ :ref:`47 <refMiyamoto92>`
63 (sec. :ref:`constraintalg`). The implementation of SETTLE in |Gromacs|
64 is a slight modification of the original algorithm, in that it completely
65 avoids the calculation of the center of mass of the water molecule.
66 Apart from saving a few operations, the main gain of this is a reduction
67 in rouding errors. For large coordinates, the floating pointing precision of constrained
68 distances is reduced, which leads to an energy drift which usually depends
69 quadratically on the coordinate. For SETTLE this dependence is now linear, which enables
70 accurate integration of systems in single precision up to 1000 nm in size.
71 But note that the drift due to SHAKE and LINCS still has a quadratic
72 dependence, which limits the size of systems with normal constraints
73 in single precision to 100 to 200 nm.
75 For velocity Verlet, an additional round of constraining must be done,
76 to constrain the velocities of the second velocity half step, removing
77 any component of the velocity parallel to the bond vector. This step is
78 called RATTLE, and is covered in more detail in the original Andersen
79 paper \ :ref:`48 <refAndersen1983a>`.
81 LINCS
82 ~~~~~
84 .. _lincs:
86 The LINCS algorithm
87 ^^^^^^^^^^^^^^^^^^^
89 LINCS is an algorithm that resets bonds to their correct lengths after
90 an unconstrained update \ :ref:`49 <refHess97>`. The method is non-iterative,
91 as it always uses two steps. Although LINCS is based on matrices, no
92 matrix-matrix multiplications are needed. The method is more stable and
93 faster than SHAKE, but it can only be used with bond constraints and
94 isolated angle constraints, such as the proton angle in OH. Because of
95 its stability, LINCS is especially useful for Brownian dynamics. LINCS
96 has two parameters, which are explained in the subsection parameters.
97 The parallel version of LINCS, P-LINCS, is described in subsection
98 :ref:`plincs`.
100 The LINCS formulas
101 ^^^^^^^^^^^^^^^^^^
103 We consider a system of :math:`N` particles, with positions given by a
104 :math:`3N` vector :math:`\mathbf{r}(t)`. For molecular
105 dynamics the equations of motion are given by Newton’s Law
107 .. math:: {{\mbox{d}}^2 \mathbf{r} \over {\mbox{d}}t^2} = {{\mathbf{M}}^{-1}}\mathbf{F},
108           :label: eqnc1
110 where :math:`\mathbf{F}` is the :math:`3N` force vector
111 and :math:`{\mathbf{M}}` is a :math:`3N \times 3N`
112 diagonal matrix, containing the masses of the particles. The system is
113 constrained by :math:`K` time-independent constraint equations
115 .. math:: g_i(\mathbf{r}) = | \mathbf{r}_{i_1}-\mathbf{r}_{i_2} | - d_i = 0 ~~~~~~i=1,\ldots,K.
116           :label: eqnc2
118 In a numerical integration scheme, LINCS is applied after an
119 unconstrained update, just like SHAKE. The algorithm works in two steps
120 (see figure :numref:`Fig. %s <fig-lincs>`). In the first step, the
121 projections of the new bonds on the old bonds are set to zero. In the
122 second step, a correction is applied for the lengthening of the bonds
123 due to rotation. The numerics for the first step and the second step are
124 very similar. A complete derivation of the algorithm can be found in
125 :ref:`49 <refHess97>`. Only a short description of the first step is given
126 here.
128 .. _fig-lincs:
130 .. figure:: plots/lincs.*
131    :height: 5.00000cm
133    The three position updates needed for one time step. The dashed
134    line is the old bond of length :math:`d`, the solid lines are the new
135    bonds. :math:`l=d \cos \theta` and :math:`p=(2 d^2 - l^2)^{1 \over 2}`.
137 A new notation is introduced for the gradient matrix of the constraint
138 equations which appears on the right hand side of this equation:
140 .. math::  B_{hi} = {{\partial}g_h \over {\partial}r_i}
141            :label: eqnc3
143 Notice that :math:`{\mathbf{B}}` is a :math:`K \times 3N`
144 matrix, it contains the directions of the constraints. The following
145 equation shows how the new constrained coordinates
146 :math:`\mathbf{r}_{n+1}` are related to the unconstrained
147 coordinates :math:`\mathbf{r}_{n+1}^{unc}` by
149 .. math::  \begin{array}{c}
150            \mathbf{r}_{n+1}=(\mathbf{I}-\mathbf{T}_n \mathbf{B}_n) \mathbf{r}_{n+1}^{unc} + {\mathbf{T}}_n \mathbf{d}=  
151            \\[2mm]
152            \mathbf{r}_{n+1}^{unc} - 
153            {{\mathbf{M}}^{-1}}\mathbf{B}_n ({\mathbf{B}}_n {{\mathbf{M}}^{-1}}{\mathbf{B}}_n^T)^{-1} ({\mathbf{B}}_n \mathbf{r}_{n+1}^{unc} - \mathbf{d}) 
154            \end{array}
155            :label: eqnm0
157 where
159 .. math:: {\mathbf{T}}= {{\mathbf{M}}^{-1}}{\mathbf{B}}^T ({\mathbf{B}}{{\mathbf{M}}^{-1}}{\mathbf{B}}^T)^{-1}
160           :label: eqnnm01
162 The derivation of this equation from :eq:`eqns. %s <eqnc1>` and
163 :eq:`%s <eqnc2>` can be found in :ref:`49 <refHess97>`.
165 This first step does not set the real bond lengths to the prescribed
166 lengths, but the projection of the new bonds onto the old directions of
167 the bonds. To correct for the rotation of bond :math:`i`, the projection
168 of the bond, :math:`p_i`, on the old direction is set to
170 .. math::  p_i=\sqrt{2 d_i^2 - l_i^2},
171            :label: eqnm1a
173 where :math:`l_i` is the bond length after the first projection. The
174 corrected positions are
176 .. math::  \mathbf{r}_{n+1}^*=(\mathbf{I}-\mathbf{T}_n \mathbf{B}_n)\mathbf{r}_{n+1} + {\mathbf{T}}_n \mathbf{p}.
177            :label: eqnm1b
179 This correction for rotational effects is actually an iterative
180 process, but during MD only one iteration is applied. The relative
181 constraint deviation after this procedure will be less than 0.0001 for
182 every constraint. In energy minimization, this might not be accurate
183 enough, so the number of iterations is equal to the order of the
184 expansion (see below).
186 Half of the CPU time goes to inverting the constraint coupling matrix
187 :math:`{\mathbf{B}}_n {{\mathbf{M}}^{-1}}{\mathbf{B}}_n^T`,
188 which has to be done every time step. This :math:`K \times K` matrix has
189 :math:`1/m_{i_1} + 1/m_{i_2}` on the diagonal. The off-diagonal elements
190 are only non-zero when two bonds are connected, then the element is
191 :math:`\cos \phi /m_c`, where :math:`m_c` is the mass of the atom
192 connecting the two bonds and :math:`\phi` is the angle between the
193 bonds.
195 The matrix :math:`\mathbf{T}` is inverted through a power
196 expansion. A :math:`K \times K` matrix :math:`\mathbf{S}`
197 is introduced which is the inverse square root of the diagonal of
198 :math:`\mathbf{B}_n {{\mathbf{M}}^{-1}}{\mathbf{B}}_n^T`.
199 This matrix is used to convert the diagonal elements of the coupling
200 matrix to one:
202 .. math:: \begin{array}{c}
203           ({\mathbf{B}}_n {{\mathbf{M}}^{-1}}{\mathbf{B}}_n^T)^{-1}
204           = {\mathbf{S}}{\mathbf{S}}^{-1} ({\mathbf{B}}_n {{\mathbf{M}}^{-1}}{\mathbf{B}}_n^T)^{-1} {\mathbf{S}}^{-1} {\mathbf{S}}\\[2mm]
205           = {\mathbf{S}}({\mathbf{S}}{\mathbf{B}}_n {{\mathbf{M}}^{-1}}{\mathbf{B}}_n^T {\mathbf{S}})^{-1} {\mathbf{S}}=
206           {\mathbf{S}}(\mathbf{I} - \mathbf{A}_n)^{-1} {\mathbf{S}}\end{array}
207           :label: eqnm2
209 The matrix :math:`\mathbf{A}_n` is symmetric and sparse
210 and has zeros on the diagonal. Thus a simple trick can be used to
211 calculate the inverse:
213 .. math:: (\mathbf{I}-\mathbf{A}_n)^{-1}= 
214           \mathbf{I} + \mathbf{A}_n + \mathbf{A}_n^2 + \mathbf{A}_n^3 + \ldots
215           :label: eqnm3
217 This inversion method is only valid if the absolute values of all the
218 eigenvalues of :math:`\mathbf{A}_n` are smaller than one.
219 In molecules with only bond constraints, the connectivity is so low that
220 this will always be true, even if ring structures are present. Problems
221 can arise in angle-constrained molecules. By constraining angles with
222 additional distance constraints, multiple small ring structures are
223 introduced. This gives a high connectivity, leading to large
224 eigenvalues. Therefore LINCS should NOT be used with coupled
225 angle-constraints.
227 For molecules with all bonds constrained the eigenvalues of :math:`A`
228 are around 0.4. This means that with each additional order in the
229 expansion :eq:`eqn. %s <eqnm3>` the deviations decrease by a factor 0.4. But for
230 relatively isolated triangles of constraints the largest eigenvalue is
231 around 0.7. Such triangles can occur when removing hydrogen angle
232 vibrations with an additional angle constraint in alcohol groups or when
233 constraining water molecules with LINCS, for instance with flexible
234 constraints. The constraints in such triangles converge twice as slow as
235 the other constraints. Therefore, starting with |Gromacs| 4, additional
236 terms are added to the expansion for such triangles
238 .. math:: (\mathbf{I}-\mathbf{A}_n)^{-1} \approx
239           \mathbf{I} + \mathbf{A}_n + \ldots + \mathbf{A}_n^{N_i} +
240           \left(\mathbf{A}^*_n + \ldots + {\mathbf{A}_n^*}^{N_i} \right) \mathbf{A}_n^{N_i}
241           :label: eqnm3ang
243 where :math:`N_i` is the normal order of the expansion and
244 :math:`\mathbf{A}^*` only contains the elements of
245 :math:`\mathbf{A}` that couple constraints within rigid
246 triangles, all other elements are zero. In this manner, the accuracy of
247 angle constraints comes close to that of the other constraints, while
248 the series of matrix vector multiplications required for determining the
249 expansion only needs to be extended for a few constraint couplings. This
250 procedure is described in the P-LINCS paper\ :ref:`50 <refHess2008a>`.
252 The LINCS Parameters
253 ^^^^^^^^^^^^^^^^^^^^
255 The accuracy of LINCS depends on the number of matrices used in the
256 expansion :eq:`eqn. %s <eqnm3>`. For MD calculations a fourth order expansion is
257 enough. For Brownian dynamics with large time steps an eighth order
258 expansion may be necessary. The order is a parameter in the :ref:`mdp` file.
259 The implementation of LINCS is done in such a way that the algorithm
260 will never crash. Even when it is impossible to to reset the constraints
261 LINCS will generate a conformation which fulfills the constraints as
262 well as possible. However, LINCS will generate a warning when in one
263 step a bond rotates over more than a predefined angle. This angle is set
264 by the user in the :ref:`mdp` file.