Move constraint and bonded filtering info into DDSystemInfo
authorBerk Hess <hess@kth.se>
Wed, 4 Sep 2019 06:45:19 +0000 (08:45 +0200)
committerBerk Hess <hess@kth.se>
Wed, 4 Sep 2019 18:56:37 +0000 (20:56 +0200)
Change-Id: Ic3cc4b5da309ab3f8e719739a5a162cd2c2d62bf

src/gromacs/domdec/domdec.cpp
src/gromacs/domdec/domdec.h
src/gromacs/domdec/domdec_constraints.cpp
src/gromacs/domdec/domdec_internal.h
src/gromacs/domdec/domdec_struct.h
src/gromacs/domdec/domdec_topology.cpp
src/gromacs/domdec/partition.cpp
src/gromacs/mdlib/constr.cpp

index 71dbfb60b5b824f7044d1bce35378d8e106161b1..37514e9ef2397e4479006c86547811a9d753a4b6 100644 (file)
@@ -667,7 +667,7 @@ real dd_cutoff_multibody(const gmx_domdec_t *dd)
             {
                 r = std::min(r, comm.cellsize_min[dd->dim[di]]);
             }
-            if (comm.bBondComm)
+            if (comm.systemInfo.filterBondedCommunication)
             {
                 r = std::max(r, comm.cutoff_mbody);
             }
@@ -1081,6 +1081,11 @@ int dd_pme_maxshift_y(const gmx_domdec_t *dd)
     }
 }
 
+bool ddHaveSplitConstraints(const gmx_domdec_t &dd)
+{
+    return dd.comm->systemInfo.haveSplitConstraints;
+}
+
 void dd_cycles_add(const gmx_domdec_t *dd, float cycles, int ddCycl)
 {
     /* Note that the cycles value can be incorrect, either 0 or some
@@ -2167,13 +2172,13 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
 
     if (systemInfo.useUpdateGroups)
     {
-        dd->splitConstraints = false;
-        dd->splitSettles     = false;
+        systemInfo.haveSplitConstraints = false;
+        systemInfo.haveSplitSettles     = false;
     }
     else
     {
-        dd->splitConstraints = gmx::inter_charge_group_constraints(*mtop);
-        dd->splitSettles     = gmx::inter_charge_group_settles(*mtop);
+        systemInfo.haveSplitConstraints = gmx::inter_charge_group_constraints(*mtop);
+        systemInfo.haveSplitSettles     = gmx::inter_charge_group_settles(*mtop);
     }
 
     if (ir->rlist == 0)
@@ -2191,8 +2196,8 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
     systemInfo.minCutoffForMultiBody = 0;
 
     /* Determine the minimum cell size limit, affected by many factors */
-    systemInfo.cellsizeLimit = 0;
-    comm->bBondComm          = FALSE;
+    systemInfo.cellsizeLimit             = 0;
+    systemInfo.filterBondedCommunication = false;
 
     /* We do not allow home atoms to move beyond the neighboring domain
      * between domain decomposition steps, which limits the cell size.
@@ -2230,7 +2235,7 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
                 atomToAtomIntoDomainToDomainCutoff(systemInfo, options.minimumCommunicationRange);
             if (options.useBondedCommunication)
             {
-                comm->bBondComm = (systemInfo.minCutoffForMultiBody > systemInfo.cutoff);
+                systemInfo.filterBondedCommunication = (systemInfo.minCutoffForMultiBody > systemInfo.cutoff);
             }
             else
             {
@@ -2268,7 +2273,8 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
                 {
                     r_bonded        = std::max(r_2b, r_mb);
                     r_bonded_limit  = tenPercentMargin*r_bonded;
-                    comm->bBondComm = TRUE;
+                    /* This is the (only) place where we turn on the filtering */
+                    systemInfo.filterBondedCommunication = true;
                 }
                 else
                 {
@@ -2296,7 +2302,7 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
     }
 
     real rconstr = 0;
