Remove run-time check from calc_shifts
authorMark Abraham <mark.j.abraham@gmail.com>
Tue, 13 Aug 2019 08:27:48 +0000 (10:27 +0200)
committerErik Lindahl <erik.lindahl@gmail.com>
Thu, 15 Aug 2019 12:21:21 +0000 (14:21 +0200)
Now that this logic is unit tested, we don't have to
check it every step when using pressure coupling.

Also fixed some minor issues with pbcenums test that
I used as a template for the new source file.

Change-Id: I8f5dfaa44d2f2a92beee61fa6ba4ec94de21f50e

src/gromacs/pbcutil/pbc.cpp
src/gromacs/pbcutil/tests/CMakeLists.txt
src/gromacs/pbcutil/tests/pbc.cpp [new file with mode: 0644]
src/gromacs/pbcutil/tests/pbcenums.cpp
src/gromacs/pbcutil/tests/refdata/PbcTest_CalcShiftsWorks.xml [new file with mode: 0644]

index 0397683415186f86a1bc23231bff768680814190..55770da1f97532ecf5fa9ec87b0e645adcb68ee7 100644 (file)
@@ -1198,21 +1198,13 @@ void pbc_dx_d(const t_pbc *pbc, const dvec x1, const dvec x2, dvec dx)
 
 void calc_shifts(const matrix box, rvec shift_vec[])
 {
-    int k, l, m, d, n, test;
-
-    n = 0;
-    for (m = -D_BOX_Z; m <= D_BOX_Z; m++)
+    for (int n = 0, m = -D_BOX_Z; m <= D_BOX_Z; m++)
     {
-        for (l = -D_BOX_Y; l <= D_BOX_Y; l++)
+        for (int l = -D_BOX_Y; l <= D_BOX_Y; l++)
         {
-            for (k = -D_BOX_X; k <= D_BOX_X; k++, n++)
+            for (int k = -D_BOX_X; k <= D_BOX_X; k++, n++)
             {
-                test = XYZ2IS(k, l, m);
-                if (n != test)
-                {
-                    gmx_incons("inconsistent shift numbering");
-                }
-                for (d = 0; d < DIM; d++)
+                for (int d = 0; d < DIM; d++)
                 {
                     shift_vec[n][d] = k*box[XX][d] + l*box[YY][d] + m*box[ZZ][d];
                 }
index 05ccb45d2456ac59f36b2ce43565fb4c27890b93..2b9cfd5e66b298f0d10c71e9f060345a96660f2e 100644 (file)
@@ -32,6 +32,7 @@
 # To help us fund GROMACS development, we humbly ask that you cite
 # the research papers on the package. Check out http://www.gromacs.org.
 
-gmx_add_unit_test(PbcUtilityTests pbcutil-test
+gmx_add_unit_test(PbcutilUnitTest pbcutil-test
+                  pbc.cpp
                   pbcenums.cpp
                   )
diff --git a/src/gromacs/pbcutil/tests/pbc.cpp b/src/gromacs/pbcutil/tests/pbc.cpp
new file mode 100644 (file)
index 0000000..eba2aa7
--- /dev/null
@@ -0,0 +1,78 @@
+/*
+ * 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.
+ */
+/*! \internal \file
+ * \brief
+ * Tests PBC code
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ * \ingroup module_pbcutil
+ */
+#include "gmxpre.h"
+
+#include "gromacs/pbcutil/pbc.h"
+
+#include <gtest/gtest.h>
+
+#include "gromacs/math/vectypes.h"
+#include "gromacs/pbcutil/ishift.h"
+
+#include "testutils/refdata.h"
+
+namespace gmx
+{
+
+namespace test
+{
+
+TEST(PbcTest, CalcShiftsWorks)
+{
+    // Choose box vector entries whose magnitudes will lead to unique
+    // shift vector values when the largest box shift in any dimension
+    // is two.
+    const matrix box = {{ 0.01,     1, -100 },
+                        {  300, -0.03,    3 },
+                        {   -6,  -600, 0.06 }};
+    rvec         shiftVectors[SHIFTS];
+
+    calc_shifts(box, shiftVectors);
+
+    gmx::test::TestReferenceData    data;
+    gmx::test::TestReferenceChecker checker(data.rootChecker());
+    checker.checkSequence(std::begin(shiftVectors), std::end(shiftVectors), "ShiftVectors");
+}
+
+} // namespace test
+
+} // namespace gmx
index 864a197471d8950b26270ad8b5ef06edb6010d30..47458fd01f97b88161e1b11aaab010101ba8fc2e 100644 (file)
  * Tests PBC enumerations
  *
  * \author Paul Bauer <paul.bauer.q@gmail.com>
- * \ingroup module_math
+ * \ingroup module_pbcutil
  */
 #include "gmxpre.h"
 
 #include "gromacs/pbcutil/pbcenums.h"
 
-#include <gmock/gmock.h>
 #include <gtest/gtest.h>
 
 namespace gmx
diff --git a/src/gromacs/pbcutil/tests/refdata/PbcTest_CalcShiftsWorks.xml b/src/gromacs/pbcutil/tests/refdata/PbcTest_CalcShiftsWorks.xml
new file mode 100644 (file)
index 0000000..79091a1
--- /dev/null
@@ -0,0 +1,232 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Sequence Name="ShiftVectors">
+    <Int Name="Length">45</Int>
+    <Vector>
+      <Real Name="X">-294.01999999999998</Real>
+      <Real Name="Y">598.02999999999997</Real>
+      <Real Name="Z">196.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-294.00999999999999</Real>
+      <Real Name="Y">599.02999999999997</Real>
+      <Real Name="Z">96.939999999999998</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-294</Real>
+      <Real Name="Y">600.02999999999997</Real>
+      <Real Name="Z">-3.0600000000000001</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-293.99000000000001</Real>
+      <Real Name="Y">601.02999999999997</Real>
+      <Real Name="Z">-103.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-293.98000000000002</Real>
+      <Real Name="Y">602.02999999999997</Real>
+      <Real Name="Z">-203.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">5.9800000000000004</Real>
+      <Real Name="Y">598</Real>
+      <Real Name="Z">199.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">5.9900000000000002</Real>
+      <Real Name="Y">599</Real>
+      <Real Name="Z">99.939999999999998</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">6</Real>
+      <Real Name="Y">600</Real>
+      <Real Name="Z">-0.059999999999999998</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">6.0099999999999998</Real>
+      <Real Name="Y">601</Real>
+      <Real Name="Z">-100.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">6.0199999999999996</Real>
+      <Real Name="Y">602</Real>
+      <Real Name="Z">-200.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">305.98000000000002</Real>
+      <Real Name="Y">597.97000000000003</Real>
+      <Real Name="Z">202.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">305.99000000000001</Real>
+      <Real Name="Y">598.97000000000003</Real>
+      <Real Name="Z">102.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">306</Real>
+      <Real Name="Y">599.97000000000003</Real>
+      <Real Name="Z">2.9399999999999999</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">306.00999999999999</Real>
+      <Real Name="Y">600.97000000000003</Real>
+      <Real Name="Z">-97.060000000000002</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">306.01999999999998</Real>
+      <Real Name="Y">601.97000000000003</Real>
+      <Real Name="Z">-197.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-300.01999999999998</Real>
+      <Real Name="Y">-1.97</Real>
+      <Real Name="Z">197</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-300.00999999999999</Real>
+      <Real Name="Y">-0.96999999999999997</Real>
+      <Real Name="Z">97</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-300</Real>
+      <Real Name="Y">0.029999999999999999</Real>
+      <Real Name="Z">-3</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-299.99000000000001</Real>
+      <Real Name="Y">1.03</Real>
+      <Real Name="Z">-103</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-299.98000000000002</Real>
+      <Real Name="Y">2.0299999999999998</Real>
+      <Real Name="Z">-203</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-0.02</Real>
+      <Real Name="Y">-2</Real>
+      <Real Name="Z">200</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-0.01</Real>
+      <Real Name="Y">-1</Real>
+      <Real Name="Z">100</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">0</Real>
+      <Real Name="Y">0</Real>
+      <Real Name="Z">0</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">0.01</Real>
+      <Real Name="Y">1</Real>
+      <Real Name="Z">-100</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">0.02</Real>
+      <Real Name="Y">2</Real>
+      <Real Name="Z">-200</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">299.98000000000002</Real>
+      <Real Name="Y">-2.0299999999999998</Real>
+      <Real Name="Z">203</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">299.99000000000001</Real>
+      <Real Name="Y">-1.03</Real>
+      <Real Name="Z">103</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">300</Real>
+      <Real Name="Y">-0.029999999999999999</Real>
+      <Real Name="Z">3</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">300.00999999999999</Real>
+      <Real Name="Y">0.96999999999999997</Real>
+      <Real Name="Z">-97</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">300.01999999999998</Real>
+      <Real Name="Y">1.97</Real>
+      <Real Name="Z">-197</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-306.01999999999998</Real>
+      <Real Name="Y">-601.97000000000003</Real>
+      <Real Name="Z">197.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-306.00999999999999</Real>
+      <Real Name="Y">-600.97000000000003</Real>
+      <Real Name="Z">97.060000000000002</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-306</Real>
+      <Real Name="Y">-599.97000000000003</Real>
+      <Real Name="Z">-2.9399999999999999</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-305.99000000000001</Real>
+      <Real Name="Y">-598.97000000000003</Real>
+      <Real Name="Z">-102.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-305.98000000000002</Real>
+      <Real Name="Y">-597.97000000000003</Real>
+      <Real Name="Z">-202.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-6.0199999999999996</Real>
+      <Real Name="Y">-602</Real>
+      <Real Name="Z">200.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-6.0099999999999998</Real>
+      <Real Name="Y">-601</Real>
+      <Real Name="Z">100.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-6</Real>
+      <Real Name="Y">-600</Real>
+      <Real Name="Z">0.059999999999999998</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-5.9900000000000002</Real>
+      <Real Name="Y">-599</Real>
+      <Real Name="Z">-99.939999999999998</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">-5.9800000000000004</Real>
+      <Real Name="Y">-598</Real>
+      <Real Name="Z">-199.94</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">293.98000000000002</Real>
+      <Real Name="Y">-602.02999999999997</Real>
+      <Real Name="Z">203.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">293.99000000000001</Real>
+      <Real Name="Y">-601.02999999999997</Real>
+      <Real Name="Z">103.06</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">294</Real>
+      <Real Name="Y">-600.02999999999997</Real>
+      <Real Name="Z">3.0600000000000001</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">294.00999999999999</Real>
+      <Real Name="Y">-599.02999999999997</Real>
+      <Real Name="Z">-96.939999999999998</Real>
+    </Vector>
+    <Vector>
+      <Real Name="X">294.01999999999998</Real>
+      <Real Name="Y">-598.02999999999997</Real>
+      <Real Name="Z">-196.94</Real>
+    </Vector>
+  </Sequence>
+</ReferenceData>