Fix UB vector usage
authorRoland Schulz <roland.schulz@intel.com>
Tue, 26 Mar 2019 00:24:39 +0000 (17:24 -0700)
committerPaul Bauer <paul.bauer.q@gmail.com>
Thu, 28 Mar 2019 10:10:17 +0000 (11:10 +0100)
It is UB to
- increment past end.
- decrement end iterator for empty vector.
- use operator[] on end iterator.

Also fixes a buffer overflow for c_simdBestPairAlignment=2.

All found with _LIBCPP_DEBUG=1.

Change-Id: Ib21ca875244673b27748a01373e7fc10252a7c44

src/gromacs/gmxpreprocess/resall.cpp
src/gromacs/gmxpreprocess/toppush.cpp
src/gromacs/mdlib/lincs.cpp
src/gromacs/nbnxm/atomdata.cpp

index fac51b47a22947911ab3df0c0b24035a819a98c9..df04c52f831ae06ef0b6a79d19d1633638bca4d7 100644 (file)
@@ -376,7 +376,7 @@ void readResidueDatabase(const std::string &rrdb, std::vector<PreprocessResidue>
         print_resall_header(stderr, gmx::arrayRefFromArray(&header_settings, 1));
     }
     /* We don't know the current size of rrtp, but simply realloc immediately */
-    auto oldArrayEnd = rtpDBEntry->end() - 1;
+    auto oldArrayEnd = rtpDBEntry->end();
     while (!feof(in))
     {
         /* Initialise rtp entry structure */
index c4d7b16420f34fe6d836cd8edeae903a668167c8..8f4c46d437557c1ce2d0d5a992e1a43f81c213e2 100644 (file)
@@ -1786,9 +1786,16 @@ defaultInteractionTypeParameters(int ftype, gmx::ArrayRef<InteractionTypeParamet
              * The rule in that case is that additional matches
              * HAVE to be on adjacent lines!
              */
-            bool bSame = true;
+            bool       bSame = true;
+            //Advance iterator (like std::advance) without incrementing past end (UB)
+            const auto safeAdvance = [](auto &it, auto n, auto end) {
+                    it = end-it > n ? it+n : end;
+                };
             /* Continue from current iterator position */
-            for (auto nextPos = prevPos + 2; (nextPos < bt[ftype].interactionTypes.end()) && bSame; nextPos += 2)
+            auto       nextPos = prevPos;
+            const auto endIter = bt[ftype].interactionTypes.end();
+            safeAdvance(nextPos, 2, endIter);
+            for (; nextPos < endIter && bSame; safeAdvance(nextPos, 2, endIter))
             {
                 bSame = (prevPos->ai() == nextPos->ai() && prevPos->aj() == nextPos->aj() && prevPos->ak() == nextPos->ak() && prevPos->al() == nextPos->al());
                 if (bSame)
index 185d968f57da4926708bcb1716ec6f8a2ad1e6b4..169591e0e8c9ffc7b57eb331ee7621c4410c1107 100644 (file)
@@ -1877,7 +1877,7 @@ static void set_matrix_indices(Lincs                *li,
         if (bSortMatrix)
         {
             /* Order the blbnb matrix to optimize memory access */
-            std::sort(&(li->blbnb[li->blnr[b]]), &(li->blbnb[li->blnr[b+1]]));
+            std::sort(li->blbnb.begin()+li->blnr[b], li->blbnb.begin()+li->blnr[b+1]);
         }
     }
 }
index ac1839e3ec29106172556d63b11ce2b0fc859e91..38c09fcbdf16efd36894b7a4ac39f3f65c44a282 100644 (file)
@@ -273,8 +273,11 @@ static void set_lj_parameter_data(nbnxn_atomdata_t::Params *params, gmx_bool bSI
             {
                 params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+0] = params->nbfp[(i*nt+j)*2+0];
                 params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+1] = params->nbfp[(i*nt+j)*2+1];
-                params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+2] = 0;
-                params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+3] = 0;
+                if (c_simdBestPairAlignment > 2)
+                {
+                    params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+2] = 0;
+                    params->nbfp_aligned[(i*nt+j)*c_simdBestPairAlignment+3] = 0;
+                }
             }
         }
 #endif