Add nblib virials function with test
authorJoe Jordan <ejjordan12@gmail.com>
Mon, 28 Jun 2021 09:18:17 +0000 (09:18 +0000)
committerJoe Jordan <ejjordan12@gmail.com>
Mon, 28 Jun 2021 09:18:17 +0000 (09:18 +0000)
api/nblib/CMakeLists.txt
api/nblib/tests/CMakeLists.txt
api/nblib/tests/virials.cpp [new file with mode: 0644]
api/nblib/virials.cpp [new file with mode: 0644]
api/nblib/virials.h [new file with mode: 0644]

index 82e406d806952503736ca23ac5b909ea9137f0db..0259b3b43f01fa2fe30a20830af7de207847b9b1 100644 (file)
@@ -99,6 +99,7 @@ target_sources(nblib
         simulationstate.cpp
         topologyhelpers.cpp
         topology.cpp
+        virials.cpp
         )
 
 gmx_target_compile_options(nblib)
index ae965706f9cbc17b6122c480a80e2fb738cc3056..6346ff29fb326058ea824e90696d6e4574b8e111 100644 (file)
@@ -61,6 +61,7 @@ gmx_add_gtest_executable(
         molecules.cpp
         nbnxmsetup.cpp
         topology.cpp
+        virials.cpp
     )
 target_link_libraries(${exename} PRIVATE mdrun_test_infrastructure)
 target_link_libraries(${exename} PRIVATE nblib_test_infrastructure nblib)
diff --git a/api/nblib/tests/virials.cpp b/api/nblib/tests/virials.cpp
new file mode 100644 (file)
index 0000000..235de30
--- /dev/null
@@ -0,0 +1,75 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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
+ * This implements virials tests
+ *
+ * \author Victor Holanda <victor.holanda@cscs.ch>
+ * \author Joe Jordan <ejjordan@kth.se>
+ * \author Prashanth Kanduri <kanduri@cscs.ch>
+ * \author Sebastian Keller <keller@cscs.ch>
+ */
+#include <gtest/gtest.h>
+
+#include "gromacs/mdtypes/forcerec.h"
+#include "nblib/box.h"
+#include "nblib/nbnxmsetuphelpers.h"
+#include "nblib/virials.h"
+
+namespace nblib
+{
+namespace test
+{
+
+TEST(VirialsTest, computeVirialTensorWorks)
+{
+    std::vector<Vec3> coords = { { 0, 1, 2 }, { 2, 3, 4 } };
+    std::vector<Vec3> forces = { { 2, 1, 2 }, { 4, 3, 4 } };
+    std::vector<Vec3> shiftForces(gmx::c_numShiftVectors, Vec3(0.0, 1.0, 0.0));
+    Box               box(1, 2, 3);
+    t_forcerec        forcerec;
+    updateForcerec(&forcerec, box.legacyMatrix());
+    std::vector<Vec3> shiftVectors(gmx::c_numShiftVectors);
+    // copy shift vectors from ForceRec
+    std::copy(forcerec.shift_vec.begin(), forcerec.shift_vec.end(), shiftVectors.begin());
+    std::vector<real> virialTest(9, 0);
+    computeVirialTensor(coords, forces, shiftVectors, shiftForces, box, virialTest);
+    std::vector<real> virialRef{ -4, -3, -4, -7, -5, -7, -10, -7, -10 };
+    EXPECT_EQ(virialRef, virialTest);
+}
+
+} // namespace test
+
+} // namespace nblib
diff --git a/api/nblib/virials.cpp b/api/nblib/virials.cpp
new file mode 100644 (file)
index 0000000..74f90a6
--- /dev/null
@@ -0,0 +1,85 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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.
+ */
+/*! \inpublicapi \file
+ * \brief
+ * Implements functionality to compute virials from force data
+ *
+ * \author Victor Holanda <victor.holanda@cscs.ch>
+ * \author Joe Jordan <ejjordan@kth.se>
+ * \author Prashanth Kanduri <kanduri@cscs.ch>
+ * \author Sebastian Keller <keller@cscs.ch>
+ */
+
+#include "gromacs/mdlib/calcvir.h"
+#include "gromacs/utility/arrayref.h"
+#include "nblib/box.h"
+#include "nblib/exception.h"
+#include "nblib/virials.h"
+
+namespace nblib
+{
+
+void computeVirialTensor(gmx::ArrayRef<const Vec3> coordinates,
+                         gmx::ArrayRef<const Vec3> forces,
+                         gmx::ArrayRef<const Vec3> shiftVectors,
+                         gmx::ArrayRef<const Vec3> shiftForces,
+                         const Box&                box,
+                         gmx::ArrayRef<real>       virialOutput)
+{
+    if (virialOutput.size() != 9)
+    {
+        throw InputException("Virial array size incorrect, should be 9");
+    }
+
+    // set virial output array to zero
+    std::fill(virialOutput.begin(), virialOutput.end(), 0.0);
+    // use legacy tensor format
+    rvec* virial = reinterpret_cast<rvec*>(virialOutput.data());
+
+    // compute the virials from surrounding boxes
+    const rvec* fshift    = as_rvec_array(shiftForces.data());
+    const rvec* shift_vec = as_rvec_array(shiftVectors.data());
+    calc_vir(numShiftVectors, shift_vec, fshift, virial, false, box.legacyMatrix());
+
+    // calculate partial virial, for local atoms only, based on short range
+    // NOTE: GROMACS uses a variable called mdatoms.homenr which is basically number of atoms in current processor
+    //       Used numAtoms from the coordinate array in its place
+    auto        numAtoms = coordinates.size();
+    const rvec* f        = as_rvec_array(forces.data());
+    const rvec* x        = as_rvec_array(coordinates.data());
+    calc_vir(numAtoms, x, f, virial, false, box.legacyMatrix());
+}
+
+} // namespace nblib
diff --git a/api/nblib/virials.h b/api/nblib/virials.h
new file mode 100644 (file)
index 0000000..bdaab73
--- /dev/null
@@ -0,0 +1,69 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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
+ * Implements functionality to compute virials from force data
+ *
+ * \author Victor Holanda <victor.holanda@cscs.ch>
+ * \author Joe Jordan <ejjordan@kth.se>
+ * \author Prashanth Kanduri <kanduri@cscs.ch>
+ * \author Sebastian Keller <keller@cscs.ch>
+ */
+#ifndef NBLIB_VIRIALS_H
+#define NBLIB_VIRIALS_H
+
+#include "nblib/basicdefinitions.h"
+#include "nblib/vector.h"
+
+namespace gmx
+{
+template<typename>
+class ArrayRef;
+}
+
+namespace nblib
+{
+class Box;
+
+//! Computes virial tensor and stores it in an array of size 9
+void computeVirialTensor(gmx::ArrayRef<const Vec3> coordinates,
+                         gmx::ArrayRef<const Vec3> forces,
+                         gmx::ArrayRef<const Vec3> shiftVectors,
+                         gmx::ArrayRef<const Vec3> shiftForces,
+                         const Box&                box,
+                         gmx::ArrayRef<real>       virialOutput);
+
+} // namespace nblib
+#endif // NBLIB_VIRIALS_H