Enable update groups
authorBerk Hess <hess@kth.se>
Wed, 19 Sep 2018 10:51:38 +0000 (12:51 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 8 Oct 2018 22:04:47 +0000 (00:04 +0200)
Also changed the naming of the decomposition unit. Everything was
named charge group, now it can also be atom or atom group.

Update groups are now mentioned in the manual and release notes

Change-Id: Ide81a759aa623be22672583dd5a8bac332c4c904

docs/reference-manual/algorithms/parallelization-domain-decomp.rst
docs/release-notes/performance.rst
src/gromacs/domdec/domdec.cpp
src/gromacs/domdec/domdec_topology.cpp
src/gromacs/domdec/partition.cpp
src/gromacs/domdec/redistribute.cpp
src/gromacs/mdlib/lincs.cpp

index 02d89c770b895e5e7e8f78c185297d1055c84135..76b862a5457b95c5d1e34fb29c601c060a13e830 100644 (file)
@@ -35,11 +35,18 @@ In the most general case of a triclinic unit cell, the space in divided
 with a 1-, 2-, or 3-D grid in parallelepipeds that we call domain
 decomposition cells. Each cell is assigned to a particle-particle rank.
 The system is partitioned over the ranks at the beginning of each MD
-step in which neighbor searching is performed. Since the neighbor
-searching is based on charge groups, charge groups are also the units
-for the domain decomposition. Charge groups are assigned to the cell
-where their center of geometry resides. Before the forces can be
-calculated, the coordinates from some neighboring cells need to be
+step in which neighbor searching is performed. The minimum unit of
+partitioning can be an atom, or a charge group with the (deprecated)
+group cut-off scheme or an update group. An update group is a group
+of atoms that has dependencies during update, which occurs when using
+constraints and/or virtual sites. Thus different update groups can be
+updated independenly. Currently update groups can only be used with at most
+two sequential constraints, which is the case when only constraining
+bonds involving hydrogen atoms. The advantages of update groups are that
+no communication is required in the update and that this allows updating part
+of the system while computing forces for other parts. Atom groups are assigned
+to the cell where their center of geometry resides. Before the forces can
+be calculated, the coordinates from some neighboring cells need to be
 communicated, and after the forces are calculated, the forces need to be
 communicated in the other direction. The communication and force
 assignment is based on zones that can cover one or multiple cells. An
index 3bbff7ba1a80c3ef5fa72cd352497a2a80b521c8..fb5e87151d859bc6155c0ce943a49075b1e384a5 100644 (file)
@@ -1,6 +1,15 @@
 Performance improvements
 ^^^^^^^^^^^^^^^^^^^^^^^^
 
+Implemented update groups
+"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
+Domain decomposition can now be based on so-called update groups. These
+are groups of atoms with dependencies during the update, which can be
+constraints and virtual sites. Update groups can typically be used when
+only bonds involving hydrogens are constrained and are enabled
+automatically when possible. This improves performance by eliminating
+MPI and OpenMP communication for constraints and virtual sites.
+
 PME on GPU when running free energy perturbations not involving charges
 """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
 PME can now be run on a GPU when doing free energy perturbations
index 5f556966e403ef7d2ed94b07c600a043e044d8d8..ec57e3ff013dc29a093693efa4e8b6287f7e49db 100644 (file)
 
 static const char *edlbs_names[int(DlbState::nr)] = { "off", "auto", "locked", "on", "on" };
 
-/* The size per charge group of the cggl_flag buffer in gmx_domdec_comm_t */
+/* The size per atom group of the cggl_flag buffer in gmx_domdec_comm_t */
 #define DD_CGIBS 2
 
 /* The flags for the cggl_flag buffer in gmx_domdec_comm_t */
@@ -2112,8 +2112,6 @@ static void setupUpdateGroups(const gmx::MDLogger &mdlog,
          * should be compatible with the box size.
          */
         comm->useUpdateGroups = (atomToAtomIntoDomainToDomainCutoff(*comm, 0) < cutoffMargin);
-        /* TODO: Enable update groups when all infrastructure is present */
-        comm->useUpdateGroups = false;
 
         if (comm->useUpdateGroups)
         {
@@ -2125,9 +2123,7 @@ static void setupUpdateGroups(const gmx::MDLogger &mdlog,
         }
         else
         {
-            /* TODO: Remove this comment when enabling update groups
-             * GMX_LOG(mdlog.info).appendTextFormatted("The combination of rlist and box size prohibits the use of update groups\n");
-             */
+            GMX_LOG(mdlog.info).appendTextFormatted("The combination of rlist and box size prohibits the use of update groups\n");
             comm->updateGroupingPerMoleculetype.clear();
             comm->updateGroupsCog.reset(nullptr);
         }
@@ -2661,7 +2657,21 @@ static void writeSettings(gmx::TextWriter       *log,
         bInterCGVsites ||
         dd->splitConstraints || dd->splitSettles)
     {
-        log->writeLine("The maximum allowed distance for charge groups involved in interactions is:");
+        std::string decompUnits;
+        if (comm->bCGs)
+        {
+            decompUnits = "charge groups";
+        }
+        else if (comm->useUpdateGroups)
+        {
+            decompUnits = "atom groups";
+        }
+        else
+        {
+            decompUnits = "atoms";
+        }
+
+        log->writeLineFormatted("The maximum allowed distance for %s involved in interactions is:", decompUnits.c_str());
         log->writeLineFormatted("%40s  %-7s %6.3f nm", "non-bonded interactions", "", comm->cutoff);
 
         if (bDynLoadBal)
index e553e87f3f3ce94d3404b8e79b2a6154e39bbe82..512860bc5451990909ce9f4460ccc5c7f1e89082 100644 (file)
@@ -2694,7 +2694,7 @@ void dd_bonded_cg_distance(const gmx::MDLogger &mdlog,
 
     if (*r_2b > 0 || *r_mb > 0)
     {
-        GMX_LOG(mdlog.info).appendText("Initial maximum inter charge-group distances:");
+        GMX_LOG(mdlog.info).appendText("Initial maximum distances in bonded interactions:");
         if (*r_2b > 0)
         {
             GMX_LOG(mdlog.info).appendTextFormatted(
index 0fe65fcfb83d87ac9bc6258d0abdda2f2135fe48..7d9e8b6cd9b7aacdd6d863c6d75b1d815641e17f 100644 (file)
@@ -2802,7 +2802,8 @@ static void dd_sort_state(gmx_domdec_t *dd, rvec *cgcm, t_forcerec *fr, t_state
     dd->ncg_home = sort->sorted.size();
     if (debug)
     {
-        fprintf(debug, "Set the new home charge group count to %d\n",
+        fprintf(debug, "Set the new home %s count to %d\n",
+                dd->comm->bCGs ? "charge group" : "atom",
                 dd->ncg_home);
     }
 
@@ -3303,7 +3304,7 @@ void dd_partition_system(FILE                    *fplog,
 
     ncg_home_old = dd->ncg_home;
 
-    /* When repartitioning we mark charge groups that will move to neighboring
+    /* When repartitioning we mark atom groups that will move to neighboring
      * DD cells, but we do not move them right away for performance reasons.
      * Thus we need to keep track of how many charge groups will move for
      * obtaining correct local charge group / atom counts.
index 3dcba3b73e779278ab19445b02f5b679b64a99cd..2e7920b9cd9e27d9d070ebf36b3e4c096e70a7b7 100644 (file)
@@ -166,25 +166,32 @@ static void print_cg_move(FILE *fplog,
                           gmx_bool bHaveCgcmOld, real limitd,
                           rvec cm_old, rvec cm_new, real pos_d)
 {
-    gmx_domdec_comm_t *comm;
-    char               buf[22];
+    const gmx_domdec_comm_t *comm = dd->comm;
+    std::string              mesg;
 
-    comm = dd->comm;
+    fprintf(fplog, "\nStep %" PRId64 ":\n", step);
 
-    fprintf(fplog, "\nStep %s:\n", gmx_step_str(step, buf));
-    if (limitd > 0)
+    if (comm->bCGs)
+    {
+        mesg += "The charge group starting at atom";
+    }
+    else if (comm->useUpdateGroups)
     {
-        fprintf(fplog, "%s %d moved more than the distance allowed by the domain decomposition (%f) in direction %c\n",
-                dd->comm->bCGs ? "The charge group starting at atom" : "Atom",
-                ddglatnr(dd, dd->atomGrouping().block(cg).begin()), limitd, dim2char(dim));
+        mesg += "The update group starting at atom";
     }
     else
     {
-        /* We don't have a limiting distance available: don't print it */
-        fprintf(fplog, "%s %d moved more than the distance allowed by the domain decomposition in direction %c\n",
-                dd->comm->bCGs ? "The charge group starting at atom" : "Atom",
-                ddglatnr(dd, dd->atomGrouping().block(cg).begin()), dim2char(dim));
+        mesg += "Atom";
     }
+    mesg += gmx::formatString(" %d moved more than the distance allowed by the domain decomposition",
+                              ddglatnr(dd, dd->atomGrouping().block(cg).begin()));
+    if (limitd > 0)
+    {
+        mesg += gmx::formatString(" (%f)", limitd);
+    }
+    mesg += gmx::formatString(" in direction %c\n", dim2char(dim));
+    fprintf(fplog, "%s", mesg.c_str());
+
     fprintf(fplog, "distance out of cell %f\n",
             dir == 1 ? pos_d - comm->cell_x1[dim] : pos_d - comm->cell_x0[dim]);
     if (bHaveCgcmOld)
@@ -217,9 +224,8 @@ static void cg_move_error(FILE *fplog,
     print_cg_move(stderr, dd, step, cg, dim, dir,
                   bHaveCgcmOld, limitd, cm_old, cm_new, pos_d);
     gmx_fatal(FARGS,
-              "%s moved too far between two domain decomposition steps\n"
-              "This usually means that your system is not well equilibrated",
-              dd->comm->bCGs ? "An atom group" : "An atom");
+              "One or more atoms moved too far between two domain decomposition steps.\n"
+              "This usually means that your system is not well equilibrated");
 }
 
 static void rotate_state_atom(t_state *state, int a)
index f57d901a034b42bf21d1c197848da938e6e014aa..cde84e0002976e07e8d5297b85d2e8496480ff24 100644 (file)
@@ -1599,7 +1599,7 @@ Lincs *init_lincs(FILE *fplog, const gmx_mtop_t &mtop,
         fprintf(fplog, "The number of constraints is %d\n", li->ncg);
         if (bPLINCS)
         {
-            fprintf(fplog, "There are inter charge-group constraints,\n"
+            fprintf(fplog, "There are constraints between atoms in different decomposition domains,\n"
                     "will communicate selected coordinates each lincs iteration\n");
         }
         if (li->ncg_triangle > 0)