Fix UB vector usage
authorRoland Schulz <roland.schulz@intel.com>
Tue, 26 Mar 2019 00:24:39 +0000 (17:24 -0700)
committerRoland Schulz <roland.schulz@intel.com>
Wed, 27 Mar 2019 15:46:02 +0000 (08:46 -0700)
Fixes a buffer overflow for c_simdBestPairAlignment=2.

Found with _LIBCPP_DEBUG=1.

Change-Id: Ib21ca875244673b27748a01373e7fc10252a7c44

src/gromacs/mdlib/nbnxn_atomdata.cpp

index 08eab170922cae82a95e5e398e27708dfdf48fa3..69052128d1d043dc1e569aad4e1732bce030d915 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2012,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.
@@ -345,8 +345,11 @@ static void set_lj_parameter_data(nbnxn_atomdata_t *nbat, gmx_bool bSIMD)
             {
                 nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+0] = nbat->nbfp[(i*nt+j)*2+0];
                 nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+1] = nbat->nbfp[(i*nt+j)*2+1];
-                nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+2] = 0;
-                nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+3] = 0;
+                if (c_simdBestPairAlignment > 2)
+                {
+                    nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+2] = 0;
+                    nbat->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+3] = 0;
+                }
             }
         }
 #endif