6fd043c49c0ad910049c7148343c6ce34e3c72d8
[alexxy/gromacs.git] / docs / reference-manual / functions / free-energy-interactions.rst
1 .. _feia:
2
3 Free energy interactions
4 ------------------------
5
6 This section describes the :math:`\lambda`-dependence of the potentials
7 used for free energy calculations (see sec. :ref:`fecalc`). All common
8 types of potentials and constraints can be interpolated smoothly from
9 state A (:math:`\lambda=0`) to state B (:math:`\lambda=1`) and vice
10 versa. All bonded interactions are interpolated by linear interpolation
11 of the interaction parameters. Non-bonded interactions can be
12 interpolated linearly or via soft-core interactions.
13
14 Starting in |Gromacs| 4.6, :math:`\lambda` is a vector, allowing different
15 components of the free energy transformation to be carried out at
16 different rates. Coulomb, Lennard-Jones, bonded, and restraint terms can
17 all be controlled independently, as described in the
18 :ref:`mdp` options.
19
20 Harmonic potentials
21 ~~~~~~~~~~~~~~~~~~~
22
23 The example given here is for the bond potential, which is harmonic in
24 |Gromacs|. However, these equations apply to the angle potential and the
25 improper dihedral potential as well.
26
27 .. math::
28
29    \begin{aligned}
30    V_b     &=&{\frac{1}{2}}\left[{(1-{\lambda})}k_b^A + 
31                    {\lambda}k_b^B\right] \left[b - {(1-{\lambda})}b_0^A - {\lambda}b_0^B\right]^2  \\
32    {\frac{\partial V_b}{\partial {\lambda}}}&=&{\frac{1}{2}}(k_b^B-k_b^A)
33                    \left[b - {(1-{\lambda})}b_0^A + {\lambda}b_0^B\right]^2 + 
34                 \nonumber\\
35            & & \phantom{{\frac{1}{2}}}(b_0^A-b_0^B) \left[b - {(1-{\lambda})}b_0^A -{\lambda}b_0^B\right]
36                 \left[{(1-{\lambda})}k_b^A + {\lambda}k_b^B \right]\end{aligned}
37
38 GROMOS-96 bonds and angles
39 ~~~~~~~~~~~~~~~~~~~~~~~~~~
40
41 Fourth-power bond stretching and cosine-based angle potentials are
42 interpolated by linear interpolation of the force constant and the
43 equilibrium position. Formulas are not given here.
44
45 Proper dihedrals
46 ~~~~~~~~~~~~~~~~
47
48 For the proper dihedrals, the equations are somewhat more complicated:
49
50 .. math::
51
52    \begin{aligned}
53    V_d     &=&\left[{(1-{\lambda})}k_d^A + {\lambda}k_d^B \right]
54            \left( 1+ \cos\left[n_{\phi} \phi - 
55                     {(1-{\lambda})}\phi_s^A - {\lambda}\phi_s^B
56                     \right]\right)\\
57    {\frac{\partial V_d}{\partial {\lambda}}}&=&(k_d^B-k_d^A) 
58             \left( 1+ \cos
59                  \left[
60                     n_{\phi} \phi- {(1-{\lambda})}\phi_s^A - {\lambda}\phi_s^B
61                  \right]
62          \right) +
63          \nonumber\\
64            &&(\phi_s^B - \phi_s^A) \left[{(1-{\lambda})}k_d^A - {\lambda}k_d^B\right] 
65            \sin\left[  n_{\phi}\phi - {(1-{\lambda})}\phi_s^A - {\lambda}\phi_s^B \right]\end{aligned}
66
67 **Note:** that the multiplicity :math:`n_{\phi}` can not be
68 parameterized because the function should remain periodic on the
69 interval :math:`[0,2\pi]`.
70
71 Tabulated bonded interactions
72 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
73
74 For tabulated bonded interactions only the force constant can
75 interpolated:
76
77 .. math::
78
79    \begin{aligned}
80          V  &=& ({(1-{\lambda})}k^A + {\lambda}k^B) \, f \\
81    {\frac{\partial V}{\partial {\lambda}}} &=& (k^B - k^A) \, f\end{aligned}
82
83 Coulomb interaction
84 ~~~~~~~~~~~~~~~~~~~
85
86 The Coulomb interaction between two particles of which the charge varies
87 with :math:`{\lambda}` is:
88
89 .. math::
90
91    \begin{aligned}
92    V_c &=& \frac{f}{{\varepsilon_{rf}}{r_{ij}}}\left[{(1-{\lambda})}q_i^A q_j^A + {\lambda}\, q_i^B q_j^B\right] \\
93    {\frac{\partial V_c}{\partial {\lambda}}}&=& \frac{f}{{\varepsilon_{rf}}{r_{ij}}}\left[- q_i^A q_j^A + q_i^B q_j^B\right]\end{aligned}
94
95 where :math:`f = \frac{1}{4\pi \varepsilon_0} = {138.935\,458}` (see
96 chapter :ref:`defunits`).
97
98 Coulomb interaction with reaction field
99 ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
100
101 The Coulomb interaction including a reaction field, between two
102 particles of which the charge varies with :math:`{\lambda}` is:
103
104 .. math:: \begin{aligned}
105           V_c     &=& f\left[\frac{1}{{r_{ij}}} + k_{rf}~ {r_{ij}}^2 -c_{rf}\right]
106           \left[{(1-{\lambda})}q_i^A q_j^A + {\lambda}\, q_i^B q_j^B\right] \\
107           {\frac{\partial V_c}{\partial {\lambda}}}&=& f\left[\frac{1}{{r_{ij}}} + k_{rf}~ {r_{ij}}^2 -c_{rf}\right]
108           \left[- q_i^A q_j^A + q_i^B q_j^B\right]
109           \end{aligned}
110           :label: eqdVcoulombdlambda
111
112 **Note** that the constants :math:`k_{rf}` and :math:`c_{rf}` are
113 defined using the dielectric constant :math:`{\varepsilon_{rf}}` of the
114 medium (see sec. :ref:`coulrf`).
115
116 Lennard-Jones interaction
117 ~~~~~~~~~~~~~~~~~~~~~~~~~
118
119 For the Lennard-Jones interaction between two particles of which the
120 *atom type* varies with :math:`{\lambda}` we can write:
121
122 .. math:: \begin{aligned}
123           V_{LJ}  &=&     \frac{{(1-{\lambda})}C_{12}^A + {\lambda}\, C_{12}^B}{{r_{ij}}^{12}} -
124           \frac{{(1-{\lambda})}C_6^A + {\lambda}\, C_6^B}{{r_{ij}}^6}   \\
125           {\frac{\partial V_{LJ}}{\partial {\lambda}}}&=&\frac{C_{12}^B - C_{12}^A}{{r_{ij}}^{12}} -
126           \frac{C_6^B - C_6^A}{{r_{ij}}^6}
127           \end{aligned}
128           :label: eqdVljdlambda
129
130 It should be noted that it is also possible to express a pathway from
131 state A to state B using :math:`\sigma` and :math:`\epsilon` (see
132 :eq:`eqn. %s <eqnsigeps>`). It may seem to make sense physically to vary the
133 force field parameters :math:`\sigma` and :math:`\epsilon` rather than
134 the derived parameters :math:`C_{12}` and :math:`C_{6}`. However, the
135 difference between the pathways in parameter space is not large, and the
136 free energy itself does not depend on the pathway, so we use the simple
137 formulation presented above.
138
139 Kinetic Energy
140 ~~~~~~~~~~~~~~
141
142 When the mass of a particle changes, there is also a contribution of the
143 kinetic energy to the free energy (note that we can not write the
144 momentum :math:`\mathbf{p}` as
145 m :math:`\mathbf{v}`, since that would result in the
146 sign of :math:`{\frac{\partial E_k}{\partial {\lambda}}}` being
147 incorrect \ :ref:`99 <refGunsteren98a>`):
148
149 .. math::
150
151    \begin{aligned}
152    E_k      &=&     {\frac{1}{2}}\frac{\mathbf{p}^2}{{(1-{\lambda})}m^A + {\lambda}m^B}        \\
153    {\frac{\partial E_k}{\partial {\lambda}}}&=&    -{\frac{1}{2}}\frac{\mathbf{p}^2(m^B-m^A)}{({(1-{\lambda})}m^A + {\lambda}m^B)^2}\end{aligned}
154
155 after taking the derivative, we *can* insert
156 :math:`\mathbf{p}` =
157 m :math:`\mathbf{v}`, such that:
158
159 .. math:: {\frac{\partial E_k}{\partial {\lambda}}}~=~    -{\frac{1}{2}}\mathbf{v}^2(m^B-m^A)
160
161 Constraints
162 ~~~~~~~~~~~
163
164 The constraints are formally part of the Hamiltonian, and therefore they
165 give a contribution to the free energy. In |Gromacs| this can be
166 calculated using the LINCS or the SHAKE algorithm. If we have
167 :math:`k = 1 \ldots K` constraint equations :math:`g_k` for LINCS, then
168
169 .. math:: g_k     =       | \mathbf{r}_{k} | - d_{k}
170
171 where :math:`\mathbf{r}_k` is the displacement vector
172 between two particles and :math:`d_k` is the constraint distance between
173 the two particles. We can express the fact that the constraint distance
174 has a :math:`{\lambda}` dependency by
175
176 .. math:: d_k     =       {(1-{\lambda})}d_{k}^A + {\lambda}d_k^B
177
178 Thus the :math:`{\lambda}`-dependent constraint equation is
179
180 .. math:: g_k     =       | \mathbf{r}_{k} | - \left({(1-{\lambda})}d_{k}^A + {\lambda}d_k^B\right).
181
182 The (zero) contribution :math:`G` to the Hamiltonian from the
183 constraints (using Lagrange multipliers :math:`\lambda_k`, which are
184 logically distinct from the free-energy :math:`{\lambda}`) is
185
186 .. math::
187
188    \begin{aligned}
189    G           &=&     \sum^K_k \lambda_k g_k    \\
190    {\frac{\partial G}{\partial {\lambda}}}    &=&     \frac{\partial G}{\partial d_k} {\frac{\partial d_k}{\partial {\lambda}}} \\
191                &=&     - \sum^K_k \lambda_k \left(d_k^B-d_k^A\right)\end{aligned}
192
193 For SHAKE, the constraint equations are
194
195 .. math:: g_k     =       \mathbf{r}_{k}^2 - d_{k}^2
196
197 with :math:`d_k` as before, so
198
199 .. math::
200
201    \begin{aligned}
202    {\frac{\partial G}{\partial {\lambda}}}    &=&     -2 \sum^K_k \lambda_k \left(d_k^B-d_k^A\right)\end{aligned}
203
204 Soft-core interactions
205 ~~~~~~~~~~~~~~~~~~~~~~
206
207 .. _fig-softcore:
208
209 .. figure:: plots/softcore.*
210    :height: 6.00000cm
211
212    Soft-core interactions at :math:`{\lambda}=0.5`, with :math:`p=2` and
213    :math:`C_6^A=C_{12}^A=C_6^B=C_{12}^B=1`.
214
215 In a free-energy calculation where particles grow out of nothing, or
216 particles disappear, using the the simple linear interpolation of the
217 Lennard-Jones and Coulomb potentials as described in
218 :eq:`Equations %s <eqdVljdlambda>` and :eq:`%s <eqdVcoulombdlambda>` may lead to poor
219 convergence. When the particles have nearly disappeared, or are close to
220 appearing (at :math:`{\lambda}` close to 0 or 1), the interaction energy
221 will be weak enough for particles to get very close to each other,
222 leading to large fluctuations in the measured values of
223 :math:`\partial V/\partial {\lambda}` (which, because of the simple
224 linear interpolation, depends on the potentials at both the endpoints of
225 :math:`{\lambda}`).
226
227 To circumvent these problems, the singularities in the potentials need
228 to be removed. This can be done by modifying the regular Lennard-Jones
229 and Coulomb potentials with “soft-core” potentials that limit the
230 energies and forces involved at :math:`{\lambda}` values between 0 and
231 1, but not *at* :math:`{\lambda}=0` or 1.
232
233 In |Gromacs| the soft-core potentials :math:`V_{sc}` are shifted versions
234 of the regular potentials, so that the singularity in the potential and
235 its derivatives at :math:`r=0` is never reached:
236
237 .. math::
238
239    \begin{aligned}
240    V_{sc}(r) &=& {(1-{\lambda})}V^A(r_A) + {\lambda}V^B(r_B)
241        \\
242    r_A &=& \left(\alpha \sigma_A^6 {\lambda}^p + r^6 \right)^\frac{1}{6}
243        \\
244    r_B &=& \left(\alpha \sigma_B^6 {(1-{\lambda})}^p + r^6 \right)^\frac{1}{6}\end{aligned}
245
246 where :math:`V^A` and :math:`V^B` are the normal “hard core” Van der
247 Waals or electrostatic potentials in state A (:math:`{\lambda}=0`) and
248 state B (:math:`{\lambda}=1`) respectively, :math:`\alpha` is the
249 soft-core parameter (set with ``sc_alpha`` in the
250 :ref:`mdp` file), :math:`p` is the soft-core :math:`{\lambda}`
251 power (set with ``sc_power``), :math:`\sigma` is the radius
252 of the interaction, which is :math:`(C_{12}/C_6)^{1/6}` or an input
253 parameter (``sc_sigma``) when :math:`C_6` or :math:`C_{12}`
254 is zero.
255
256 For intermediate :math:`{\lambda}`, :math:`r_A` and :math:`r_B` alter
257 the interactions very little for :math:`r > \alpha^{1/6} \sigma` and
258 quickly switch the soft-core interaction to an almost constant value for
259 smaller :math:`r` (:numref:`Fig. %s <fig-softcore>`). The force is:
260
261 .. math::
262
263    F_{sc}(r) = -\frac{\partial V_{sc}(r)}{\partial r} =
264     {(1-{\lambda})}F^A(r_A) \left(\frac{r}{r_A}\right)^5 +
265    {\lambda}F^B(r_B) \left(\frac{r}{r_B}\right)^5
266
267 where :math:`F^A` and :math:`F^B` are the “hard core” forces. The
268 contribution to the derivative of the free energy is:
269
270 .. math::
271
272    \begin{aligned}
273    {\frac{\partial V_{sc}(r)}{\partial {\lambda}}} & = &
274     V^B(r_B) -V^A(r_A)  + 
275         {(1-{\lambda})}\frac{\partial V^A(r_A)}{\partial r_A}
276                    \frac{\partial r_A}{\partial {\lambda}} + 
277         {\lambda}\frac{\partial V^B(r_B)}{\partial r_B}
278                    \frac{\partial r_B}{\partial {\lambda}}
279    \nonumber\\
280    &=&
281     V^B(r_B) -V^A(r_A)  + \nonumber \\
282     & &
283     \frac{p \alpha}{6}
284           \left[ {\lambda}F^B(r_B) r^{-5}_B \sigma_B^6 {(1-{\lambda})}^{p-1} -
285                {(1-{\lambda})}F^A(r_A) r^{-5}_A \sigma_A^6 {\lambda}^{p-1} \right]\end{aligned}
286
287 The original GROMOS Lennard-Jones soft-core
288 function\ :ref:`100 <refBeutler94>` uses :math:`p=2`, but :math:`p=1` gives a smoother
289 :math:`\partial H/\partial{\lambda}` curve. Another issue that should be
290 considered is the soft-core effect of hydrogens without Lennard-Jones
291 interaction. Their soft-core :math:`\sigma` is set with
292 ``sc_sigma`` in the :ref:`mdp` file. These
293 hydrogens produce peaks in :math:`\partial H/\partial{\lambda}` at
294 :math:`{\lambda}` is 0 and/or 1 for :math:`p=1` and close to 0 and/or 1
295 with :math:`p=2`. Lowering ``sc_sigma``
296 will decrease this effect, but it will also increase the interactions
297 with hydrogens relative to the other interactions in the soft-core
298 state.
299
300 When soft-core potentials are selected (by setting
301 ``sc_alpha >0``), and the Coulomb and Lennard-Jones
302 potentials are turned on or off sequentially, then the Coulombic
303 interaction is turned off linearly, rather than using soft-core
304 interactions, which should be less statistically noisy in most cases.
305 This behavior can be overwritten by using the :ref:`mdp` option
306 ``sc-coul`` to ``yes``. Note that the
307 ``sc-coul`` is only taken into account when lambda states
308 are used, not with ``couple-lambda0``  /
309 ``couple-lambda1``, and you can still turn off soft-core
310 interactions by setting ``sc-alpha=0``. Additionally, the
311 soft-core interaction potential is only applied when either the A or B
312 state has zero interaction potential. If both A and B states have
313 nonzero interaction potential, default linear scaling described above is
314 used. When both Coulombic and Lennard-Jones interactions are turned off
315 simultaneously, a soft-core potential is used, and a hydrogen is being
316 introduced or deleted, the sigma is set to ``sc-sigma-min``,
317 which itself defaults to ``sc-sigma-default``.
318
319 Recently, a new formulation of the soft-core approach has been derived
320 that in most cases gives lower and more even statistical variance than
321 the standard soft-core path described above \ :ref:`101 <refPham2011>`,
322 :ref:`102 <refPham2012>`. Specifically, we have:
323
324 .. math::
325
326    \begin{aligned}
327    V_{sc}(r) &=& {(1-{\lambda})}V^A(r_A) + {\lambda}V^B(r_B)
328        \\
329    r_A &=& \left(\alpha \sigma_A^{48} {\lambda}^p + r^{48} \right)^\frac{1}{48}
330        \\
331    r_B &=& \left(\alpha \sigma_B^{48} {(1-{\lambda})}^p + r^{48} \right)^\frac{1}{48}\end{aligned}
332
333 This “1-1-48” path is also implemented in |Gromacs|. Note that for this
334 path the soft core :math:`\alpha` should satisfy
335 :math:`0.001 < \alpha < 0.003`, rather than :math:`\alpha \approx
336 0.5`.