Added missing DD cycle counting
authorBerk Hess <hess@kth.se>
Tue, 16 Jun 2015 18:16:56 +0000 (20:16 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 16 Jun 2015 20:28:31 +0000 (22:28 +0200)
Cycle counting was missing for DD repartitioning after replica
exchange or coord swap. Removed DD cycle counting for initial DD.
Moved DD cycle counting into dd_partition_system and added
subract_cycles function with assertion to detect cycle wrapping.

Fixes #1677.

Change-Id: I7f1b19397b36456f1d120dbc0080146a384def5a

src/gromacs/domdec/domdec.cpp
src/gromacs/mdlib/minimize.cpp
src/gromacs/timing/wallcycle.c
src/programs/mdrun/md.cpp

index 28ffa7ab2dc0cf7b35ddf6fcdf93d3155debe774..0b6001e6676d202d83eb81e759ade83500a44c32 100644 (file)
@@ -9329,6 +9329,8 @@ void dd_partition_system(FILE                *fplog,
     real               grid_density;
     char               sbuf[22];
 
+    wallcycle_start(wcycle, ewcDOMDEC);
+
     dd   = cr->dd;
     comm = dd->comm;
 
@@ -9915,4 +9917,6 @@ void dd_partition_system(FILE                *fplog,
         check_index_consistency(dd, top_global->natoms, ncg_mtop(top_global),
                                 "after partitioning");
     }
+
+    wallcycle_stop(wcycle, ewcDOMDEC);
 }
index b29a8dc7e2e715f2793341fc9fb0007452c1d3f6..e70ceb8be27528e1ec7eb924bde71630c57d51fb 100644 (file)
@@ -673,14 +673,12 @@ static void em_dd_partition_system(FILE *fplog, int step, t_commrec *cr,
                                    t_nrnb *nrnb, gmx_wallcycle_t wcycle)
 {
     /* Repartition the domain decomposition */
-    wallcycle_start(wcycle, ewcDOMDEC);
     dd_partition_system(fplog, step, cr, FALSE, 1,
                         NULL, top_global, ir,
                         &ems->s, &ems->f,
                         mdatoms, top, fr, vsite, NULL, constr,
                         nrnb, wcycle, FALSE);
     dd_store_state(cr->dd, &ems->s);
-    wallcycle_stop(wcycle, ewcDOMDEC);
 }
 
 static void evaluate_energy(FILE *fplog, t_commrec *cr,
index b0b6a70852ac3b6062152b532256c682027b66c7..4961ce58a3303529327e44df7fd005a9c644f14a 100644 (file)
@@ -388,6 +388,17 @@ static gmx_bool is_pme_subcounter(int ewc)
     return (ewc >= ewcPME_REDISTXF && ewc < ewcPMEWAITCOMM);
 }
 
+/* Subtract counter ewc_sub timed inside a timing block for ewc_main */
+static void subtract_cycles(wallcc_t *wcc, int ewc_main, int ewc_sub)
+{
+    if (wcc[ewc_sub].n > 0)
+    {
+        assert(wcc[ewc_main].c >= wcc[ewc_sub].c);
+
+        wcc[ewc_main].c -= wcc[ewc_sub].c;
+    }
+}
+
 void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
 {
     wallcc_t *wcc;
@@ -439,26 +450,15 @@ void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
         }
     }
 
-    if (wcc[ewcDDCOMMLOAD].n > 0)
-    {
-        wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMLOAD].c;
-    }
-    if (wcc[ewcDDCOMMBOUND].n > 0)
-    {
-        wcc[ewcDOMDEC].c -= wcc[ewcDDCOMMBOUND].c;
-    }
-    if (wcc[ewcPME_FFTCOMM].n > 0)
-    {
-        wcc[ewcPME_FFT].c -= wcc[ewcPME_FFTCOMM].c;
-    }
+    subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMLOAD);
+    subtract_cycles(wcc, ewcDOMDEC, ewcDDCOMMBOUND);
+
+    subtract_cycles(wcc, ewcPME_FFT, ewcPME_FFTCOMM);
 
     if (cr->npmenodes == 0)
     {
         /* All nodes do PME (or no PME at all) */
-        if (wcc[ewcPMEMESH].n > 0)
-        {
-            wcc[ewcFORCE].c -= wcc[ewcPMEMESH].c;
-        }
+        subtract_cycles(wcc, ewcFORCE, ewcPMEMESH);
     }
     else
     {
@@ -466,6 +466,7 @@ void wallcycle_sum(t_commrec *cr, gmx_wallcycle_t wc)
         if (wcc[ewcPMEMESH].n > 0)
         {
             /* This must be a PME only node, calculate the Wait + Comm. time */
+            assert(wcc[ewcRUN].c >= wcc[ewcPMEMESH].c);
             wcc[ewcPMEWAITCOMM].c = wcc[ewcRUN].c - wcc[ewcPMEMESH].c;
         }
     }
index 4c34f399c159880c1400d98148b085f3c8f5923e..c377038774f94a683c65df23c39995bd31dd395f 100644 (file)
@@ -414,7 +414,7 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
                             state_global, top_global, ir,
                             state, &f, mdatoms, top, fr,
                             vsite, shellfc, constr,
-                            nrnb, wcycle, FALSE);
+                            nrnb, NULL, FALSE);
 
     }
 
@@ -917,7 +917,6 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
             if (DOMAINDECOMP(cr))
             {
                 /* Repartition the domain decomposition */
-                wallcycle_start(wcycle, ewcDOMDEC);
                 dd_partition_system(fplog, step, cr,
                                     bMasterState, nstglobalcomm,
                                     state_global, top_global, ir,
@@ -925,7 +924,6 @@ double do_md(FILE *fplog, t_commrec *cr, int nfile, const t_filenm fnm[],
                                     vsite, shellfc, constr,
                                     nrnb, wcycle,
                                     do_verbose && !bPMETunePrinting);
-                wallcycle_stop(wcycle, ewcDOMDEC);
             }
         }