Fix incorrect AWH free-energy output with multiple walkers
authorMagnus Lundborg <magnus.lundborg@scilifelab.se>
Mon, 14 Dec 2020 11:22:48 +0000 (11:22 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 14 Dec 2020 11:22:48 +0000 (11:22 +0000)
The free-energy was averaged exponentially over walkers using the wrong
sign for the exponent. The meant that the free-energy output was off
for large update sizes. This error was quadratic in the update interval
and the number of walkers. This likely has no measurable effects with
the default AWH free-energy update interval. But the effect is noticeable
with the not uncommon update interval of 1000 steps.

Fixes: #3828

docs/release-notes/2020/2020.5.rst
src/gromacs/awh/biasstate.cpp

index b7966f176e782a171fdc7968c2f1ee7d3df62c66..27e334a923a1219e45ba82998b8731acfc5de04d 100644 (file)
@@ -37,6 +37,20 @@ results.
 
 :issue:`3750`
 
+Fix incorrect AWH free-energies when multiple walkers share a bias
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+The AWH free-energy output was incorrect when multiple walkers shared
+an AWH bias. The error went up quadratically with the free-energy update
+interval, as well as with the number of walkers. The error decreases as
+update size decreases with time. This meant that with default AWH settings
+the error was negligible. With a free-energy update interval of 2 ps,
+we observed an error about equal to the statistical error with 32 walkers
+for a rather fast reaction coordinate. For slower coordinates the error
+will be smaller than the statistical error.
+
+:issue:`3828`
+
 Fixed conserved energy for MTTK
 """""""""""""""""""""""""""""""
 
index 081fb7f13f0a0e40b394ceb917d33b97c45d6cc9..98835d6e6756256f9d6b02527e1d4849e3d57051 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2015,2016,2017,2018,2019, by the GROMACS development team, led by
+ * Copyright (c) 2015,2016,2017,2018,2019,2020, 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.
@@ -163,7 +163,7 @@ void sumPmf(gmx::ArrayRef<PointState> pointState,
     /* Need to temporarily exponentiate the log weights to sum over simulations */
     for (size_t i = 0; i < buffer.size(); i++)
     {
-        buffer[i] = pointState[i].inTargetRegion() ? std::exp(-pointState[i].logPmfSum()) : 0;
+        buffer[i] = pointState[i].inTargetRegion() ? std::exp(pointState[i].logPmfSum()) : 0;
     }
 
     sumOverSimulations(gmx::ArrayRef<double>(buffer), commRecord, multiSimComm);
@@ -174,7 +174,7 @@ void sumPmf(gmx::ArrayRef<PointState> pointState,
     {
         if (pointState[i].inTargetRegion())
         {
-            pointState[i].setLogPmfSum(-std::log(buffer[i] * normFac));
+            pointState[i].setLogPmfSum(std::log(buffer[i] * normFac));
         }
     }
 }