Made gromacs work again without MPI
authorBerk Hess <hess@kth.se>
Wed, 13 Oct 2021 19:43:53 +0000 (19:43 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 13 Oct 2021 19:43:53 +0000 (19:43 +0000)
Made the communication functions in network.h that do nothing with a
single rank also work without MPI and thread-MPI. One conditional
in mdrunner() needed to be fixed.

src/gromacs/gmxlib/network.cpp
src/gromacs/gmxlib/network.h
src/gromacs/mdrun/runner.cpp

index c5edd4eb51eea884d5dc6f9ff190f4a0a5abae60..5bac7f1e0f973059cea997da96fd6d7232f429f9 100644 (file)
@@ -235,23 +235,21 @@ void gmx_barrier(MPI_Comm gmx_unused communicator)
 
 void gmx_bcast(int gmx_unused nbytes, void gmx_unused* b, MPI_Comm gmx_unused communicator)
 {
-#if !GMX_MPI
-    GMX_RELEASE_ASSERT(false, "Invalid call to gmx_bcast");
-#else
+    // Without MPI we have a single rank, so bcast is a no-op
+#if GMX_MPI
     MPI_Bcast(b, nbytes, MPI_BYTE, 0, communicator);
 #endif
 }
 
 void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unused* cr)
 {
+    // Without MPI we have a single rank, so sum is a no-op
+#if GMX_MPI
     if (cr->sizeOfMyGroupCommunicator == 1)
     {
         return;
     }
 
-#if !GMX_MPI
-    GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumd");
-#else
     if (cr->nc.bUse)
     {
         if (cr->nc.rank_intra == 0)
@@ -278,14 +276,13 @@ void gmx_sumd(int gmx_unused nr, double gmx_unused r[], const t_commrec gmx_unus
 
 void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unused* cr)
 {
+    // Without MPI we have a single rank, so sum is a no-op
+#if GMX_MPI
     if (cr->sizeOfMyGroupCommunicator == 1)
     {
         return;
     }
 
-#if !GMX_MPI
-    GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumf");
-#else
     if (cr->nc.bUse)
     {
         /* Use two step summing.  */
@@ -312,14 +309,13 @@ void gmx_sumf(int gmx_unused nr, float gmx_unused r[], const t_commrec gmx_unuse
 
 void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused* cr)
 {
+    // Without MPI we have a single rank, so sum is a no-op
+#if GMX_MPI
     if (cr->sizeOfMyGroupCommunicator == 1)
     {
         return;
     }
 
-#if !GMX_MPI
-    GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumi");
-#else
     if (cr->nc.bUse)
     {
         /* Use two step summing */
@@ -344,35 +340,6 @@ void gmx_sumi(int gmx_unused nr, int gmx_unused r[], const t_commrec gmx_unused*
 #endif
 }
 
-void gmx_sumli(int gmx_unused nr, int64_t gmx_unused r[], const t_commrec gmx_unused* cr)
-{
-#if !GMX_MPI
-    GMX_RELEASE_ASSERT(false, "Invalid call to gmx_sumli");
-#else
-    if (cr->nc.bUse)
-    {
-        /* Use two step summing */
-        if (cr->nc.rank_intra == 0)
-        {
-            MPI_Reduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
-            /* Sum with the buffers reversed */
-            MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->nc.comm_inter);
-        }
-        else
-        {
-            /* This is here because of the silly MPI specification
-                that MPI_IN_PLACE should be put in sendbuf instead of recvbuf */
-            MPI_Reduce(r, nullptr, nr, MPI_INT64_T, MPI_SUM, 0, cr->nc.comm_intra);
-        }
-        MPI_Bcast(r, nr, MPI_INT64_T, 0, cr->nc.comm_intra);
-    }
-    else
-    {
-        MPI_Allreduce(MPI_IN_PLACE, r, nr, MPI_INT64_T, MPI_SUM, cr->mpi_comm_mygroup);
-    }
-#endif
-}
-
 const char* opt2fn_master(const char* opt, int nfile, const t_filenm fnm[], t_commrec* cr)
 {
     return SIMMASTER(cr) ? opt2fn(opt, nfile, fnm) : nullptr;
index 3bc41df9b3d7722db17929ab55904c2429585da3..3ac723afd7b7c9834c2eb78f1428bea3b308952d 100644 (file)
@@ -4,7 +4,7 @@
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
  * Copyright (c) 2013,2014,2015,2016,2017 by the GROMACS development team.
- * 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.
@@ -78,20 +78,29 @@ void gmx_setup_nodecomm(FILE* fplog, struct t_commrec* cr);
 //! Wait until all processes in communicator have reached the barrier
 void gmx_barrier(MPI_Comm communicator);
 
-//! Broadcast nbytes bytes from the master to communicator
+/*! Broadcast nbytes bytes from the master to communicator
+ *
+ * Can be called with a single rank or without MPI
+ */
 void gmx_bcast(int nbytes, void* b, MPI_Comm communicator);
 
+/*! Calculate the global sum of an array of ints
+ *
+ * Can be called with a single rank or without MPI
+ */
 void gmx_sumi(int nr, int r[], const struct t_commrec* cr);
-/* Calculate the global sum of an array of ints */
-
-void gmx_sumli(int nr, int64_t r[], const struct t_commrec* cr);
-/* Calculate the global sum of an array of large ints */
 
+/*! Calculate the global sum of an array of floats
+ *
+ * Can be called with a single rank or without MPI
+ */
 void gmx_sumf(int nr, float r[], const struct t_commrec* cr);
-/* Calculate the global sum of an array of floats */
 
+/*! Calculate the global sum of an array of doubles
+ *
+ * Can be called with a single rank or without MPI
+ */
 void gmx_sumd(int nr, double r[], const struct t_commrec* cr);
-/* Calculate the global sum of an array of doubles */
 
 #if GMX_DOUBLE
 #    define gmx_sum gmx_sumd
index 8c3fbdf517855a22e7f5e9689a515619c193e7c5..3bf1d5071cd6d9d21bb56dd62ad182f178488c99 100644 (file)
@@ -937,7 +937,7 @@ int Mdrunner::mdrunner()
         // master and spawned threads joins at the end of this block.
     }
 
-    GMX_RELEASE_ASSERT(ms || simulationCommunicator != MPI_COMM_NULL,
+    GMX_RELEASE_ASSERT(!GMX_MPI || ms || simulationCommunicator != MPI_COMM_NULL,
                        "Must have valid communicator unless running a multi-simulation");
     CommrecHandle crHandle = init_commrec(simulationCommunicator);
     t_commrec*    cr       = crHandle.get();