Fixed parallel distribution and thread count reporting
authorBerk Hess <hess@kth.se>
Tue, 24 Feb 2015 14:52:06 +0000 (15:52 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 24 Feb 2015 22:31:50 +0000 (23:31 +0100)
The domain decomposition no longer prints the atom count for all ranks
(which gets far too long at high parallelization), but rather av,
stddev, min and max.
Corrected the OpenMP thread count print, which often printed 0 thread
for (non-existent) PME ranks.
Fixes #1681, #1685.

Change-Id: Ia0e01c1c458e464ba2469d45c94d9e8c6b85d128

src/gromacs/domdec/domdec.cpp
src/gromacs/gmxlib/gmx_omp_nthreads.c

index c9bfb95cf985f28f58cc5a0f6acf58df747ebc78..43a2d6d8541c536fd6915db41d7bb72adc728145 100644 (file)
@@ -3911,7 +3911,7 @@ static void check_screw_box(matrix box)
     }
 }
 
-static void distribute_cg(FILE *fplog, gmx_int64_t step,
+static void distribute_cg(FILE *fplog,
                           matrix box, ivec tric_dir, t_block *cgs, rvec pos[],
                           gmx_domdec_t *dd)
 {
@@ -4065,18 +4065,33 @@ static void distribute_cg(FILE *fplog, gmx_int64_t step,
 
     if (fplog)
     {
-        char buf[22];
-        fprintf(fplog, "Charge group distribution at step %s:",
-                gmx_step_str(step, buf));
+        /* Here we avoid int overflows due to #atoms^2: use double, dsqr */
+        int    nat_sum, nat_min, nat_max;
+        double nat2_sum;
+
+        nat_sum  = 0;
+        nat2_sum = 0;
+        nat_min  = ma->nat[0];
+        nat_max  = ma->nat[0];
         for (i = 0; i < dd->nnodes; i++)
         {
-            fprintf(fplog, " %d", ma->ncg[i]);
+            nat_sum  += ma->nat[i];
+            nat2_sum += dsqr(ma->nat[i]);
+            nat_min   = std::min(nat_min, ma->nat[i]);
+            nat_max   = std::max(nat_max, ma->nat[i]);
         }
-        fprintf(fplog, "\n");
+        nat_sum  /= dd->nnodes;
+        nat2_sum /= dd->nnodes;
+
+        fprintf(fplog, "Atom distribution over %d domains: av %d stddev %d min %d max %d\n",
+                dd->nnodes,
+                nat_sum,
+                static_cast<int>(sqrt(nat2_sum - dsqr(nat_sum) + 0.5)),
+                nat_min, nat_max);
     }
 }
 
-static void get_cg_distribution(FILE *fplog, gmx_int64_t step, gmx_domdec_t *dd,
+static void get_cg_distribution(FILE *fplog, gmx_domdec_t *dd,
                                 t_block *cgs, matrix box, gmx_ddbox_t *ddbox,
                                 rvec pos[])
 {
@@ -4097,7 +4112,7 @@ static void get_cg_distribution(FILE *fplog, gmx_int64_t step, gmx_domdec_t *dd,
 
         set_dd_cell_sizes_slb(dd, ddbox, setcellsizeslbMASTER, npulse);
 
-        distribute_cg(fplog, step, box, ddbox->tric_dir, cgs, pos, dd);
+        distribute_cg(fplog, box, ddbox->tric_dir, cgs, pos, dd);
         for (i = 0; i < dd->nnodes; i++)
         {
             ma->ibuf[2*i]   = ma->ncg[i];
@@ -9466,7 +9481,7 @@ void dd_partition_system(FILE                *fplog,
         set_ddbox(dd, bMasterState, cr, ir, state_global->box,
                   TRUE, cgs_gl, state_global->x, &ddbox);
 
-        get_cg_distribution(fplog, step, dd, cgs_gl,
+        get_cg_distribution(fplog, dd, cgs_gl,
                             state_global->box, &ddbox, state_global->x);
 
         dd_distribute_state(dd, cgs_gl,
index ec6a1f6f354f3081a7305cb3368307a22b94bfbe..dc1688d36e829d11d2e7609e86cf780056c185b7 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014,2015, 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.
@@ -509,7 +509,7 @@ void gmx_omp_nthreads_init(FILE *fplog, t_commrec *cr,
     }
 #endif
 
-    reportOpenmpSettings(fplog, cr, bOMP, bSepPME, bFullOmpSupport);
+    reportOpenmpSettings(fplog, cr, bOMP, bFullOmpSupport, bSepPME);
     issueOversubscriptionWarning(fplog, cr, bSepPME, nthreads_hw_avail, nppn);
 }