Fix torsion angle calculation
authorPaul Bauer <paul.bauer.q@gmail.com>
Fri, 6 Dec 2019 11:18:23 +0000 (12:18 +0100)
committerChristian Blau <cblau@gerrit.gromacs.org>
Mon, 16 Dec 2019 09:30:11 +0000 (10:30 +0100)
Fix originally contributed by Boris Timofeev.

Fixes #3225

Change-Id: I4862679e4aeae514736df53f0eff44bf3b85f928

docs/release-notes/2019/2019.5.rst
src/gromacs/gmxana/anadih.cpp
src/gromacs/gmxana/angle_correction.cpp [new file with mode: 0644]
src/gromacs/gmxana/angle_correction.h [new file with mode: 0644]
src/gromacs/gmxana/gmx_angle.cpp

index 5f2530099e43a29ab40e2421df7722c2fc8d05f3..7bc2700155878760e18948fbe07a2364610268b3 100644 (file)
@@ -87,6 +87,15 @@ Output gave number of events included in histogram bar as *a.u.*,
 which was not clear for users.
 
 
+Fix dihedral angle calculation near 180 degree boundary
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+The analysis tools could incorrectly calculate properties of torsion angles and their averages
+when close to the -180 or 180 degree boundary.
+
+:issue:`3225`
+
+
 Fixes that affect portability
 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
 
index cf7869ced0ee85cec290c293dfa34323cc697dc8..c266e8961dc0078d39c2c4a2417b5fd0b8ac4c95 100644 (file)
@@ -3,7 +3,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,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -45,6 +45,7 @@
 #include "gromacs/fileio/confio.h"
 #include "gromacs/fileio/trxio.h"
 #include "gromacs/fileio/xvgr.h"
+#include "gromacs/gmxana/angle_correction.h"
 #include "gromacs/gmxana/gstat.h"
 #include "gromacs/listed-forces/bonded.h"
 #include "gromacs/math/functions.h"
