Files
lammps/doc/src/Howto_bioFF.rst

269 lines
11 KiB
ReStructuredText

CHARMM, AMBER, COMPASS, and DREIDING force fields
=================================================
A compact summary of the concepts, definitions, and properties of
force fields with explicit bonded interactions (like the ones discussed
in this HowTo) is given in :ref:`(Gissinger) <Typelabel2>`.
A force field has 2 parts: the formulas that define it and the
coefficients used for a particular system. Here we only discuss
formulas implemented in LAMMPS that correspond to formulas commonly used
in the CHARMM, AMBER, COMPASS, and DREIDING force fields. Setting
coefficients is done either from special sections in an input data file
via the :doc:`read_data <read_data>` command or in the input script with
commands like :doc:`pair_coeff <pair_coeff>` or :doc:`bond_coeff
<bond_coeff>` and so on. See the :doc:`Tools <Tools>` doc page for
additional tools that can use CHARMM, AMBER, or Materials Studio
generated files to assign force field coefficients and convert their
output into LAMMPS input. LAMMPS input scripts can also be generated by
`charmm-gui.org <https://charmm-gui.org/>`_.
CHARMM and AMBER
----------------
The `CHARMM force field
<https://mackerell.umaryland.edu/charmm_ff.shtml>`_ :ref:`(MacKerell)
<howto-MacKerell>` and `AMBER force field
<https://ambermd.org/AmberModels.php>`_ :ref:`(Cornell) <howto-Cornell>`
have potential energy function of the form
.. math::
V & = \sum_{bonds} E_b + \sum_{angles} \!E_a + \!\overbrace{\sum_{dihedral} \!\!E_d}^{\substack{
\text{charmm} \\
\text{charmmfsw}
}} +\!\!\! \sum_{impropers} \!\!\!E_i \\[.6em]
& \quad + \!\!\!\!\!\!\!\!\!\!\underbrace{~\sum_{pairs} \left(E_{LJ}+E_{coul}\right)}_{\substack{
\text{lj/charmm/coul/charmm} \\
\text{lj/charmm/coul/charmm/implicit} \\
\text{lj/charmm/coul/long} \\
\text{lj/charmm/coul/msm} \\
\text{lj/charmmfsw/coul/charmmfsh} \\
\text{lj/charmmfsw/coul/long}
}} \!\!\!\!\!\!\!\!+ \!\!\sum_{special}\! E_s + \!\!\!\!\sum_{residues} \!\!\!{\scriptstyle\mathrm{CMAP}(\phi,\psi)}
The terms are computed by bond styles (relationship between two atoms),
angle styles (between 3 atoms) , dihedral/improper styles (between 4
atoms), pair styles (non-covalently bonded pair interactions) and
special bonds. The CMAP term (see :doc:`fix cmap <fix_cmap>` command for
details) corrects for pairs of dihedral angles ("Correction MAP") to
significantly improve the structural and dynamic properties of proteins
in crystalline and solution environments :ref:`(Brooks)
<howto-Brooks>`. The AMBER force field does not include the CMAP term.
The interaction styles listed below compute force field formulas that
are consistent with common options in CHARMM or AMBER. See each
command's documentation for the formula it computes.
* :doc:`bond_style <bond_harmonic>` harmonic
* :doc:`angle_style <angle_charmm>` charmm
* :doc:`dihedral_style <dihedral_charmm>` charmmfsh
* :doc:`dihedral_style <dihedral_charmm>` charmm
* :doc:`pair_style <pair_charmm>` lj/charmmfsw/coul/charmmfsh
* :doc:`pair_style <pair_charmm>` lj/charmmfsw/coul/long
* :doc:`pair_style <pair_charmm>` lj/charmm/coul/charmm
* :doc:`pair_style <pair_charmm>` lj/charmm/coul/charmm/implicit
* :doc:`pair_style <pair_charmm>` lj/charmm/coul/long
* :doc:`special_bonds <special_bonds>` charmm
* :doc:`special_bonds <special_bonds>` amber
The pair styles compute Lennard Jones (LJ) and Coulombic interactions
with additional switching or shifting functions that ramp the energy
and/or force smoothly to zero between an inner :math:`(a)` and outer
:math:`(b)` cutoff. The older styles with *charmm* (not *charmmfsw* or
*charmmfsh*\ ) in their name compute the LJ and Coulombic interactions
with an energy switching function (esw) :math:`S(r)` which ramps the energy
smoothly to zero between the inner and outer cutoff. This can cause
irregularities in pairwise forces (due to the discontinuous second
derivative of energy at the boundaries of the switching region), which
in some cases can result in complications in energy minimization and
detectable artifacts in MD simulations.
.. grid:: 1 1 2 2
.. grid-item::
.. math::
LJ(r) &= 4 \epsilon \left[ \left(\frac{\sigma}{r}\right)^{12} -
\left(\frac{\sigma}{r}\right)^6 \right]\\[.6em]
C(r) &= \frac{C q_i q_j}{ \epsilon r}\\[.6em]
S(r) &= \frac{ \left(b^2 - r^2\right)^2 \left(b^2 + 2r^2 - 3{a^2}\right)}
{ \left(b^2 - a^2\right)^3 }\\[.6em]
E_{LJ}(r) &= \begin{cases}
LJ(r), & r \leq a \\
LJ(r) S(r), & a < r \leq b \\
0, &r > b
\end{cases} \\[.6em]
E_{coul}(r) &= \begin{cases}
C(r), & r \leq a \\
C(r) S(r), & a < r \leq b \\
0, & r > b
\end{cases}
.. grid-item::
.. image:: img/howto_charmm_ELJ.png
:align: center
The newer styles with *charmmfsw* or *charmmfsh* in their name replace
energy switching with force switching (fsw) for LJ interactions and
force shifting (fsh) functions for Coulombic interactions
:ref:`(Steinbach) <howto-Steinbach>`
.. grid:: 1 1 2 2
.. grid-item::
.. math::
E_{LJ}(r) = & \begin{cases}
4 \epsilon \sigma^6 \left(\frac{\displaystyle\sigma
^6-r^6}{\displaystyle r^{12}}-\frac{\displaystyle\sigma ^6}{\displaystyle a^6
b^6}+\frac{\displaystyle 1}{\displaystyle a^3 b^3}\right) & r\leq a \\
\frac{\displaystyle 4 \epsilon \sigma^6 \left(\sigma ^6
\left(b^6-r^6\right)^2-b^3 r^6 \left(a^3+b^3\right)
\left(b^3-r^3\right)^2\right)}{\displaystyle b^6 r^{12}
\left(b^6-a^6\right)} & a<r \leq b\\
0, & r>b
\end{cases}\\[.6em]
E_{coul}(r) & = \begin{cases}
C(r) \frac{\displaystyle (b-r)^2}{\displaystyle r b^2}, & r \leq b \\
0, & r > b
\end{cases}
.. grid-item::
.. image:: img/howto_charmmfsw_ELJ.png
:align: center
These styles are used by LAMMPS input scripts generated by
https://charmm-gui.org/ :ref:`(Brooks) <howto-Brooks>`.
.. note::
For CHARMM, newer *charmmfsw* or *charmmfsh* styles were released in
March 2017. We recommend they be used instead of the older *charmm*
styles. See discussion of the differences on the :doc:`pair charmm
<pair_charmm>` and :doc:`dihedral charmm <dihedral_charmm>` doc
pages.
.. note::
The TIP3P water model is strongly recommended for use with the CHARMM
force field. In fact, `"using the SPC model with CHARMM parameters is
a bad idea"
<https://matsci.org/t/using-spc-water-with-charmm-ff/24715>`_ and `"to
enable TIP4P style water in CHARMM, you would have to write a new pair
style"
<https://matsci.org/t/hybrid-pair-styles-for-charmm-and-tip4p-ew/32609>`_
. LAMMPS input scripts generated by Solution Builder on https://charmm-gui.org
use TIP3P molecules for solvation. Any other water model can and
probably will lead to false conclusions.
COMPASS
-------
COMPASS is a general force field for atomistic simulation of common
organic molecules, inorganic small molecules, and polymers which was
developed using ab initio and empirical parameterization techniques
:ref:`(Sun) <howto-Sun>`. See the :doc:`Tools <Tools>` page for the
msi2lmp tool for creating LAMMPS template input and data files from
BIOVIA's Materials Studio files. Please note that the msi2lmp tool is
very old and largely unmaintained, so it does not support all features
of Materials Studio provided force field files, especially additions
during the last decade. You should watch the output carefully and
compare results, where possible. See :ref:`(Sun) <howto-Sun>` for a
description of the COMPASS force field.
These interaction styles listed below compute force field formulas that
are consistent with the COMPASS force field. See each command's
documentation for the formula it computes.
* :doc:`bond_style <bond_class2>` class2
* :doc:`angle_style <angle_class2>` class2
* :doc:`dihedral_style <dihedral_class2>` class2
* :doc:`improper_style <improper_class2>` class2
* :doc:`pair_style <pair_class2>` lj/class2
* :doc:`pair_style <pair_class2>` lj/class2/coul/cut
* :doc:`pair_style <pair_class2>` lj/class2/coul/long
* :doc:`special_bonds <special_bonds>` lj/coul 0 0 1
DREIDING
--------
DREIDING is a generic force field developed by the `Goddard group
<http://www.wag.caltech.edu>`_ at Caltech and is useful for predicting
structures and dynamics of organic, biological and main-group inorganic
molecules. The philosophy in DREIDING is to use general force constants
and geometry parameters based on simple hybridization considerations,
rather than individual force constants and geometric parameters that
depend on the particular combinations of atoms involved in the bond,
angle, or torsion terms. DREIDING has an :doc:`explicit hydrogen bond
term <pair_hbond_dreiding>` to describe interactions involving a
hydrogen atom on very electronegative atoms (N, O, F). Unlike CHARMM
or AMBER, the DREIDING force field has not been parameterized for
considering solvents (like water).
See :ref:`(Mayo) <howto-Mayo>` for a description of the DREIDING force field
The interaction styles listed below compute force field formulas that
are consistent with the DREIDING force field. See each command's
documentation for the formula it computes.
* :doc:`bond_style <bond_harmonic>` harmonic
* :doc:`bond_style <bond_morse>` morse
* :doc:`angle_style <angle_cosine_squared>` cosine/squared
* :doc:`angle_style <angle_harmonic>` harmonic
* :doc:`angle_style <angle_cosine>` cosine
* :doc:`angle_style <angle_cosine_periodic>` cosine/periodic
* :doc:`dihedral_style <dihedral_charmm>` charmm
* :doc:`improper_style <improper_umbrella>` umbrella
* :doc:`pair_style <pair_buck>` buck
* :doc:`pair_style <pair_buck>` buck/coul/cut
* :doc:`pair_style <pair_buck>` buck/coul/long
* :doc:`pair_style <pair_lj>` lj/cut
* :doc:`pair_style <pair_lj_cut_coul>` lj/cut/coul/cut
* :doc:`pair_style <pair_lj_cut_coul>` lj/cut/coul/long
* :doc:`pair_style <pair_hbond_dreiding>` hbond/dreiding/lj
* :doc:`pair_style <pair_hbond_dreiding>` hbond/dreiding/morse
* :doc:`special_bonds <special_bonds>` dreiding
----------
.. _Typelabel2:
**(Gissinger)** J. R. Gissinger, I. Nikiforov, Y. Afshar, B. Waters, M. Choi, D. S. Karls, A. Stukowski, W. Im, H. Heinz, A. Kohlmeyer, and E. B. Tadmor, J Phys Chem B, 128, 3282-3297 (2024).
.. _howto-MacKerell:
**(MacKerell)** MacKerell, Bashford, Bellott, Dunbrack, Evanseck, Field, Fischer, Gao, Guo, Ha, et al (1998). J Phys Chem, 102, 3586 . https://doi.org/10.1021/jp973084f
.. _howto-Cornell:
**(Cornell)** Cornell, Cieplak, Bayly, Gould, Merz, Ferguson, Spellmeyer, Fox, Caldwell, Kollman (1995). JACS 117, 5179-5197. https://doi.org/10.1021/ja00124a002
.. _howto-Steinbach:
**(Steinbach)** Steinbach, Brooks (1994). J Comput Chem, 15, 667. https://doi.org/10.1002/jcc.540150702
.. _howto-Brooks:
**(Brooks)** Brooks, et al (2009). J Comput Chem, 30, 1545. https://onlinelibrary.wiley.com/doi/10.1002/jcc.21287
.. _howto-Sun:
**(Sun)** Sun (1998). J. Phys. Chem. B, 102, 7338-7364. https://doi.org/10.1021/jp980939v
.. _howto-Mayo:
**(Mayo)** Mayo, Olfason, Goddard III (1990). J Phys Chem, 94, 8897-8909. https://doi.org/10.1021/j100389a010