Fix MPI deadlock in affinity setting
authorBerk Hess <hess@kth.se>
Mon, 24 Jun 2019 09:50:35 +0000 (11:50 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 1 Jul 2019 16:19:52 +0000 (18:19 +0200)
Commit 96d28d6b introduced the possiblity for deadlocks when
the affinity masks changed over time.

Fixes #2989

Change-Id: I454b38deb9ff11c90cdf4a19aaa80e79e8898df5

src/gromacs/mdrun/runner.cpp
src/gromacs/mdrunutility/threadaffinity.cpp
src/gromacs/mdrunutility/threadaffinity.h

index e6835bfc0c24e49c024ff3afec5e667c02775af4..2fb115f9698286c34baad0bc08dcdad22b1c7c73 100644 (file)
@@ -994,7 +994,7 @@ int Mdrunner::mdrunner()
     // that existing affinity setting was from OpenMP or something
     // else, so we run this code both before and after we initialize
     // the OpenMP support.
-    gmx_check_thread_affinity_set(mdlog, cr,
+    gmx_check_thread_affinity_set(mdlog,
                                   &hw_opt, hwinfo->nthreads_hw_avail, FALSE);
     /* Check and update the number of OpenMP threads requested */
     checkAndUpdateRequestedNumOpenmpThreads(&hw_opt, *hwinfo, cr, ms, physicalNodeComm.size_,
@@ -1190,7 +1190,7 @@ int Mdrunner::mdrunner()
          * - which indicates that probably the OpenMP library has changed it
          * since we first checked).
          */
-        gmx_check_thread_affinity_set(mdlog, cr,
+        gmx_check_thread_affinity_set(mdlog,
                                       &hw_opt, hwinfo->nthreads_hw_avail, TRUE);
 
         int numThreadsOnThisNode, intraNodeThreadOffset;
index 80ec82bd4a4367829f1956614230d93ba5249286..7d5bbebdd4399b8c5f9de6d3d88eb5f64c390340 100644 (file)
@@ -464,7 +464,6 @@ gmx_set_thread_affinity(const gmx::MDLogger         &mdlog,
  */
 void
 gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog,
-                              const t_commrec     *cr,
                               gmx_hw_opt_t        *hw_opt,
                               int  gmx_unused      nthreads_hw_avail,
                               gmx_bool             bAfterOpenmpInit)
@@ -494,16 +493,8 @@ gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog,
             }
         }
 
-        /* With thread-MPI this is needed as pinning might get turned off,
-         * which needs to be known before starting thread-MPI.
-         * With thread-MPI hw_opt is processed here on the master rank
-         * and passed to the other ranks later, so we only do this on master.
-         */
-        if (!SIMMASTER(cr))
-        {
-            return;
-        }
 #if !GMX_THREAD_MPI
+        // TODO: Remove this return so we get early reporting without thread-MPI as well
         return;
 #endif
     }
@@ -526,6 +517,8 @@ gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog,
         {
             fprintf(debug, "Failed to query affinity mask (error %d)", ret);
         }
+        // TODO: This might cause divergence between ranks on the same node.
+        //       Replace this return by a boolean and combine with bAllSet.
         return;
     }
 
@@ -540,18 +533,29 @@ gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog,
             fprintf(debug, "%d hardware threads detected, but %d was returned by CPU_COUNT",
                     nthreads_hw_avail, CPU_COUNT(&mask_current));
         }
+        // TODO: This might cause divergence between ranks on the same node.
+        //       Replace this return by a boolean and combine with bAllSet.
         return;
     }
 #endif /* CPU_COUNT */
 
+    /* Here we check whether all bits of the affinity mask are set.
+     * Note that this mask can change while the program is executing,
+     * e.g. by the system dynamically enabling or disabling cores or
+     * by some scheduler changing the affinities. Thus we can not assume
+     * that the result is the same on all ranks within a node
+     * when running this code at different times.
+     */
     gmx_bool bAllSet = TRUE;
     for (int i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
     {
         bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
     }
 
-    // With thread-MPI, bAllSet must already be the same on all ranks.
-#if GMX_LIB_MPI
+#if GMX_MPI
+    int mpiIsInitialized;
+    MPI_Initialized(&mpiIsInitialized);
+    if (mpiIsInitialized)
     {
         bool bAllSet_All;
         MPI_Allreduce(&bAllSet, &bAllSet_All, 1, MPI_C_BOOL, MPI_LAND, MPI_COMM_WORLD);
index 384842b79c9d65b6bb05b730fe5448f3590947e7..f5b890660f56ad501b3b8a82154031a8340804a8 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,2015,2016,2017,2018,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.
@@ -116,7 +116,7 @@ gmx_set_thread_affinity(const gmx::MDLogger         &mdlog,
  * variables for setting the affinity are set.
  */
 void
-gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog, const t_commrec *cr,
+gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog,
                               gmx_hw_opt_t *hw_opt, int ncpus,
                               gmx_bool bAfterOpenmpInit);