Merge origin/release-2019 into release-2020
authorPaul Bauer <paul.bauer.q@gmail.com>
Fri, 27 Dec 2019 12:58:41 +0000 (13:58 +0100)
committerPaul Bauer <paul.bauer.q@gmail.com>
Fri, 27 Dec 2019 12:58:41 +0000 (13:58 +0100)
Resolved Conflicts:
cmake/FindHwloc.cmake
cmake/gmxVersionInfo.cmake
docs/CMakeLists.txt
src/gromacs/gmxana/anadih.cpp
src/gromacs/gmxana/gmx_angle.cpp
src/gromacs/mdlib/forcerec.cpp

Change-Id: Ie4408c932a712d4f9d1f42117a512e0474535df5

docs/CMakeLists.txt
docs/release-notes/2019/2019.5.rst
docs/release-notes/2019/2019.6.rst [new file with mode: 0644]
docs/release-notes/index.rst
src/gromacs/fileio/pdbio.cpp
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/cmat.cpp
src/gromacs/gmxana/gmx_angle.cpp
src/gromacs/mdlib/forcerec.cpp

index b13a72c54dffd727f7ca89f71dd39ece142ff06d..32636d3dd1e1aa1305379bd4adf267b86d2abc2f 100644 (file)
@@ -370,6 +370,7 @@ if (SPHINX_FOUND)
         release-notes/2020/major/deprecated-functionality.rst
         release-notes/2020/major/portability.rst
         release-notes/2020/major/miscellaneous.rst
+        release-notes/2019/2019.6.rst
         release-notes/2019/2019.5.rst
         release-notes/2019/2019.4.rst
         release-notes/2019/2019.3.rst
index 7a757cc08eba888fd77e4b8f8ea464c89ab81e70..f526dc3ee1ecd0aad6c83978865e9c0eb196fa4b 100644 (file)
@@ -1,10 +1,10 @@
 GROMACS 2019.5 release notes
 ----------------------------
 
-This version was released on TODO, 2019. These release notes
+This version was released on December 23rd, 2019. These release notes
 document the changes that have taken place in GROMACS since the
 previous 2019.4 version, to fix known issues. It also incorporates all
-fixes made in version 2018.7 and earlier, which you can find described
+fixes made in version 2018.8 and earlier, which you can find described
 in the :ref:`release-notes`.
 
 .. Note to developers!
@@ -80,9 +80,44 @@ Remove assertion failure with AWH when not using the initial stage
 Fixes for ``gmx`` tools
 ^^^^^^^^^^^^^^^^^^^^^^^
 
