Make init_dires independent of t_commrec
authorPascal Merz <pascal.merz@me.com>
Wed, 22 Apr 2020 17:21:14 +0000 (17:21 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 22 Apr 2020 17:21:14 +0000 (17:21 +0000)
init_disres was requesting a full pointer to the commrec,
but only uses a single communicator and checks for master rank
and whether the run is parallel. This information is now passed
in explicitly, simplifying the planned splitting of t_commrec.

Note that passing a nullptr for commrec was (mis)used by
gmx_disre only - effectively signalling that init_disres was
called from an analysis tool and not from mdrun. This has
been made explicit.

Refs #2395

src/gromacs/gmxana/gmx_disre.cpp
src/gromacs/listed_forces/disre.cpp
src/gromacs/listed_forces/disre.h
src/gromacs/mdrun/runner.cpp
src/gromacs/mdtypes/commrec.h

index baee36ec07d6d377adaf13f770c551105e146a6e..6a84cc8a81991e37e0f65af822decbee3217282a 100644 (file)
@@ -60,6 +60,7 @@
 #include "gromacs/math/vec.h"
 #include "gromacs/mdlib/force.h"
 #include "gromacs/mdlib/mdatoms.h"
+#include "gromacs/mdtypes/commrec.h"
 #include "gromacs/mdtypes/fcdata.h"
 #include "gromacs/mdtypes/inputrec.h"
 #include "gromacs/mdtypes/md_enums.h"
@@ -799,7 +800,8 @@ int gmx_disre(int argc, char* argv[])
     }
 
     ir->dr_tau = 0.0;
-    init_disres(fplog, topInfo.mtop(), ir, nullptr, nullptr, &fcd, nullptr, FALSE);
+    init_disres(fplog, topInfo.mtop(), ir, DisResRunMode::AnalysisTool, DDRole::Master,
+                NumRanks::Single, MPI_COMM_NULL, nullptr, &fcd, nullptr, FALSE);
 
     int natoms = read_first_x(oenv, &status, ftp2fn(efTRX, NFILE, fnm), &t, &x, box);
     snew(f, 5 * natoms);
