Change default MTS force group selection
authorBerk Hess <hess@kth.se>
Tue, 19 Jan 2021 16:38:19 +0000 (16:38 +0000)
committerArtem Zhmurov <zhmurov@gmail.com>
Tue, 19 Jan 2021 16:38:19 +0000 (16:38 +0000)
Because of instabilities of aliphatic protons with the Amber
force fields, the default force group selection for multiple
time stepping has been changed to include longrange-nonbonded only.

Refs #3804

docs/user-guide/mdp-options.rst
src/gromacs/gmxpreprocess/readir.cpp

index 6a6fd20750a62b2e24a85e3d5668369e9aa0ff5f..74f140063075e9c2fe814edd9750eeada5960c83 100644 (file)
@@ -252,18 +252,15 @@ Run control
 
 .. mdp:: mts-level2-forces
 
-   (longrange-nonbonded nonbonded pair dihedral)
-   A list of force groups that will be evaluated only every
+   (longrange-nonbonded)
+   A list of one or more force groups that will be evaluated only every
    :mdp:`mts-level2-factor` steps. Supported entries are:
    ``longrange-nonbonded``, ``nonbonded``, ``pair``, ``dihedral``, ``angle``,
    ``pull`` and ``awh``. With ``pair`` the listed pair forces (such as 1-4)
    are selected. With ``dihedral`` all dihedrals are selected, including cmap.
    All other forces, including all restraints, are evaluated and
    integrated every step. When PME or Ewald is used for electrostatics
-   and/or LJ interactions, ``longrange-nonbonded`` has to be entered here.
-   The default value should work well for most standard atomistic simulations
-   and in particular for replacing virtual site treatment for increasing
-   the time step.
+   and/or LJ interactions, ``longrange-nonbonded`` can not be omitted here.
 
 .. mdp:: mts-level2-factor
 
index ded12b74afe43836b64be410296aa1103856d289..f7cb14caecc42010d5d40e86f288caa5fbd33207 100644 (file)
@@ -1996,8 +1996,7 @@ void get_ir(const char*     mdparin,
         opts->numMtsLevels = get_eint(&inp, "mts-levels", 2, wi);
         ir->mtsLevels.resize(2);
         gmx::MtsLevel& mtsLevel = ir->mtsLevels[1];
-        opts->mtsLevel2Forces   = setStringEntry(&inp, "mts-level2-forces",
-                                               "longrange-nonbonded nonbonded pair dihedral");
+        opts->mtsLevel2Forces   = setStringEntry(&inp, "mts-level2-forces", "longrange-nonbonded");
         mtsLevel.stepFactor     = get_eint(&inp, "mts-level2-factor", 2, wi);
 
         // We clear after reading without dynamics to not force the user to remove MTS mdp options