Add AWH references
authorBerk Hess <hess@kth.se>
Thu, 27 May 2021 18:45:19 +0000 (18:45 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 27 May 2021 18:45:19 +0000 (18:45 +0000)
docs/reference-manual/references.rst
docs/reference-manual/special/awh.rst
docs/reference-manual/special/free-energy-implementation.rst
src/gromacs/applied_forces/awh/awh.cpp
src/gromacs/utility/pleasecite.cpp

index eb9099d164c147512c505e9aa96788c945ee2474..0d2aa5f93a7f9245f72362d1c0d8b4ce924f5f46 100644 (file)
@@ -2600,9 +2600,42 @@ structures into cryoelectron microscopy maps using biased molecular dynamics sim
 
 .. _refBernetti2020:
 
-:sup:`184` Bernetti, M. and Bussi, G.,
+:sup:`184` Bernetti, M. and Bussi G.,
 "Pressure control using stochastic cell rescaling", *J. Chem. Phys.*, **153**, 114107 (2020).
 
+.. raw:: html
+
+   </div>
+
+  <div id="ref-Lidmar2012">
+
+.. _refLidmar2012:
+
+:sup:`185` Lidmar J.,
+"Improving the efficiency of extended ensemble simulations: The accelerated weight histogram method", *Phys. Rev. E*, **85**, 0256708 (2012).
+
+.. raw:: html
+
+   </div>
+
+   <div id="ref-Lindahl2018">
+
+.. _refLindahl2018:
+
+:sup:`186` Lindahl V., Lidmar J. and Hess B.,
+"Riemann metric approach to optimal sampling of multidimensional free-energy landscapes", *Phys. Rev. E*, **98**, 023312 (2018).
+
+.. raw:: html
+
+   </div>
+
+   <div id="ref-LundBorg2021">
+
+.. _refLundborg2021:
+
+:sup:`187` Lundborg M., Lidmar J. and Hess B.,
+"The accelerated weight histogram method for alchemical free energy calculations", *J. Chem. Phys.*, **154**, 204103 (2021).
+
 .. raw:: html
 
    </div>
index f5b2fbf225608561a66d6e6aecbf94d2b4a20315..f9e30ea273aa056c9762b5d53fed21eb996c8937 100644 (file)
@@ -3,7 +3,7 @@
 Adaptive biasing with AWH
 -------------------------
 
-The accelerated weight histogram method
+The accelerated weight histogram method :ref:`185 <refLidmar2012>`
 :ref:`137 <reflindahl2014accelerated>` calculates the PMF along a reaction coordinate by adding
 an adaptively determined biasing potential. AWH flattens free energy
 barriers along the reaction coordinate by applying a history-dependent
@@ -64,7 +64,7 @@ determined accurately. Thus, AWH adaptively calculates
 toward :math:`\rho(\lambda)`.
 
 It is also possible to directly control the :math:`\lambda` state
-of, e.g., alchemical free energy perturbations. In that case there is no harmonic
+of, e.g., alchemical free energy perturbations :ref:`187 <reflundborg2021>`. In that case there is no harmonic
 potential and :math:`\lambda` changes in discrete steps along the reaction coordinate
 depending on the biased free energy difference between the :math:`\lambda` states.
 N.b., it is not yet possible to use AWH in combination with perturbed masses or
@@ -566,7 +566,7 @@ centered at :math:`\lambda` and
 is the deviation of the force. The factors :math:`\omega(\lambda|x(t))`,
 see :eq:`Eq %s <eqawhomega>`, reweight the samples.
 :math:`\eta_{\mu\nu}(\lambda)` is a friction
-tensor \ :ref:`144 <refsivak2012thermodynamic>`. Its matrix elements are inversely proportional to local
+tensor :ref:`186 <reflindahl2018>` and :ref:`144 <refsivak2012thermodynamic>`. Its matrix elements are inversely proportional to local
 diffusion coefficients. A measure of sampling (in)efficiency at each
 :math:`\lambda` is given by
 
index 3233975edd29961ddf89df0d8aa7279844df6733..6de27ad64ff2845aa82c1ece51168be2b6861d29 100644 (file)
@@ -119,7 +119,7 @@ force. Each method has its limitations, which are listed below.
    molecules.
 
 -  **AWH code:** currently acts on coordinates provided by the pull
-   code.
+   code or the free-energy lambda parameter.
 
 -  **free-energy code with harmonic bonds or constraints:** between
    single atoms.
index 0f3c8c28b4c55c6e5436eaebb0bde15c06e469b5..02f030b08c08c5d84ce59390ac535c36887ec51d 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2015,2016,2017,2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2015,2016,2017,2018,2019,2020,2021, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -107,17 +107,18 @@ struct BiasCoupledToSystem
     /* Here AWH can be extended to work on other coordinates than pull. */
 };
 