index 9eebef862dc5433231c83a1159b909608716e5f6..1d66f2ba52ac77a3bb5fd3bb164e2823a56230e6 100644 (file)
 void init_disres(FILE*                 fplog,
                  const gmx_mtop_t*     mtop,
                  t_inputrec*           ir,
-                 const t_commrec*      cr,
+                 DisResRunMode         disResRunMode,
+                 DDRole                ddRole,
+                 NumRanks              numRanks,
+                 MPI_Comm              communicator,
                  const gmx_multisim_t* ms,
                  t_fcdata*             fcd,
                  t_state*              state,
@@ -116,10 +119,10 @@ void init_disres(FILE*                 fplog,
     {
         /* We store the r^-6 time averages in an array that is indexed
          * with the local disres iatom index, so this doesn't work with DD.
-         * Note that DD is not initialized yet here, so we check for PAR(cr),
+         * Note that DD is not initialized yet here, so we check that we are on multiple ranks,
          * but there are probably also issues with e.g. NM MPI parallelization.
          */
-        if (cr && PAR(cr))
+        if ((disResRunMode == DisResRunMode::MDRun) && (numRanks == NumRanks::Multiple))
         {
             gmx_fatal(FARGS,
                       "Time-averaged distance restraints are not supported with MPI "
@@ -168,7 +171,7 @@ void init_disres(FILE*                 fplog,
         }
     }
 
-    if (cr && PAR(cr) && ir->nstdisreout > 0)
+    if ((disResRunMode == DisResRunMode::MDRun) && (numRanks == NumRanks::Multiple) && ir->nstdisreout > 0)
     {
         /* With DD we currently only have local pair information available */
         gmx_fatal(FARGS,
@@ -218,7 +221,7 @@ void init_disres(FILE*                 fplog,
     dd->Rtav_6 = &(dd->Rt_6[dd->nres]);
 
     ptr = getenv("GMX_DISRE_ENSEMBLE_SIZE");
-    if (cr && ms != nullptr && ptr != nullptr && !bIsREMD)
+    if ((disResRunMode == DisResRunMode::MDRun) && ms != nullptr && ptr != nullptr && !bIsREMD)
     {
 #if GMX_MPI
         dd->nsystems = 0;
@@ -230,11 +233,11 @@ void init_disres(FILE*                 fplog,
         /* This check is only valid on MASTER(cr), so probably
          * ensemble-averaged distance restraints are broken on more
          * than one processor per simulation system. */
-        if (MASTER(cr))
+        if (ddRole == DDRole::Master)
         {
             check_multi_int(fplog, ms, dd->nsystems, "the number of systems per ensemble", FALSE);
         }
-        gmx_bcast(sizeof(int), &dd->nsystems, cr->mpi_comm_mysim);
+        gmx_bcast(sizeof(int), &dd->nsystems, communicator);
 
         /* We use to allow any value of nsystems which was a divisor
          * of ms->nsim. But this required an extra communicator which
@@ -285,7 +288,7 @@ void init_disres(FILE*                 fplog,
          * checks from appropriate processes (since check_multi_int is
          * too broken to check whether the communication will
          * succeed...) */
-        if (cr && ms && dd->nsystems > 1 && MASTER(cr))
+        if ((disResRunMode == DisResRunMode::MDRun) && ms && dd->nsystems > 1 && (ddRole == DDRole::Master))
         {
             check_multi_int(fplog, ms, fcd->disres.nres, "the number of distance restraints", FALSE);
         }
index 2b4adc3e01cad3665fc4ee96f055c03df0efa527..43983bfb8271ce73fe183b9f553e3aa86cc6bcce 100644 (file)
@@ -49,6 +49,7 @@
 
 #include "gromacs/topology/ifunc.h"
 #include "gromacs/utility/basedefinitions.h"
+#include "gromacs/utility/gmxmpi.h"
 
 struct gmx_mtop_t;
 struct gmx_multisim_t;
@@ -57,6 +58,15 @@ struct t_commrec;
 struct t_inputrec;
 struct t_pbc;
 class t_state;
+enum class DDRole;
+enum class NumRanks;
+
+//! Whether distance restraints are called from mdrun or from an analysis tool
+enum class DisResRunMode
+{
+    MDRun,
+    AnalysisTool
+};
 
 /*! \brief
  * Initiates *fcd data.
@@ -71,7 +81,10 @@ class t_state;
 void init_disres(FILE*                 fplog,
                  const gmx_mtop_t*     mtop,
                  t_inputrec*           ir,
-                 const t_commrec*      cr,
+                 DisResRunMode         disResRunMode,
+                 DDRole                ddRole,
+                 NumRanks              numRanks,
+                 MPI_Comm              communicator,
                  const gmx_multisim_t* ms,
                  t_fcdata*             fcd,
                  t_state*              state,
index 0656ac7bfe3924a70be89f28ac5677df892abc2d..648d220aa684718d470f1a8884a87e4d10b36f0f 100644 (file)
@@ -1037,7 +1037,9 @@ int Mdrunner::mdrunner()
     snew(fcd, 1);
 
     /* This needs to be called before read_checkpoint to extend the state */
-    init_disres(fplog, &mtop, inputrec, cr, ms, fcd, globalState.get(), replExParams.exchangeInterval > 0);
+    init_disres(fplog, &mtop, inputrec, DisResRunMode::MDRun, MASTER(cr) ? DDRole::Master : DDRole::Agent,
+                PAR(cr) ? NumRanks::Multiple : NumRanks::Single, cr->mpi_comm_mysim, ms, fcd,
+                globalState.get(), replExParams.exchangeInterval > 0);
 
     init_orires(fplog, &mtop, inputrec, cr, ms, globalState.get(), &(fcd->orires));
 
index 0f00d52528004975d30bc02172d06186a158cf7c..804b9043d17c856bf5e28749397c3f011caad458 100644 (file)
@@ -50,6 +50,20 @@ struct gmx_domdec_t;
 #define DUTY_PP (1U << 0U)
 #define DUTY_PME (1U << 1U)
 
+//! Whether the current DD role is master or slave
+enum class DDRole
+{
+    Master,
+    Agent
+};
+
+//! Whether one or more ranks are used
+enum class NumRanks
+{
+    Single,
+    Multiple
+};
+
 typedef struct
 {
     int      bUse;