Collect settings for DD in DDSettings
authorBerk Hess <hess@kth.se>
Wed, 4 Sep 2019 11:49:10 +0000 (13:49 +0200)
committerSzilárd Páll <pall.szilard@gmail.com>
Wed, 4 Sep 2019 22:09:26 +0000 (00:09 +0200)
This change is only refactoring.

Change-Id: Icb32f8f5211ff3033c4d5c70afbfd24d4a8ffbf3

src/gromacs/domdec/cellsizes.cpp
src/gromacs/domdec/dlbtiming.cpp
src/gromacs/domdec/domdec.cpp
src/gromacs/domdec/domdec_internal.h
src/gromacs/domdec/domdec_network.cpp
src/gromacs/domdec/domdec_struct.h
src/gromacs/domdec/partition.cpp
src/gromacs/domdec/utility.h

index 22c9e6f4bdc2d7b9e1273573bf63ea73058722f4..05b01b16d22568dcb0e4d934029c37b9685b17ed 100644 (file)
@@ -564,7 +564,7 @@ static void set_dd_cell_sizes_dlb_root(gmx_domdec_t *dd,
     int                range[] = { 0, 0 };
 
     /* Convert the maximum change from the input percentage to a fraction */
-    const real          change_limit = comm->dlb_scale_lim*0.01;
+    const real          change_limit = comm->ddSettings.dlb_scale_lim*0.01;
 
     const int           ncd          = dd->nc[dim];
 
index 8466bfc69c192289f96082d6a220d252bfe0112d..6e920ea3b4ba2c31b7b758cb4189265a2d2b7f4d 100644 (file)
@@ -86,7 +86,7 @@ static BalanceRegion *getBalanceRegion(const gmx_domdec_t *dd)
 void DDBalanceRegionHandler::openRegionCpuImpl(DdAllowBalanceRegionReopen gmx_unused allowReopen) const
 {
     BalanceRegion *reg = getBalanceRegion(dd_);
-    if (dd_->comm->bRecordLoad)
+    if (dd_->comm->ddSettings.recordLoad)
     {
         GMX_ASSERT(allowReopen == DdAllowBalanceRegionReopen::yes || !reg->isOpen, "Should not open an already opened region");
 
@@ -211,7 +211,7 @@ static double force_flop_count(const t_nrnb *nrnb)
 
 void dd_force_flop_start(gmx_domdec_t *dd, t_nrnb *nrnb)
 {
-    if (dd->comm->eFlop)
+    if (dd->comm->ddSettings.eFlop)
     {
         dd->comm->flop -= force_flop_count(nrnb);
     }
@@ -219,7 +219,7 @@ void dd_force_flop_start(gmx_domdec_t *dd, t_nrnb *nrnb)
 
 void dd_force_flop_stop(gmx_domdec_t *dd, t_nrnb *nrnb)
 {
-    if (dd->comm->eFlop)
+    if (dd->comm->ddSettings.eFlop)
     {
         dd->comm->flop += force_flop_count(nrnb);
         dd->comm->flop_n++;
index 65076a4423ade9456916f9176b2ac92390ab455d..05d365bd9d3fb79458e18df9f253686aa39cebe4 100644 (file)
@@ -1377,7 +1377,7 @@ static void setup_neighbor_relations(gmx_domdec_t *dd)
         dd->comm->cellsizesWithDlb.resize(dd->ndim);
     }
 
-    if (dd->comm->bRecordLoad)
+    if (dd->comm->ddSettings.recordLoad)
     {
         make_load_communicators(dd);
     }
@@ -2298,19 +2298,20 @@ getSystemInfo(const gmx::MDLogger           &mdlog,
 static void set_dd_limits_and_grid(const gmx::MDLogger &mdlog,
                                    t_commrec *cr, gmx_domdec_t *dd,
                                    const DomdecOptions &options,
-                                   const gmx::MdrunOptions &mdrunOptions,
+                                   const DDSettings &ddSettings,
                                    const gmx_mtop_t *mtop,
                                    const t_inputrec *ir,
                                    const matrix box,
                                    gmx::ArrayRef<const gmx::RVec> xGlobal,
                                    gmx_ddbox_t *ddbox)
 {
-    gmx_domdec_comm_t *comm             = dd->comm;
+    gmx_domdec_comm_t *comm = dd->comm;
+    comm->ddSettings        = ddSettings;
 
     /* Initialize to GPU share count to 0, might change later */
     comm->nrank_gpu_shared = 0;
 
-    comm->dlbState         = determineInitialDlbState(mdlog, options.dlbOption, comm->bRecordLoad, mdrunOptions, ir);
+    comm->dlbState         = comm->ddSettings.initialDlbState;
     dd_dlb_set_should_check_whether_to_turn_dlb_on(dd, TRUE);
     /* To consider turning DLB on after 2*nstlist steps we need to check
      * at partitioning count 3. Thus we need to increase the first count by 2.
@@ -2926,38 +2927,42 @@ static void set_ddgrid_parameters(const gmx::MDLogger &mdlog,
                                  static_cast<int>(vol_frac*natoms_tot));
 }
 
-/*! \brief Set some important DD parameters that can be modified by env.vars */
-static void set_dd_envvar_options(const gmx::MDLogger &mdlog,
-                                  gmx_domdec_t *dd, int rank_mysim)
+/*! \brief Get some important DD parameters which can be modified by env.vars */
+static DDSettings
+getDDSettings(const gmx::MDLogger     &mdlog,
+              const DomdecOptions     &options,
+              const gmx::MdrunOptions &mdrunOptions,
+              const t_inputrec        &ir)
 {
-    gmx_domdec_comm_t *comm = dd->comm;
+    DDSettings ddSettings;
 
-    dd->bSendRecv2      = (dd_getenv(mdlog, "GMX_DD_USE_SENDRECV2", 0) != 0);
-    comm->dlb_scale_lim = dd_getenv(mdlog, "GMX_DLB_MAX_BOX_SCALING", 10);
-    comm->eFlop         = dd_getenv(mdlog, "GMX_DLB_BASED_ON_FLOPS", 0);
-    int recload         = dd_getenv(mdlog, "GMX_DD_RECORD_LOAD", 1);
-    comm->nstDDDump     = dd_getenv(mdlog, "GMX_DD_NST_DUMP", 0);
-    comm->nstDDDumpGrid = dd_getenv(mdlog, "GMX_DD_NST_DUMP_GRID", 0);
-    comm->DD_debug      = dd_getenv(mdlog, "GMX_DD_DEBUG", 0);
+    ddSettings.useSendRecv2  = (dd_getenv(mdlog, "GMX_DD_USE_SENDRECV2", 0) != 0);
+    ddSettings.dlb_scale_lim = dd_getenv(mdlog, "GMX_DLB_MAX_BOX_SCALING", 10);
+    ddSettings.eFlop         = dd_getenv(mdlog, "GMX_DLB_BASED_ON_FLOPS", 0);
+    const int recload        = dd_getenv(mdlog, "GMX_DD_RECORD_LOAD", 1);
+    ddSettings.nstDDDump     = dd_getenv(mdlog, "GMX_DD_NST_DUMP", 0);
+    ddSettings.nstDDDumpGrid = dd_getenv(mdlog, "GMX_DD_NST_DUMP_GRID", 0);
+    ddSettings.DD_debug      = dd_getenv(mdlog, "GMX_DD_DEBUG", 0);
 
-    if (dd->bSendRecv2)
+    if (ddSettings.useSendRecv2)
     {
         GMX_LOG(mdlog.info).appendText("Will use two sequential MPI_Sendrecv calls instead of two simultaneous non-blocking MPI_Irecv and MPI_Isend pairs for constraint and vsite communication");
     }
 
-    if (comm->eFlop)
+    if (ddSettings.eFlop)
     {
         GMX_LOG(mdlog.info).appendText("Will load balance based on FLOP count");
-        if (comm->eFlop > 1)
-        {
-            srand(1 + rank_mysim);
-        }
-        comm->bRecordLoad = TRUE;
+        ddSettings.recordLoad = true;
     }
     else
     {
-        comm->bRecordLoad = (wallcycle_have_counter() && recload > 0);
+        ddSettings.recordLoad = (wallcycle_have_counter() && recload > 0);
     }
+
+    ddSettings.initialDlbState =
+        determineInitialDlbState(mdlog, options.dlbOption, ddSettings.recordLoad, mdrunOptions, &ir);
+
+    return ddSettings;
 }
 
 gmx_domdec_t::gmx_domdec_t(const t_inputrec &ir) :
@@ -2984,10 +2989,16 @@ gmx_domdec_t *init_domain_decomposition(const gmx::MDLogger           &mdlog,
 
     dd->comm = init_dd_comm();
 
-    set_dd_envvar_options(mdlog, dd, cr->nodeid);
+    DDSettings  ddSettings = getDDSettings(mdlog, options, mdrunOptions, *ir);
+    if (ddSettings.eFlop > 1)
+    {
+        /* Ensure that we have different random flop counts on different ranks */
+        srand(1 + cr->nodeid);
+    }
 
     gmx_ddbox_t ddbox = {0};
-    set_dd_limits_and_grid(mdlog, cr, dd, options, mdrunOptions,
+    set_dd_limits_and_grid(mdlog, cr, dd, options,
+                           ddSettings,
                            mtop, ir,
                            box, xGlobal,
                            &ddbox);
index 6464035791f24ca80d86669d8b806eebf82c0ae0..539b852b2600157d4f5e6ad12b7c7fce8bb7071b 100644 (file)
@@ -470,6 +470,40 @@ struct DDSystemInfo
     bool increaseMultiBodyCutoff = false;
 };
 
+/*! \brief Settings that affect the behavior of the domain decomposition
+ *
+ * These settings depend on options chosen by the user, set by enviroment
+ * variables, as well as hardware support. The initial DLB state also
+ * depends on the integrator.
+ *
+ * Note: Settings that depend on the simulated system are in DDSystemInfo.
+ */
+struct DDSettings
+{
+    //! Use MPI_Sendrecv communication instead of non-blocking calls
+    bool useSendRecv2 = false;
+
+    /* Information for managing the dynamic load balancing */
+    //! Maximum DLB scaling per load balancing step in percent
+    int  dlb_scale_lim = 0;
+    //! Flop counter (0=no,1=yes,2=with (eFlop-1)*5% noise
+    int  eFlop = 0;
+
+    //! Whether we should record the load
+    bool recordLoad = false;
+
+    /* Debugging */
+    //! Step interval for dumping the local+non-local atoms to pdb
+    int  nstDDDump = 0;
+    //! Step interval for duming the DD grid to pdb
+    int  nstDDDumpGrid = 0;
+    //! DD debug print level: 0, 1, 2
+    int  DD_debug = 0;
+
+    //! The DLB state at the start of the run
+    DlbState initialDlbState = DlbState::offCanTurnOn;
+};
+
 /*! \brief Struct for domain decomposition communication
  *
  * This struct contains most information about domain decomposition
@@ -482,6 +516,9 @@ struct DDSystemInfo
  */
 struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
 {
+    /**< Constant parameters that control DD behavior */
+    DDSettings ddSettings;
+
     /* PME and Cartesian communicator stuff */
     /**< The number of decomposition dimensions for PME, 0: no PME */
     int         npmedecompdim = 0;
@@ -631,8 +668,6 @@ struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
     std::vector<DDCellsizesWithDlb> cellsizesWithDlb;
 
     /* Stuff for load communication */
-    /**< Should we record the load */
-    gmx_bool        bRecordLoad = false;
     /**< The recorded load data */
     domdec_load_t  *load = nullptr;
     /**< The number of MPI ranks sharing the GPU our rank is using */
@@ -644,10 +679,6 @@ struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
     MPI_Comm        mpi_comm_gpu_shared;
 #endif
 
-    /* Information for managing the dynamic load balancing */
-    /**< Maximum DLB scaling per load balancing step in percent */
-    int            dlb_scale_lim = 0;
-
     /**< Struct for timing the force load balancing region */
     BalanceRegion *balanceRegion = nullptr;
 
@@ -658,8 +689,6 @@ struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
     int    cycl_n[ddCyclNr] = { };
     /**< The maximum cycle count */
     float  cycl_max[ddCyclNr] = { };
-    /** Flop counter (0=no,1=yes,2=with (eFlop-1)*5% noise */
-    int    eFlop = 0;
     /**< Total flops counted */
     double flop = 0.0;
     /**< The number of flop recordings */
@@ -703,14 +732,6 @@ struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
 
     /** The last partition step */
     int64_t partition_step = INT_MIN;
-
-    /* Debugging */
-    /**< Step interval for dumping the local+non-local atoms to pdb */
-    int  nstDDDump = 0;
-    /**< Step interval for duming the DD grid to pdb */
-    int  nstDDDumpGrid = 0;
-    /**< DD debug print level: 0, 1, 2 */
-    int  DD_debug = 0;
 };
 
 /*! \brief DD zone permutation
index 31bacdc185a6ce60f9006a144f008b3b9ed38a52..28f71eaf8b58589d62b260d8de343895cf0611ea 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2008,2009,2010,2012,2013,2014,2015,2016,2018, by the GROMACS development team, led by
+ * Copyright (c) 2008,2009,2010,2012,2013,2014,2015,2016,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.
@@ -53,6 +53,8 @@
 #include "gromacs/domdec/domdec_struct.h"
 #include "gromacs/utility/gmxmpi.h"
 
+#include "domdec_internal.h"
+
 
 /*! \brief Returns the MPI rank of the domain decomposition master rank */
 #define DDMASTERRANK(dd)   ((dd)->masterrank)
@@ -161,7 +163,7 @@ void dd_sendrecv2_rvec(const struct gmx_domdec_t gmx_unused *dd,
     rank_fw = dd->neighbor[ddimind][0];
     rank_bw = dd->neighbor[ddimind][1];
 
-    if (!dd->bSendRecv2)
+    if (!dd->comm->ddSettings.useSendRecv2)
     {
         /* Try to send and receive in two directions simultaneously.
          * Should be faster, especially on machines
index 4a43576e33067ff78ada75ff25d62ebd2d1f7907..88bc93be50d35f026301a1b0c507c73dff1987c4 100644 (file)
@@ -161,8 +161,6 @@ struct gmx_domdec_t { //NOLINT(clang-analyzer-optin.performance.Padding)
      */
     int                    nnodes       = 0;
     MPI_Comm               mpi_comm_all = MPI_COMM_NULL;
-    /* Use MPI_Sendrecv communication instead of non-blocking calls */
-    gmx_bool               bSendRecv2 = FALSE;
     /* The local DD cell index and rank */
     ivec                   ci         = { 0, 0, 0 };
     int                    rank       = 0;
index cf1adc3acf312965725af4ed3ad096a8ee135098..f8b37d908769dc32de6d40dec7aff4d07ab119e3 100644 (file)
@@ -582,7 +582,7 @@ static void check_index_consistency(gmx_domdec_t *dd,
 
     const int numAtomsInZones = dd->comm->atomRanges.end(DDAtomRanges::Type::Zones);
 
-    if (dd->comm->DD_debug > 1)
+    if (dd->comm->ddSettings.DD_debug > 1)
     {
         std::vector<int> have(natoms_sys);
         for (int a = 0; a < numAtomsInZones; a++)
@@ -732,12 +732,12 @@ static float dd_force_load(gmx_domdec_comm_t *comm)
 {
     float load;
 
-    if (comm->eFlop)
+    if (comm->ddSettings.eFlop)
     {
         load = comm->flop;
-        if (comm->eFlop > 1)
+        if (comm->ddSettings.eFlop > 1)
         {
-            load *= 1.0 + (comm->eFlop - 1)*(0.1*rand()/RAND_MAX - 0.05);
+            load *= 1.0 + (comm->ddSettings.eFlop - 1)*(0.1*rand()/RAND_MAX - 0.05);
         }
     }
     else
@@ -2715,7 +2715,7 @@ void print_dd_statistics(const t_commrec *cr, const t_inputrec *ir, FILE *fplog)
     }
     fprintf(fplog, "\n");
 
-    if (comm->bRecordLoad && EI_DYNAMICS(ir->eI))
+    if (comm->ddSettings.recordLoad && EI_DYNAMICS(ir->eI))
     {
         print_dd_load_av(fplog, cr->dd);
     }
@@ -2813,7 +2813,7 @@ void dd_partition_system(FILE                    *fplog,
     }
 
     /* Check if we have recorded loads on the nodes */
-    if (comm->bRecordLoad && dd_load_count(comm) > 0)
+    if (comm->ddSettings.recordLoad && dd_load_count(comm) > 0)
     {
         bCheckWhetherToTurnDlbOn = dd_dlb_get_should_check_whether_to_turn_dlb_on(dd);
 
@@ -3037,7 +3037,7 @@ void dd_partition_system(FILE                    *fplog,
     set_dd_cell_sizes(dd, &ddbox, dd->unitCellInfo.ddBoxIsDynamic, bMasterState, bDoDLB,
                       step, wcycle);
 
-    if (comm->nstDDDumpGrid > 0 && step % comm->nstDDDumpGrid == 0)
+    if (comm->ddSettings.nstDDDumpGrid > 0 && step % comm->ddSettings.nstDDDumpGrid == 0)
     {
         write_dd_grid_pdb("dd_grid", step, dd, state_local->box, &ddbox);
     }
@@ -3316,7 +3316,7 @@ void dd_partition_system(FILE                    *fplog,
 
     wallcycle_sub_stop(wcycle, ewcsDD_TOPOTHER);
 
-    if (comm->nstDDDump > 0 && step % comm->nstDDDump == 0)
+    if (comm->ddSettings.nstDDDump > 0 && step % comm->ddSettings.nstDDDump == 0)
     {
         dd_move_x(dd, state_local->box, state_local->x, nullWallcycle);
         write_dd_pdb("dd_dump", step, "dump", &top_global, cr,
@@ -3338,7 +3338,7 @@ void dd_partition_system(FILE                    *fplog,
         comm->master_cg_ddp_count = (bSortCG ? 0 : dd->ddp_count);
     }
 
-    if (comm->DD_debug > 0)
+    if (comm->ddSettings.DD_debug > 0)
     {
         /* Set the env var GMX_DD_DEBUG if you suspect corrupted indices */
         check_index_consistency(dd, top_global.natoms, ncg_mtop(&top_global),
index 51f3f34f282988aca906dacfd011e25afbbaea64..151cd3d833aa12034104620269612ac9a18d890c 100644 (file)
@@ -86,7 +86,7 @@ static inline int ddcginfo(const cginfo_mb_t *cginfo_mb,
 /*! \brief Returns the number of MD steps for which load has been recorded */
 static inline int dd_load_count(const gmx_domdec_comm_t *comm)
 {
-    return (comm->eFlop ? comm->flop_n : comm->cycl_n[ddCyclF]);
+    return (comm->ddSettings.eFlop ? comm->flop_n : comm->cycl_n[ddCyclF]);
 }
 
 /*! \brief Resize the state and f, if !=nullptr, to natoms */