-    if (dd->splitConstraints && options.constraintCommunicationRange <= 0)
+    if (systemInfo.haveSplitConstraints && options.constraintCommunicationRange <= 0)
     {
         /* There is a cell size limit due to the constraints (P-LINCS) */
         rconstr = gmx::constr_r_max(mdlog, mtop, ir);
@@ -2357,7 +2363,7 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
         if (ddSetup.numDomains[XX] == 0)
         {
             char     buf[STRLEN];
-            gmx_bool bC = (dd->splitConstraints && rconstr > r_bonded_limit);
+            gmx_bool bC = (systemInfo.haveSplitConstraints && rconstr > r_bonded_limit);
             sprintf(buf, "Change the number of ranks or mdrun option %s%s%s",
                     !bC ? "-rdd" : "-rcon",
                     comm->dlbState != DlbState::offUser ? " or -dds" : "",
@@ -2477,7 +2483,7 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
     comm->cellsize_limit = systemInfo.cellsizeLimit;
     if (systemInfo.haveInterDomainBondeds && comm->cutoff_mbody == 0)
     {
-        if (comm->bBondComm || !isDlbDisabled(comm))
+        if (systemInfo.filterBondedCommunication || !isDlbDisabled(comm))
         {
             /* Set the bonded communication distance to halfway
              * the minimum and the maximum,
@@ -2491,7 +2497,7 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
                 comm->cutoff_mbody = std::min(comm->cutoff_mbody,
                                               options.dlbScaling*acs);
             }
-            if (!comm->bBondComm)
+            if (!systemInfo.filterBondedCommunication)
             {
                 /* Without bBondComm do not go beyond the n.b. cut-off */
                 comm->cutoff_mbody = std::min(comm->cutoff_mbody, systemInfo.cutoff);
@@ -2520,7 +2526,8 @@ static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
     {
         fprintf(debug, "Bonded atom communication beyond the cut-off: %s\n"
                 "cellsize limit %f\n",
-                gmx::boolToString(comm->bBondComm), comm->cellsize_limit);
+                gmx::boolToString(systemInfo.filterBondedCommunication),
+                comm->cellsize_limit);
     }
 
     if (MASTER(cr))
@@ -2557,7 +2564,7 @@ void dd_init_bondeds(FILE *fplog,
 
     comm = dd->comm;
 
-    if (comm->bBondComm)
+    if (comm->systemInfo.filterBondedCommunication)
     {
         /* Communicate atoms beyond the cut-off for bonded interactions */
         comm = dd->comm;
@@ -2646,7 +2653,8 @@ static void writeSettings(gmx::TextWriter       *log,
 
     if (comm->systemInfo.haveInterDomainBondeds ||
         haveInterDomainVsites ||
-        dd->splitConstraints || dd->splitSettles)
+        comm->systemInfo.haveSplitConstraints ||
+        comm->systemInfo.haveSplitSettles)
     {
         std::string decompUnits;
         if (comm->systemInfo.useUpdateGroups)
@@ -2685,14 +2693,14 @@ static void writeSettings(gmx::TextWriter       *log,
                                     std::max(comm->systemInfo.cutoff, comm->cutoff_mbody));
             log->writeLineFormatted("%40s  %-7s %6.3f nm",
                                     "multi-body bonded interactions", "(-rdd)",
-                                    (comm->bBondComm || isDlbOn(dd->comm)) ? comm->cutoff_mbody : std::min(comm->systemInfo.cutoff, limit));
+                                    (comm->systemInfo.filterBondedCommunication || isDlbOn(dd->comm)) ? comm->cutoff_mbody : std::min(comm->systemInfo.cutoff, limit));
         }
         if (haveInterDomainVsites)
         {
             log->writeLineFormatted("%40s  %-7s %6.3f nm",
                                     "virtual site constructions", "(-rcon)", limit);
         }
-        if (dd->splitConstraints || dd->splitSettles)
+        if (comm->systemInfo.haveSplitConstraints || comm->systemInfo.haveSplitSettles)
         {
             std::string separation = gmx::formatString("atoms separated by up to %d constraints",
                                                        1+ir->nProjOrder);
index 6f8fa8c8ef42ff0a174175b26966a353f2401783..948d6ec55660560d801570207d59ec7c99542851 100644 (file)
@@ -153,6 +153,9 @@ int dd_pme_maxshift_x(const gmx_domdec_t *dd);
 /*! \brief Returns the maximum shift for coordinate communication in PME, dim y */
 int dd_pme_maxshift_y(const gmx_domdec_t *dd);
 
+/*! \brief Return whether constraints, not including settles, cross domain boundaries */
+bool ddHaveSplitConstraints(const gmx_domdec_t &dd);
+
 /*! \brief Initialized the domain decomposition, chooses the DD grid and PME ranks, return the DD struct */
 gmx_domdec_t *
 init_domain_decomposition(const gmx::MDLogger            &mdlog,
index 8954364bb040f9e711a4834a2e9b1889ee0af59d..5e3d382ca113b5218fc7253aa13c7a6ce4ae212c 100644 (file)
@@ -69,6 +69,7 @@
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/smalloc.h"
 
+#include "domdec_internal.h"
 #include "domdec_specatomcomm.h"
 
 /*! \brief Struct used during constraint setup with domain decomposition */
@@ -421,7 +422,9 @@ int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
     // This code should not be called unless this condition is true,
     // because that's the only time init_domdec_constraints is
     // called...
-    GMX_RELEASE_ASSERT(dd->splitConstraints || dd->splitSettles, "dd_make_local_constraints called when there are no local constraints");
+    GMX_RELEASE_ASSERT(dd->comm->systemInfo.haveSplitConstraints ||
+                       dd->comm->systemInfo.haveSplitSettles,
+                       "dd_make_local_constraints called when there are no local constraints");
     // ... and init_domdec_constraints always sets
     // dd->constraint_comm...
     GMX_RELEASE_ASSERT(dd->constraint_comm, "Invalid use of dd_make_local_constraints before construction of constraint_comm");
@@ -454,7 +457,7 @@ int dd_make_local_constraints(gmx_domdec_t *dd, int at_start,
 
     gmx::ArrayRef < const std::vector < int>> at2settle_mt;
     /* When settle works inside charge groups, we assigned them already */
-    if (dd->splitSettles)
+    if (dd->comm->systemInfo.haveSplitSettles)
     {
         // TODO Perhaps gmx_domdec_constraints_t should keep a valid constr?
         GMX_RELEASE_ASSERT(constr != nullptr, "Must have valid constraints object");
index 812ab3aa567fa28228891780e438745c0a75bf1c..bb6c37fc6387b7a59029f7f30b4629affd010250 100644 (file)
@@ -456,6 +456,14 @@ struct DDSystemInfo
     real cutoff = 0;
     //! The lower limit for the DD cell size
     real cellsizeLimit = 0;
+
+    //! Can atoms connected by constraints be assigned to different domains?
+    bool haveSplitConstraints = false;
+    //! Can atoms connected by settles be assigned to different domains?
+    bool haveSplitSettles = false;
+
+    //! Whether to only communicate atoms beyond the non-bonded cut-off when they are involved in bonded interactions with non-local atoms
+    bool filterBondedCommunication = false;
 };
 
 /*! \brief Struct for domain decomposition communication
@@ -509,9 +517,7 @@ struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
     //! Centers of mass of local update groups
     std::unique_ptr<gmx::UpdateGroupsCog> updateGroupsCog;
 
-    /* Data for the optional bonded interaction atom communication range */
-    /**< Only communicate atoms beyond the non-bonded cut-off when they are involved in bonded interactions with non-local atoms */
-    gmx_bool  bBondComm = false;
+    /* Data for the optional filtering of communication of atoms for bonded interactions */
     /**< Links between cg's through bonded interactions */
     t_blocka *cglink = nullptr;
     /**< Local cg availability, TODO: remove when group scheme is removed */
index 0d4cd47e95c0ca1016ace649ae425cecfd4ea60b..4a43576e33067ff78ada75ff25d62ebd2d1f7907 100644 (file)
@@ -189,11 +189,6 @@ struct gmx_domdec_t { //NOLINT(clang-analyzer-optin.performance.Padding)
     /* Only available on the master node */
     std::unique_ptr<AtomDistribution> ma;
 
-    /* Can atoms connected by constraints be assigned to different domains? */
-    bool splitConstraints = false;
-    /* Can atoms connected by settles be assigned to different domains? */
-    bool splitSettles = false;
-
     /* Global atom number to interaction list */
     gmx_reverse_top_t  *reverse_top    = nullptr;
     int                 nbonded_global = 0;
index 5009a18fafb23a2fcffd6eec7678ea5723d29e23..5acc05c167ce8cfd88c777605b19f6c30edf633e 100644 (file)
@@ -746,7 +746,8 @@ void dd_make_reverse_top(FILE *fplog,
     dd->reverse_top  = new gmx_reverse_top_t;
     *dd->reverse_top =
         make_reverse_top(mtop, ir->efep != efepNO,
-                         !dd->splitConstraints, !dd->splitSettles,
+                         !dd->comm->systemInfo.haveSplitConstraints,
+                         !dd->comm->systemInfo.haveSplitSettles,
                          bBCheck, &dd->nbonded_global);
 
     gmx_reverse_top_t *rt = dd->reverse_top;
@@ -796,7 +797,7 @@ void dd_make_reverse_top(FILE *fplog,
         init_domdec_vsites(dd, vsite->numInterUpdategroupVsites);
     }
 
-    if (dd->splitConstraints || dd->splitSettles)
+    if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
     {
         init_domdec_constraints(dd, mtop);
     }
index 15a3a8809c18e6c796bf471658b8315d25cf5846..cf1adc3acf312965725af4ed3ad096a8ee135098 100644 (file)
@@ -1913,7 +1913,7 @@ static void setup_dd_communication(gmx_domdec_t *dd,
         comm->dth.resize(numThreads);
     }
 
-    bBondComm = comm->bBondComm;
+    bBondComm = comm->systemInfo.filterBondedCommunication;
 
     /* Do we need to determine extra distances for multi-body bondeds? */
     bDistMB = (comm->systemInfo.haveInterDomainMultiBodyBondeds && isDlbOn(dd->comm) && dd->ndim > 1);
@@ -3212,7 +3212,7 @@ void dd_partition_system(FILE                    *fplog,
                 }
                 break;
             case DDAtomRanges::Type::Constraints:
-                if (dd->splitConstraints || dd->splitSettles)
+                if (dd->comm->systemInfo.haveSplitConstraints || dd->comm->systemInfo.haveSplitSettles)
                 {
                     /* Only for inter-cg constraints we need special code */
                     n = dd_make_local_constraints(dd, n, &top_global, fr->cginfo.data(),
index 11c7bb7d6d477aa0f53a5bd8252da534aae72199..2be5fd3c3a553e9dc360114002286ba2b5660e15 100644 (file)
@@ -1069,13 +1069,13 @@ Constraints::Impl::Impl(const gmx_mtop_t     &mtop_p,
         {
             lincsd = init_lincs(log, mtop,
                                 nflexcon, at2con_mt,
-                                DOMAINDECOMP(cr) && cr->dd->splitConstraints,
+                                DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd),
                                 ir.nLincsIter, ir.nProjOrder);
         }
 
         if (ir.eConstrAlg == econtSHAKE)
         {
-            if (DOMAINDECOMP(cr) && cr->dd->splitConstraints)
+            if (DOMAINDECOMP(cr) && ddHaveSplitConstraints(*cr->dd))
             {
                 gmx_fatal(FARGS, "SHAKE is not supported with domain decomposition and constraint that cross charge group boundaries, use LINCS");
             }