Eliminate GMX_MIMIC symbol from source code
authorMark Abraham <mark.j.abraham@gmail.com>
Thu, 16 Sep 2021 10:45:17 +0000 (10:45 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 16 Sep 2021 10:45:17 +0000 (10:45 +0000)
src/config.h.cmakein
src/gromacs/mdrun/mimic.cpp
src/gromacs/mimic/CMakeLists.txt
src/gromacs/mimic/communicator.cpp
src/gromacs/mimic/communicator.h
src/gromacs/mimic/communicator_stub.cpp [new file with mode: 0644]

index 298d9f40b6103e41fe7c7bbbb2b6d31ce668ad8b..b082d5e702a24b591b13e63e40c8d447db9b2b9b 100644 (file)
 /* Build using clang analyzer */
 #cmakedefine01 GMX_CLANG_ANALYZER
 
-/* Use MiMiC QM/MM interface */
-#cmakedefine01 GMX_MIMIC
-
 /* Use Interactive Molecular Dynamics */
 #cmakedefine01 GMX_IMD
 
index 29b74e0d8bdd4fc5c0e37326a4fa232794ec8375..995336e1574d2bea863c19a813e39e54c1681345 100644 (file)
@@ -452,7 +452,7 @@ void gmx::LegacySimulator::do_mimic()
 
         if (MASTER(cr))
         {
-            MimicCommunicator::getCoords(&state_global->x, state_global->natoms);
+            MimicCommunicator::getCoords(state_global->x, state_global->natoms);
         }
 
         if (ir->efep != FreeEnergyPerturbationType::No)
index 84040bd4c7173dfcd59cac6ea45d8a41984a0572..d4832ab1c5dfb7c864613f716a3703b12e27ea70 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+# Copyright (c) 2018,2019,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.
 
 # Set up the module library
 add_library(mimic INTERFACE)
-file(GLOB MIMIC_SOURCES *.cpp)
+if (GMX_MIMIC)
+    file(GLOB MIMIC_SOURCES
+        communicator.cpp
+        utilities.cpp
+        )
+else()
+    file(GLOB MIMIC_SOURCES
+        communicator_stub.cpp
+        utilities.cpp
+        )
+endif()
 set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${MIMIC_SOURCES} PARENT_SCOPE)
 
 # Source files have the following dependencies on library infrastructure.
index 177dc8813b0f883a2ecd2d98cc24e349da3eaee9..9b3f8a9fc1b2b6df8f91ac5abf66bc6bbccae862 100644 (file)
 #include <unordered_set>
 
 #include "gromacs/math/units.h"
-#include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/arrayref.h"
 
-#if GMX_MIMIC
-#    include <DataTypes.h>
-#    include <MessageApi.h>
-#endif
+// Include headers from MiMiC library
+#include <DataTypes.h>
+#include <MessageApi.h>
 