-/*! \brief Checks whether any dimension uses pulling as a coordinate provider.
+/*! \brief Checks whether any dimension uses the given coordinate provider type.
  *
  * \param[in] awhBiasParams The bias params to check.
- * \returns true if any dimension of the bias is linked to pulling.
+ * \param[in] awhCoordProvider The type of coordinate provider
+ * \returns true if any dimension of the bias is linked to the given provider
  */
-static bool anyDimUsesPull(const AwhBiasParams& awhBiasParams)
+static bool anyDimUsesProvider(const AwhBiasParams& awhBiasParams, const int awhCoordProvider)
 {
     for (int d = 0; d < awhBiasParams.ndim; d++)
     {
         const AwhDimParams& awhDimParams = awhBiasParams.dimParams[d];
-        if (awhDimParams.eCoordProvider == eawhcoordproviderPULL)
+        if (awhDimParams.eCoordProvider == awhCoordProvider)
         {
             return true;
         }
@@ -125,17 +126,18 @@ static bool anyDimUsesPull(const AwhBiasParams& awhBiasParams)
     return false;
 }
 
-/*! \brief Checks whether any dimension uses pulling as a coordinate provider.
+/*! \brief Checks whether any dimension uses the given coordinate provider type.
  *
  * \param[in] awhParams The AWH params to check.
- * \returns true if any dimension of awh is linked to pulling.
+ * \param[in] awhCoordProvider The type of coordinate provider
+ * \returns true if any dimension of awh is linked to the given provider type.
  */
-static bool anyDimUsesPull(const AwhParams& awhParams)
+static bool anyDimUsesProvider(const AwhParams& awhParams, const int awhCoordProvider)
 {
     for (int k = 0; k < awhParams.numBias; k++)
     {
         const AwhBiasParams& awhBiasParams = awhParams.awhBiasParams[k];
-        if (anyDimUsesPull(awhBiasParams))
+        if (anyDimUsesProvider(awhBiasParams, awhCoordProvider))
         {
             return true;
         }
@@ -188,7 +190,7 @@ Awh::Awh(FILE*                 fplog,
     numFepLambdaStates_(numFepLambdaStates),
     fepLambdaState_(fepLambdaState)
 {
-    if (anyDimUsesPull(awhParams))
+    if (anyDimUsesProvider(awhParams, eawhcoordproviderPULL))
     {
         GMX_RELEASE_ASSERT(inputRecord.pull != nullptr, "With AWH we should have pull parameters");
         GMX_RELEASE_ASSERT(pull_work != nullptr,
@@ -198,6 +200,11 @@ Awh::Awh(FILE*                 fplog,
     if (fplog != nullptr)
     {
         please_cite(fplog, "Lindahl2014");
+
+        if (anyDimUsesProvider(awhParams, eawhcoordproviderFREE_ENERGY_LAMBDA))
+        {
+            please_cite(fplog, "Lundborg2021");
+        }
     }
 
     if (haveBiasSharingWithinSimulation(awhParams))
@@ -461,7 +468,8 @@ const char* Awh::externalPotentialString()
 
 void Awh::registerAwhWithPull(const AwhParams& awhParams, pull_t* pull_work)
 {
-    GMX_RELEASE_ASSERT(!anyDimUsesPull(awhParams) || pull_work, "Need a valid pull object");
+    GMX_RELEASE_ASSERT(!anyDimUsesProvider(awhParams, eawhcoordproviderPULL) || pull_work,
+                       "Need a valid pull object");
 
     for (int k = 0; k < awhParams.numBias; k++)
     {
index c5985c48fdbc7e1dd791a349233537951f17a1e4..e333ca5d729731430acb1fbe507493832e356e51 100644 (file)
@@ -4,7 +4,7 @@
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,2020,2021, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -296,6 +296,9 @@ void please_cite(FILE* fp, const char* key)
         { "Bernetti2020", "M. Bernetti, G. Bussi",
           "Pressure control using stochastic cell rescaling", "J. Chem. Phys.", 153, 2020,
           "114107" },
+        { "Lundborg2021", "M. Lundborg, J. Lidmar, B. Hess",
+          "The accelerated weight histogram method for alchemical free energy calculations",
+          "J. Chem. Phys.", 154, 2021, "204103" },
     };
 #define NSTR static_cast<int>(asize(citedb))