Avoid integer overflow in dispersion correction
authorPaul Bauer <paul.bauer.q@gmail.com>
Fri, 21 Feb 2020 14:49:21 +0000 (15:49 +0100)
committerPaul Bauer <paul.bauer.q@gmail.com>
Tue, 25 Feb 2020 10:34:35 +0000 (11:34 +0100)
Index could become out of range and turn negative.

Fixes #3391

Change-Id: I91a8ea0ab01370fb972c02994e2b996284e1d06a

docs/release-notes/2020/2020.1.rst
src/gromacs/mdlib/dispersioncorrection.cpp

index 278fab5f2d144eb75584b66a5023c35fc19172cd..2c38caf0ecde4bdcff92296393cf82e9948c7d3b 100644 (file)
@@ -112,6 +112,14 @@ steps.
 
 :issue:`3388`
 
+Avoid integer overflow when using dispersioncorrection
+"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+A change in the integer type storing the index meant that the value could overflow
+and turn negative, leading to wrong lookup and unphysical values.
+
+:issue:`3391`
+
 Fixes for ``gmx`` tools
 ^^^^^^^^^^^^^^^^^^^^^^^
 
index 7cbd9e41c7c76cc1cf368a1033193d4840fc4b8c..2d06168cbedf1fee0be36f28aeda75b5811aee51 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2019, by the GROMACS development team, led by
+ * Copyright (c) 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.
@@ -136,10 +136,10 @@ DispersionCorrection::TopologyParams::TopologyParams(const gmx_mtop_t&         m
 
     for (int q = 0; q < (inputrec.efep == efepNO ? 1 : 2); q++)
     {
-        double csix    = 0;
-        double ctwelve = 0;
-        int    npair   = 0;
-        int    nexcl   = 0;
+        double  csix    = 0;
+        double  ctwelve = 0;
+        int64_t npair   = 0;
+        int64_t nexcl   = 0;
         if (!EI_TPI(inputrec.eI))
         {
             numAtomsForDensity_ = mtop.natoms;
@@ -153,9 +153,9 @@ DispersionCorrection::TopologyParams::TopologyParams(const gmx_mtop_t&         m
             {
                 for (int tpj = tpi; tpj < ntp; tpj++)
                 {
-                    const int iCount = typecount[tpi];
-                    const int jCount = typecount[tpj];
-                    int       npair_ij;
+                    const int64_t iCount = typecount[tpi];
+                    const int64_t jCount = typecount[tpj];
+                    int64_t       npair_ij;
                     if (tpi != tpj)
                     {
                         npair_ij = iCount * jCount;
@@ -280,7 +280,7 @@ DispersionCorrection::TopologyParams::TopologyParams(const gmx_mtop_t&         m
         }
         if (debug)
         {
-            fprintf(debug, "Counted %d exclusions\n", nexcl);
+            fprintf(debug, "Counted %" PRId64 " exclusions\n", nexcl);
             fprintf(debug, "Average C6 parameter is: %10g\n", csix);
             fprintf(debug, "Average C12 parameter is: %10g\n", ctwelve);
         }