Fix illegal memory access in FE calculations
authorBerk Hess <hess@kth.se>
Wed, 23 Oct 2019 09:50:02 +0000 (11:50 +0200)
committerBerk Hess <hess@kth.se>
Wed, 23 Oct 2019 09:50:02 +0000 (11:50 +0200)
With free-energy calculations not using lambda states, an output
buffer would be accessed one element beyond it's allocated size.

Note that this code should be completely refactored, but not
in a release branch.

Fixes #3173

Change-Id: I677e602ba96c9f64fbf79a626e43c9e590c18bea

docs/release-notes/2019/2019.5.rst
src/gromacs/mdlib/mdebin_bar.cpp

index 90488e4c151ca6e1f962e407c9f347e603dedaf8..a34b88e85b88aea2bb82dfce0139d22a808088c8 100644 (file)
@@ -33,6 +33,16 @@ with gcc as the host side compiler.
 
 :issue:`3120`
 
+Fix out of range memory access with free-energy calculations
+""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+
+With free-energy calculations not using lambda states, an output
+buffer would be accessed one element beyond it's allocated size.
+We don't expect this to have caused incorrect results, but
+a memory checker would complain.
+
+:issue:`3173`
 Fixes for ``gmx`` tools
 ^^^^^^^^^^^^^^^^^^^^^^^
 
index 7b2c030e7978a33437097010bedd3f5c55a07e0b..dcd9b7c1150e69e0af101952f5987d66d8ce934b 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.
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/smalloc.h"
 
+/* Number of entries in subblock_d preceding the lambda components */
+constexpr int c_subblockDNumPreEntries = 5;
+/* Number of entries in subblock_i preceding the lambda components */
+constexpr int c_subblockINumPreEntries = 2;
+
 /* reset the delta_h list to prepare it for new values */
 static void mde_delta_h_reset(t_mde_delta_h *dh)
 {
@@ -428,8 +433,8 @@ void mde_delta_h_coll_init(t_mde_delta_h_coll *dhc, const t_inputrec *ir)
         dhc->lambda_index             = -1;
     }
     /* allocate metadata subblocks */
-    snew(dhc->subblock_d, 5 + dhc->n_lambda_vec);
-    snew(dhc->subblock_i, 1 + dhc->n_lambda_vec);
+    snew(dhc->subblock_d, c_subblockDNumPreEntries + dhc->n_lambda_vec);
+    snew(dhc->subblock_i, c_subblockINumPreEntries + dhc->n_lambda_vec);
 
     /* now decide which data to write out */
     dhc->nlambda     = 0;
@@ -663,11 +668,11 @@ void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
     {
         for (i = 0; i < dhc->n_lambda_vec; i++)
         {
-            dhc->subblock_d[5+i] = dhc->native_lambda_vec[i];
+            dhc->subblock_d[c_subblockDNumPreEntries + i] = dhc->native_lambda_vec[i];
         }
     }
     blk->id          = enxDHCOLL;
-    blk->sub[0].nr   = 5 + dhc->n_lambda_vec;
+    blk->sub[0].nr   = c_subblockDNumPreEntries + dhc->n_lambda_vec;
     blk->sub[0].type = xdr_datatype_double;
     blk->sub[0].dval = dhc->subblock_d;
 
@@ -678,9 +683,9 @@ void mde_delta_h_coll_handle_block(t_mde_delta_h_coll *dhc,
         dhc->subblock_i[1] = dhc->n_lambda_vec;
         for (i = 0; i < dhc->n_lambda_vec; i++)
         {
-            dhc->subblock_i[i+2] = dhc->native_lambda_components[i];
+            dhc->subblock_i[c_subblockINumPreEntries + i] = dhc->native_lambda_components[i];
         }
-        blk->sub[1].nr   = 2 + dhc->n_lambda_vec;
+        blk->sub[1].nr   = c_subblockINumPreEntries + dhc->n_lambda_vec;
         blk->sub[1].type = xdr_datatype_int;
         blk->sub[1].ival = dhc->subblock_i;
     }