Add tests for calcvir
authorJoe Jordan <ejjordan12@gmail.com>
Mon, 29 Mar 2021 13:03:56 +0000 (13:03 +0000)
committerArtem Zhmurov <zhmurov@gmail.com>
Mon, 29 Mar 2021 13:03:56 +0000 (13:03 +0000)
src/gromacs/mdlib/tests/CMakeLists.txt
src/gromacs/mdlib/tests/calcvir.cpp [new file with mode: 0644]
src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAllAtomsInBox.xml [new file with mode: 0644]
src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAllAtomsInBoxScrew.xml [new file with mode: 0644]
src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewX.xml [new file with mode: 0644]
src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewXYZ.xml [new file with mode: 0644]
src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewY.xml [new file with mode: 0644]
src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewZ.xml [new file with mode: 0644]

index 96e4760622302a47ffec0d8ba8f5617da0387204..e50f85ab4f847cccf15a8f7b2ce03176dfdf5adb 100644 (file)
@@ -35,6 +35,7 @@
 gmx_add_unit_test(MdlibUnitTest mdlib-test HARDWARE_DETECTION
     CPP_SOURCE_FILES
         calc_verletbuf.cpp
+       calcvir.cpp
         constr.cpp
         constrtestdata.cpp
         constrtestrunners.cpp
diff --git a/src/gromacs/mdlib/tests/calcvir.cpp b/src/gromacs/mdlib/tests/calcvir.cpp
new file mode 100644 (file)
index 0000000..5504396
--- /dev/null
@@ -0,0 +1,153 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2020,2021, 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 for virial calculation.
+ *
+ * \author Joe Jordan <ejjordan@kth.se>
+ */
+#include "gmxpre.h"
+
+#include "gromacs/mdlib/calcvir.h"
+
+#include <gtest/gtest.h>
+
+#include "testutils/refdata.h"
+#include "testutils/testasserts.h"
+
+namespace gmx
+{
+namespace test
+{
+namespace
+{
+
+
+class CalcvirTest : public ::testing::Test
+{
+public:
+    TestReferenceData      refData_;
+    TestReferenceChecker   checker_;
+    std::vector<gmx::RVec> coordinates_;
+    std::vector<gmx::RVec> forces_;
+    int                    numVirialAtoms_;
+    tensor                 virial_{ { 0 } };
+    FloatingPointTolerance tolerances =
+            FloatingPointTolerance(1e-16, 1.0e-16, 1e-16, 1.0e-7, 1000, 100, false);
+
+    CalcvirTest() :
+        checker_(refData_.rootChecker()),
+        coordinates_{ { 1, 2, 3 },
+                      {
+                              4,
+                              5,
+                              6,
+                      },
+                      { 7, 8, 9 } },
+        forces_{ { 0.9, 0.1, 0.3 }, { 0.4, 0.7, 0.2 }, { 0.5, 1, 0.6 } },
+        numVirialAtoms_(coordinates_.size())
+    {
+        checker_.setDefaultTolerance(tolerances);
+    }
+
+private:
+};
+
+TEST_F(CalcvirTest, CanCalculateVirialAllAtomsInBox)
+{
+    calc_vir(numVirialAtoms_, as_rvec_array(coordinates_.data()), as_rvec_array(forces_.data()), virial_, false, nullptr);
+
+    checker_.checkVector(virial_[0], "Virial x");
+    checker_.checkVector(virial_[1], "Virial y");
+    checker_.checkVector(virial_[2], "Virial z");
+}
+
+TEST_F(CalcvirTest, CanCalculateVirialAllAtomsInBoxScrew)
+{
+
+    const matrix box = { { 10, 0, 0 }, { 0, 12, 0 }, { 0, 0, 14 } };
+    calc_vir(numVirialAtoms_, as_rvec_array(coordinates_.data()), as_rvec_array(forces_.data()), virial_, true, box);
+
+    checker_.checkVector(virial_[0], "Virial x");
+    checker_.checkVector(virial_[1], "Virial y");
+    checker_.checkVector(virial_[2], "Virial z");
+}
+
+TEST_F(CalcvirTest, CanCalculateVirialAtomsOutOfBoxScrewX)
+{
+
+    const matrix box = { { 5, 0, 0 }, { 0, 10, 0 }, { 0, 0, 12 } };
+    calc_vir(numVirialAtoms_, as_rvec_array(coordinates_.data()), as_rvec_array(forces_.data()), virial_, true, box);
+
+    checker_.checkVector(virial_[0], "Virial x");
+    checker_.checkVector(virial_[1], "Virial y");
+    checker_.checkVector(virial_[2], "Virial z");
+}
+
+TEST_F(CalcvirTest, CanCalculateVirialAtomsOutOfBoxScrewY)
+{
+
+    const matrix box = { { 10, 0, 0 }, { 0, 5, 0 }, { 0, 0, 12 } };
+    calc_vir(numVirialAtoms_, as_rvec_array(coordinates_.data()), as_rvec_array(forces_.data()), virial_, true, box);
+
+    checker_.checkVector(virial_[0], "Virial x");
+    checker_.checkVector(virial_[1], "Virial y");
+    checker_.checkVector(virial_[2], "Virial z");
+}
+
+TEST_F(CalcvirTest, CanCalculateVirialAtomsOutOfBoxScrewZ)
+{
+
+    const matrix box = { { 10, 0, 0 }, { 0, 12, 0 }, { 0, 0, 5 } };
+    calc_vir(numVirialAtoms_, as_rvec_array(coordinates_.data()), as_rvec_array(forces_.data()), virial_, true, box);
+
+    checker_.checkVector(virial_[0], "Virial x");
+    checker_.checkVector(virial_[1], "Virial y");
+    checker_.checkVector(virial_[2], "Virial z");
+}
+
+TEST_F(CalcvirTest, CanCalculateVirialAtomsOutOfBoxScrewXYZ)
+{
+
+    const matrix box = { { 4, 0, 0 }, { 0, 5, 0 }, { 0, 0, 6 } };
+    calc_vir(numVirialAtoms_, as_rvec_array(coordinates_.data()), as_rvec_array(forces_.data()), virial_, true, box);
+
+    checker_.checkVector(virial_[0], "Virial x");
+    checker_.checkVector(virial_[1], "Virial y");
+    checker_.checkVector(virial_[2], "Virial z");
+}
+
+} // namespace
+} // namespace test
+} // namespace gmx
diff --git a/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAllAtomsInBox.xml b/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAllAtomsInBox.xml
new file mode 100644 (file)
index 0000000..68e13d9
--- /dev/null
@@ -0,0 +1,19 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Vector Name="Virial x">
+    <Real Name="X">-3</Real>
+    <Real Name="Y">-4.9500000000000002</Real>
+    <Real Name="Z">-2.6499999999999999</Real>
+  </Vector>
+  <Vector Name="Virial y">
+    <Real Name="X">-3.9000000000000004</Real>
+    <Real Name="Y">-5.8499999999999996</Real>
+    <Real Name="Z">-3.2000000000000002</Real>
+  </Vector>
+  <Vector Name="Virial z">
+    <Real Name="X">-4.8000000000000007</Real>
+    <Real Name="Y">-6.75</Real>
+    <Real Name="Z">-3.75</Real>
+  </Vector>
+</ReferenceData>
diff --git a/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAllAtomsInBoxScrew.xml b/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAllAtomsInBoxScrew.xml
new file mode 100644 (file)
index 0000000..bf2f133
--- /dev/null
@@ -0,0 +1,19 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Vector Name="Virial x">
+    <Real Name="X">-5</Real>
+    <Real Name="Y">-8.4499999999999993</Real>
+    <Real Name="Z">-3.6499999999999999</Real>
+  </Vector>
+  <Vector Name="Virial y">
+    <Real Name="X">-6.3000000000000007</Real>
+    <Real Name="Y">-10.050000000000001</Real>
+    <Real Name="Z">-4.4000000000000004</Real>
+  </Vector>
+  <Vector Name="Virial z">
+    <Real Name="X">-7.6000000000000005</Real>
+    <Real Name="Y">-11.649999999999999</Real>
+    <Real Name="Z">-5.1500000000000004</Real>
+  </Vector>
+</ReferenceData>
diff --git a/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewX.xml b/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewX.xml
new file mode 100644 (file)
index 0000000..b95987f
--- /dev/null
@@ -0,0 +1,19 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Vector Name="Virial x">
+    <Real Name="X">-4</Real>
+    <Real Name="Y">-6.6999999999999993</Real>
+    <Real Name="Z">-3.1499999999999999</Real>
+  </Vector>
+  <Vector Name="Virial y">
+    <Real Name="X">-5.9000000000000004</Real>
+    <Real Name="Y">-9.3499999999999996</Real>
+    <Real Name="Z">-4.2000000000000002</Real>
+  </Vector>
+  <Vector Name="Virial z">
+    <Real Name="X">-7.2000000000000002</Real>
+    <Real Name="Y">-10.949999999999999</Real>
+    <Real Name="Z">-4.9500000000000002</Real>
+  </Vector>
+</ReferenceData>
diff --git a/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewXYZ.xml b/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewXYZ.xml
new file mode 100644 (file)
index 0000000..42b61b5
--- /dev/null
@@ -0,0 +1,19 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Vector Name="Virial x">
+    <Real Name="X">-3.7999999999999998</Real>
+    <Real Name="Y">-6.3499999999999996</Real>
+    <Real Name="Z">-3.0499999999999998</Real>
+  </Vector>
+  <Vector Name="Virial y">
+    <Real Name="X">-4.9000000000000004</Real>
+    <Real Name="Y">-7.5999999999999996</Real>
+    <Real Name="Z">-3.7000000000000002</Real>
+  </Vector>
+  <Vector Name="Virial z">
+    <Real Name="X">-6</Real>
+    <Real Name="Y">-8.8499999999999996</Real>
+    <Real Name="Z">-4.3499999999999996</Real>
+  </Vector>
+</ReferenceData>
diff --git a/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewY.xml b/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewY.xml
new file mode 100644 (file)
index 0000000..e7b0dc7
--- /dev/null
@@ -0,0 +1,19 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Vector Name="Virial x">
+    <Real Name="X">-5</Real>
+    <Real Name="Y">-8.4499999999999993</Real>
+    <Real Name="Z">-3.6499999999999999</Real>
+  </Vector>
+  <Vector Name="Virial y">
+    <Real Name="X">-4.9000000000000004</Real>
+    <Real Name="Y">-7.5999999999999996</Real>
+    <Real Name="Z">-3.7000000000000002</Real>
+  </Vector>
+  <Vector Name="Virial z">
+    <Real Name="X">-7.2000000000000002</Real>
+    <Real Name="Y">-10.949999999999999</Real>
+    <Real Name="Z">-4.9500000000000002</Real>
+  </Vector>
+</ReferenceData>
diff --git a/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewZ.xml b/src/gromacs/mdlib/tests/refdata/CalcvirTest_CanCalculateVirialAtomsOutOfBoxScrewZ.xml
new file mode 100644 (file)
index 0000000..b619e94
--- /dev/null
@@ -0,0 +1,19 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Vector Name="Virial x">
+    <Real Name="X">-5</Real>
+    <Real Name="Y">-8.4499999999999993</Real>
+    <Real Name="Z">-3.6499999999999999</Real>
+  </Vector>
+  <Vector Name="Virial y">
+    <Real Name="X">-6.3000000000000007</Real>
+    <Real Name="Y">-10.050000000000001</Real>
+    <Real Name="Z">-4.4000000000000004</Real>
+  </Vector>
+  <Vector Name="Virial z">
+    <Real Name="X">-5.8000000000000007</Real>
+    <Real Name="Y">-8.5</Real>
+    <Real Name="Z">-4.25</Real>
+  </Vector>
+</ReferenceData>