+Make histogram output clearer
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+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`
+
+
+Remove problematic output of gmx angle tool
+"""""""""""""""""""""""""""""""""""""""""""
+
+It could happen that the calculation of the standard deviation
+for angles caused a divide by zero error for empty populations.
+Because this standard deviation was meaningless, it has been
+removed.
+
+:issue:`3206`
+
 Fixes that affect portability
 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
 
+Check that libhwloc headers and runtime match
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+It could happen that the libhwloc headers and library detection would
+lead to a mismatch at compile or runtime that could cause cryptic
+crashes while using mdrun.
+
+:issue:`3200`
+
 Miscellaneous
 ^^^^^^^^^^^^^
 
@@ -101,10 +136,16 @@ and expect whitespace separation will continue to work in all cases.
 :issue:`3176`
 
 Fix duplicate PDB CONECT record output
-^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
+""""""""""""""""""""""""""""""""""""""
 
 PDB "CONECT" record output was duplicated in some instances. Since |Gromacs| does
 not use these anywhere, analysis was not affected. The behavior is now fixed.
 
 :issue:`3206`
 
+Fix performance issue with bonded interactions in wrong GPU stream
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+This could lead to a significant loss in performance.
+
+:issue:`3241`
diff --git a/docs/release-notes/2019/2019.6.rst b/docs/release-notes/2019/2019.6.rst
new file mode 100644 (file)
index 0000000..68cb9a8
--- /dev/null
@@ -0,0 +1,28 @@
+GROMACS 2019.6 release notes
+----------------------------
+
+This version was released on TODO, 2020. These release notes
+document the changes that have taken place in GROMACS since the
+previous 2019.5 version, to fix known issues. It also incorporates all
+fixes made in version 2018.8 and earlier, which you can find described
+in the :ref:`release-notes`.
+
+.. Note to developers!
+   Please use """"""" to underline the individual entries for fixed issues in the subfolders,
+   otherwise the formatting on the webpage is messed up.
+   Also, please use the syntax :issue:`number` to reference issues on redmine, without the
+   a space between the colon and number!
+
+Fixes where mdrun could behave incorrectly
+^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
+
+Fixes for ``gmx`` tools
+^^^^^^^^^^^^^^^^^^^^^^^
+
+Fixes that affect portability
+^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
+
+
+Miscellaneous
+^^^^^^^^^^^^^
+
index bfb9ab19afa60d3101646d7a99b290f10595e0ea..16cb7ea8ac75c86ddba6de2e545b61ec5b65f3e8 100644 (file)
@@ -49,6 +49,7 @@ Patch releases
 .. toctree::
    :maxdepth: 1
 
+   2019/2019.6
    2019/2019.5
    2019/2019.4
    2019/2019.3
index 3466d8d4f1114644918dbe3c9fba228062ccc701..f431614cc8a512e39d767b6fd97c9a9d9726a755 100644 (file)
@@ -754,14 +754,14 @@ static void gmx_conect_addline(gmx_conect_t* con, char* line)
 {
     int n, ai, aj;
 
-    std::string form2  = "%%*s";
-    std::string format = form2 + "%%d";
+    std::string form2  = "%*s";
+    std::string format = form2 + "%d";
     if (sscanf(line, format.c_str(), &ai) == 1)
     {
         do
         {
             form2 += "%*s";
-            format = form2 + "%%d";
+            format = form2 + "%d";
             n      = sscanf(line, format.c_str(), &aj);
             if (n == 1)
             {
index b76f99dc1aed9e2bf31e3bcd5ccedf66b7595668..bfe16bede134eead6646469ce52e79c2d1819237 100644 (file)
@@ -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"
@@ -833,7 +834,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;
@@ -946,9 +947,16 @@ void read_ang_dih(const char*             trj_fn,
         }
 
         /* Average angles */
-        aa = 0;
+        double aa = 0;
         for (i = 0; (i < nangles); 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,
@@ -960,15 +968,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;
             }
 
@@ -994,14 +994,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..3ee9a99
--- /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 cde327414c01eff198371fcab431a11d423d2111..07368d9d41b2c5ddff94daa131bd7edcd54e5d96 100644 (file)
@@ -227,7 +227,7 @@ void low_rmsd_dist(const char* fn, real maxrms, int nn, real** mat, const gmx_ou
         }
     }
 
-    fp = xvgropen(fn, "RMS Distribution", "RMS (nm)", "a.u.", oenv);
+    fp = xvgropen(fn, "RMS Distribution", "RMS (nm)", "counts", oenv);
     for (i = 0; (i < 101); i++)
     {
         fprintf(fp, "%10g  %10d\n", i / fac, histo[i]);
index 0a92c634d86bb4a2eb91a7e44024e3d772db76c5..e99bee5628e958d5ed9574e09c120adb03d3d753 100644 (file)
@@ -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"
@@ -174,7 +175,6 @@ int gmx_g_angle(int argc, char* argv[])
     gmx_bool      bAver, bRb, bPeriodic, 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*/
@@ -403,19 +403,40 @@ int gmx_g_angle(int argc, char* argv[])
     for (first = 0; (first < maxangstat - 1) && (angstat[first + 1] == 0); first++) {}
     for (last = maxangstat - 1; (last > 0) && (angstat[last - 1] == 0); last--) {}
 
-    aver = aver2 = 0;
-    for (i = 0; (i < nframes); i++)
-    {
-        aver += RAD2DEG * aver_angle[i];
-        aver2 += gmx::square(RAD2DEG * aver_angle[i]);
+    double aver = 0;
+    printf("Found points in the range from %d to %d (max %d)\n", first, last, maxangstat);
+    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;
+            }
+        }
+    }
+    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;
+        }
     }
     aver /= nframes;
-    aver2 /= nframes;
-    aversig = std::sqrt(aver2 - gmx::square(aver));
-    printf("Found points in the range from %d to %d (max %d)\n", first, last, maxangstat);
-    printf(" < angle >  = %g\n", aver);
-    printf("< angle^2 > = %g\n", aver2);
-    printf("Std. Dev.   = %g\n", aversig);
+    double aversig = correctRadianAngleRange(aver);
+    aversig *= RAD2DEG;
+    aver *= RAD2DEG;
+    printf(" < angle >  = %g\n", aversig);
 
     if (mult == 3)
     {
index b702a4a2552b91a025493e0720474fe4fc578451..c8a42aac5896ea3b95e02b2d792eb07133607c95 100644 (file)
@@ -1438,7 +1438,7 @@ void init_forcerec(FILE*                            fp,
 
         if (useGpuForBonded)
         {
-            auto stream = DOMAINDECOMP(cr)
+            auto stream = havePPDomainDecomposition(cr)
                                   ? Nbnxm::gpu_get_command_stream(
                                             fr->nbv->gpu_nbv, gmx::InteractionLocality::NonLocal)
                                   : Nbnxm::gpu_get_command_stream(fr->nbv->gpu_nbv,