-// When not built in a configuration with QMMM support, much of this
-// code is unreachable by design. Tell clang not to warn about it.
-#pragma GCC diagnostic push
-#pragma GCC diagnostic ignored "-Wmissing-noreturn"
-
-#if !GMX_MIMIC
-//! \brief Definitions to stub the ones defined in DataTypes.h
-constexpr int TYPE_INT = 0, TYPE_DOUBLE = 0;
-
-/*! \brief Stub communication library function to call in case if
- * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
- */
-static void MCL_init_client(const char*) // NOLINT(readability-named-parameter)
+namespace gmx
 {
-    GMX_RELEASE_ASSERT(
-            GMX_MIMIC,
-            "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
-}
 
-/*! \brief Stub communication library function to call in case if
- * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
- */
-static void MCL_send(void*, int, int, int) // NOLINT(readability-named-parameter)
-{
-    GMX_RELEASE_ASSERT(
-            GMX_MIMIC,
-            "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
-}
-
-/*! \brief Stub communication library function to call in case if
- * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
- */
-static void MCL_receive(void*, int, int, int) // NOLINT(readability-named-parameter)
-{
-    GMX_RELEASE_ASSERT(
-            GMX_MIMIC,
-            "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
-}
-
-/*! \brief Stub communication library function to call in case if
- * GROMACS is compiled without MiMiC. Calling causes GROMACS to exit!
- */
-static void MCL_destroy()
-{
-    GMX_RELEASE_ASSERT(
-            GMX_MIMIC,
-            "GROMACS is compiled without MiMiC support! Please, recompile with -DGMX_MIMIC=ON");
-}
-#endif
-
-void gmx::MimicCommunicator::init()
+void MimicCommunicator::init()
 {
     char path[GMX_PATH_MAX];
     gmx_getcwd(path, GMX_PATH_MAX);
     MCL_init_client(path);
 }
 
-void gmx::MimicCommunicator::sendInitData(gmx_mtop_t* mtop, PaddedHostVector<gmx::RVec> coords)
+void MimicCommunicator::sendInitData(gmx_mtop_t* mtop, ArrayRef<const RVec> coords)
 {
     MCL_send(&mtop->natoms, 1, TYPE_INT, 0);
     MCL_send(&mtop->atomtypes.nr, 1, TYPE_INT, 0);
@@ -141,7 +93,7 @@ void gmx::MimicCommunicator::sendInitData(gmx_mtop_t* mtop, PaddedHostVector<gmx
                 bonds.push_back(offset + at1 + 1);
                 bonds.push_back(offset + at2 + 1);
                 bondLengths.push_back(static_cast<double>(mtop->ffparams.iparams[contype].constr.dA)
-                                      / gmx::c_bohr2Nm);
+                                      / c_bohr2Nm);
             }
 
             for (int ncon = 0; ncon < nsettle; ++ncon)
@@ -243,32 +195,32 @@ void gmx::MimicCommunicator::sendInitData(gmx_mtop_t* mtop, PaddedHostVector<gmx
     MCL_send(&*convertedCoords.begin(), 3 * mtop->natoms, TYPE_DOUBLE, 0);
 }
 
-int64_t gmx::MimicCommunicator::getStepNumber()
+int64_t MimicCommunicator::getStepNumber()
 {
     int steps;
     MCL_receive(&steps, 1, TYPE_INT, 0);
     return steps;
 }
 
-void gmx::MimicCommunicator::getCoords(PaddedHostVector<RVec>* x, const int natoms)
+void MimicCommunicator::getCoords(ArrayRef<RVec> x, const int natoms)
 {
     std::vector<double> coords(natoms * 3);
     MCL_receive(&*coords.begin(), 3 * natoms, TYPE_DOUBLE, 0);
     for (int j = 0; j < natoms; ++j)
     {
-        (*x)[j][0] = static_cast<real>(coords[j * 3] * gmx::c_bohr2Nm);
-        (*x)[j][1] = static_cast<real>(coords[j * 3 + 1] * gmx::c_bohr2Nm);
-        (*x)[j][2] = static_cast<real>(coords[j * 3 + 2] * gmx::c_bohr2Nm);
+        x[j][0] = static_cast<real>(coords[j * 3] * gmx::c_bohr2Nm);
+        x[j][1] = static_cast<real>(coords[j * 3 + 1] * gmx::c_bohr2Nm);
+        x[j][2] = static_cast<real>(coords[j * 3 + 2] * gmx::c_bohr2Nm);
     }
 }
 
-void gmx::MimicCommunicator::sendEnergies(real energy)
+void MimicCommunicator::sendEnergies(real energy)
 {
     double convertedEnergy = energy / (gmx::c_hartree2Kj * gmx::c_avogadro);
     MCL_send(&convertedEnergy, 1, TYPE_DOUBLE, 0);
 }
 
-void gmx::MimicCommunicator::sendForces(gmx::ArrayRef<gmx::RVec> forces, int natoms)
+void MimicCommunicator::sendForces(ArrayRef<RVec> forces, int natoms)
 {
     std::vector<double> convertedForce;
     for (int j = 0; j < natoms; ++j)
@@ -280,9 +232,9 @@ void gmx::MimicCommunicator::sendForces(gmx::ArrayRef<gmx::RVec> forces, int nat
     MCL_send(&*convertedForce.begin(), convertedForce.size(), TYPE_DOUBLE, 0);
 }
 
-void gmx::MimicCommunicator::finalize()
+void MimicCommunicator::finalize()
 {
     MCL_destroy();
 }
 
-#pragma GCC diagnostic pop
+} // namespace gmx
index 6d2b9c276139420de2eeb25d2e539e841d7821d8..ea7f18b93fb4990d0714e9011ca62408ea817e6c 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,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.
 #ifndef GMX_MIMIC_COMMUNICATOR_H
 #define GMX_MIMIC_COMMUNICATOR_H
 
-#include "gromacs/gpu_utils/hostallocator.h"
-#include "gromacs/math/paddedvector.h"
 #include "gromacs/mdlib/constr.h"
 #include "gromacs/topology/topology.h"
 #include "gromacs/utility/futil.h"
 
 namespace gmx
 {
+
+template<class T>
+class ArrayRef;
+
 /**
  * \inlibraryapi
  * \internal \brief
@@ -73,7 +75,7 @@ public:
      * @param mtop global topology data
      * @param coords coordinates of all atoms
      */
-    static void sendInitData(gmx_mtop_t* mtop, PaddedHostVector<gmx::RVec> coords);
+    static void sendInitData(gmx_mtop_t* mtop, ArrayRef<const RVec> coords);
 
     /*! \brief
      * Gets the number of MD steps to perform from MiMiC
@@ -88,7 +90,7 @@ public:
      * @param x array of coordinates to fill
      * @param natoms number of atoms in the system
      */
-    static void getCoords(PaddedHostVector<RVec>* x, int natoms);
+    static void getCoords(ArrayRef<RVec> x, int natoms);
 
     /*! \brief
      * Send the potential energy value to MiMiC
diff --git a/src/gromacs/mimic/communicator_stub.cpp b/src/gromacs/mimic/communicator_stub.cpp
new file mode 100644 (file)
index 0000000..c417942
--- /dev/null
@@ -0,0 +1,97 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018,2019,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.
+ */
+#include "gmxpre.h"
+
+#include "communicator.h"
+
+#include "gromacs/math/units.h"
+#include "gromacs/utility/arrayref.h"
+#include "gromacs/utility/exceptions.h"
+
+namespace gmx
+{
+
+#ifdef __clang__
+#    pragma clang diagnostic push
+#    pragma clang diagnostic ignored "-Wmissing-noreturn"
+#endif
+
+void MimicCommunicator::init()
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+void MimicCommunicator::sendInitData(gmx_mtop_t* /*mtop*/, ArrayRef<const RVec> /*coords*/)
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+int64_t MimicCommunicator::getStepNumber()
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+void MimicCommunicator::getCoords(ArrayRef<RVec> /*x*/, const int /*natoms*/)
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+void MimicCommunicator::sendEnergies(real /*energy*/)
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+void MimicCommunicator::sendForces(ArrayRef<RVec> /*forces*/, int /*natoms*/)
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+void MimicCommunicator::finalize()
+{
+    GMX_THROW(InternalError(
+            "GROMACS is compiled without MiMiC support! Please, reconfigure with -DGMX_MIMIC=ON"));
+}
+
+#ifdef __clang__
+#    pragma clang diagnostic pop
+#endif
+
+} // namespace gmx