@@ -821,7 +822,7 @@ void read_ang_dih(const char *trj_fn,
     t_trxstatus  *status;
     int           i, angind, total, teller;
     int           nangles, n_alloc;
-    real          t, fraction, pifac, aa, angle;
+    real          t, fraction, pifac, angle;
     real         *angles[2];
     matrix        box;
     rvec         *x;
@@ -934,10 +935,17 @@ void read_ang_dih(const char *trj_fn,
         }
 
         /* Average angles */
-        aa = 0;
+        double aa = 0;
         for (i = 0; (i < nangles); i++)
         {
-            aa = aa+angles[cur][i];
+            if (!bAngles && i > 0)
+            {
+                real diffa = angles[cur][i] - angles[cur][i-1];
+                diffa          = correctRadianAngleRange(diffa);
+                angles[cur][i] = angles[cur][i-1] + diffa;
+            }
+
+            aa = aa + angles[cur][i];
 
             /* angle in rad / 2Pi * max determines bin. bins go from 0 to maxangstat,
                even though scale goes from -pi to pi (dihedral) or -pi/2 to pi/2
@@ -948,15 +956,7 @@ void read_ang_dih(const char *trj_fn,
             angle = angles[cur][i];
             if (!bAngles)
             {
-                while (angle < -M_PI)
-                {
-                    angle += 2*M_PI;
-                }
-                while (angle >= M_PI)
-                {
-                    angle -= 2*M_PI;
-                }
-
+                angle  = correctRadianAngleRange(angle);
                 angle += M_PI;
             }
 
@@ -983,14 +983,22 @@ void read_ang_dih(const char *trj_fn,
         }
 
         /* average over all angles */
-        (*aver_angle)[teller] = (aa/nangles);
+        aa                    = correctRadianAngleRange(aa/nangles);
+        (*aver_angle)[teller] = (aa);
 
         /* this copies all current dih. angles to dih[i], teller is frame */
         if (bSaveAll)
         {
             for (i = 0; i < nangles; i++)
             {
-                dih[i][teller] = angles[cur][i];
+                if (!bAngles)
+                {
+                    dih[i][teller] = correctRadianAngleRange(angles[cur][i]);
+                }
+                else
+                {
+                    dih[i][teller] = angles[cur][i];
+                }
             }
         }
 
diff --git a/src/gromacs/gmxana/angle_correction.cpp b/src/gromacs/gmxana/angle_correction.cpp
new file mode 100644 (file)
index 0000000..aaa687f
--- /dev/null
@@ -0,0 +1,56 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2019, 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+#include "gmxpre.h"
+
+#include "angle_correction.h"
+
+#include <algorithm>
+
+#include "gromacs/math/units.h"
+#include "gromacs/math/vec.h"
+
+real correctRadianAngleRange(const real angle)
+{
+    real correctedAngle = angle;
+    while (correctedAngle < -M_PI)
+    {
+        correctedAngle += 2*M_PI;
+    }
+    while (correctedAngle >= M_PI)
+    {
+        correctedAngle -= 2*M_PI;
+    }
+    return correctedAngle;
+}
diff --git a/src/gromacs/gmxana/angle_correction.h b/src/gromacs/gmxana/angle_correction.h
new file mode 100644 (file)
index 0000000..5c68c4d
--- /dev/null
@@ -0,0 +1,48 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2019, 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+#ifndef GMX_GMXANA_ANGLE_CORRECTION_H
+#define GMX_GMXANA_ANGLE_CORRECTION_H
+
+#include "gromacs/math/units.h"
+
+/*! \brief
+ * Return the angle in radians after correcting to be within range of -PI < \p angle < PI.
+ *
+ * \param[in] angle The angle in radians.
+ * \returns Angle within range.
+ */
+real correctRadianAngleRange(real angle);
+
+#endif
index 1b56598e0fd06351050efc0a969f86d1e73d6302..4fa39ad4f129962be4121bc4b7ee152c71669575 100644 (file)
@@ -3,7 +3,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,2018, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2017,2018,2019, 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.
@@ -46,6 +46,7 @@
 #include "gromacs/correlationfunctions/autocorr.h"
 #include "gromacs/fileio/trrio.h"
 #include "gromacs/fileio/xvgr.h"
+#include "gromacs/gmxana/angle_correction.h"
 #include "gromacs/gmxana/gmx_ana.h"
 #include "gromacs/gmxana/gstat.h"
 #include "gromacs/math/functions.h"
@@ -164,7 +165,6 @@ int gmx_g_angle(int argc, char *argv[])
                        bFrac,               /* calculate fraction too?  */
                        bTrans,              /* worry about transtions too? */
                        bCorr;               /* correlation function ? */
-    real              aver, aver2, aversig; /* fraction trans dihedrals */
     double            tfrac = 0;
     char              title[256];
     real            **dih = nullptr; /* mega array with all dih. angles at all times*/
@@ -411,20 +411,55 @@ int gmx_g_angle(int argc, char *argv[])
         ;
     }
 
-    aver = aver2 = 0;
-    for (i = 0; (i < nframes); i++)
-    {
-        aver  += RAD2DEG*aver_angle[i];
-        aver2 += gmx::square(RAD2DEG*aver_angle[i]);
-    }
-    aver   /= nframes;
-    aver2  /= nframes;
-    aversig = std::sqrt(aver2-gmx::square(aver));
+    double aver  = 0;
+    double aver2 = 0;
     printf("Found points in the range from %d to %d (max %d)\n",
            first, last, maxangstat);
-    printf(" < angle >  = %g\n", aver);
+    if (bTrans || bCorr  || bALL || opt2bSet("-or", NFILE, fnm))
+    {   /* It's better to re-calculate Std. Dev per sample */
+        real b_aver = aver_angle[0];
+        real b      = dih[0][0];
+        real delta;
+        for (int i = 0; (i < nframes); i++)
+        {
+            delta   = correctRadianAngleRange(aver_angle[i] - b_aver);
+            b_aver += delta;
+            aver   += b_aver;
+            for (int j = 0; (j < nangles); j++)
+            {
+                delta  = correctRadianAngleRange(dih[j][i] - b);
+                b     += delta;
+                aver2 += gmx::square(b);
+            }
+        }
+        aver2 /= nangles;
+    }
+    else
+    {   /* Incorrect  for Std. Dev. */
+        real delta, b_aver = aver_angle[0];
+        for (i = 0; (i < nframes); i++)
+        {
+            delta   = correctRadianAngleRange(aver_angle[i] - b_aver);
+            b_aver += delta;
+            aver   += b_aver;
+            aver2  += gmx::square(b_aver);
+        }
+    }
+    aver  /= nframes;
+    aver2 /= nframes/(RAD2DEG * RAD2DEG);
+    double aversig = correctRadianAngleRange(aver);
+    aversig *= RAD2DEG;
+    aver    *= RAD2DEG;
+    printf(" < angle >  = %g\n", aversig);
     printf("< angle^2 > = %g\n", aver2);
-    printf("Std. Dev.   = %g\n", aversig);
+    if ((aversig = aver2-gmx::square(aver)) > 0)
+    {
+        printf("Std. Dev.   = %g\n", std::sqrt(aversig));
+    }
+    else
+    {
+        printf("Std. Dev. - can't be calculated\n");
+    }
 
     if (mult == 3)
     {