Removed mdrun -gcom
authorMark Abraham <mark.j.abraham@gmail.com>
Sat, 13 Apr 2019 14:55:39 +0000 (16:55 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 16 Apr 2019 11:27:17 +0000 (13:27 +0200)
This was previously deprecated, and is now removed to make the
behaviour of mdrun simpler to understand and implement.

Renamed a function whose job was previously not to check a thing,
and is now clearly to compute something

Noted several TODOs to clean up behaviour related to nstcomm.

Refs #1925

Change-Id: I0b3a803fb209148e865957f796c871caef2f1fea

12 files changed:
docs/release-notes/2020/major/deprecated-functionality.rst
docs/release-notes/2020/major/removed-functionality.rst
docs/user-guide/mdrun-performance.rst
src/gromacs/domdec/domdec_internal.h
src/gromacs/domdec/partition.cpp
src/gromacs/gmxpreprocess/readir.cpp
src/gromacs/mdlib/md_support.cpp
src/gromacs/mdlib/md_support.h
src/gromacs/mdrun/legacymdrunoptions.h
src/gromacs/mdrun/md.cpp
src/gromacs/mdtypes/mdrunoptions.h
src/gromacs/timing/wallcycle.cpp

index d32f9be08cb1998aea3b1935c8d7e10c179c8f32..8dbaa88481da54b7d56a9922a6157177873c46ed 100644 (file)
@@ -21,13 +21,6 @@ These are thought to produce artefacts under some circumstances
 (unpublished results), were never well tested, are not widely used,
 and we need to simplify pdb2gmx.
 
-``gmx mdrun -gcom``
-""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
-This feature sometimes overrides the effects of various .mdp settings
-in a way that is difficult to understand and report. A user who wants
-to do communication between PP ranks less often should choose their
-``nst*`` mdp options accordingly.
-
 Benchmarking options only available with ``gmx benchmark``
 """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
 Options such as ``-confout``, ``-resethway``, ``-resetstep`` are not
index 264f6345e0aa84a69cb7ce11d5c7602301d20a14..7ceca019d4080e6edc63fede12b9e3caad0e7ead 100644 (file)
@@ -31,3 +31,11 @@ gmx morph
 """""""""
 The gmx morph tool was removed since it yields non-physical structures
 that can easily be done by a script.
+
+gmx mdrun -gcom
+"""""""""""""""
+
+This feature sometimes overrode the effects of various .mdp settings
+in a way that was difficult to understand and report. A user who wants
+to do communication between PP ranks less often should choose their
+``nst*`` mdp options accordingly.
index 5dc37bf6544cff25f72c11f83cc3eef94edb424b..20311c0476d264e604358c02b49660e7ca06303b 100644 (file)
@@ -758,22 +758,13 @@ cases.
     performance. If available, using ``-bonded gpu`` is expected
     to improve the ability of DLB to maximize performance.
 
-``-gcom``
-    During the simulation :ref:`gmx mdrun` must communicate between all ranks to
-    compute quantities such as kinetic energy. By default, this
-    happens whenever plausible, and is influenced by a lot of
-    :ref:`mdp options. <mdp-general>` The period between communication phases
-    must be a multiple of :mdp:`nstlist`, and defaults to
-    the minimum of :mdp:`nstcalcenergy` and :mdp:`nstlist`.
-    ``mdrun -gcom`` sets the number of steps that must elapse between
-    such communication phases, which can improve performance when
-    running on a lot of ranks. Note that this means that _e.g._
-    temperature coupling algorithms will
-    effectively remain at constant energy until the next
-    communication phase. :ref:`gmx mdrun` will always honor the
-    setting of ``mdrun -gcom``, by changing :mdp:`nstcalcenergy`,
-    :mdp:`nstenergy`, :mdp:`nstlog`, :mdp:`nsttcouple` and/or
-    :mdp:`nstpcouple` if necessary.
+During the simulation :ref:`gmx mdrun` must communicate between all
+PP ranks to compute quantities such as kinetic energy for log file
+reporting, or perhaps temperature coupling. By default, this happens
+whenever necessary to honor several :ref:`mdp options <mdp-general>`,
+so that the period between communication phases is the least common
+denominator of :mdp:`nstlist`, :mdp:`nstcalcenergy`,
+:mdp:`nsttcouple`, and :mdp:`nstpcouple`.
 
 Note that ``-tunepme`` has more effect when there is more than one
 :term:`node`, because the cost of communication for the PP and PME
@@ -1325,8 +1316,8 @@ Run setup
   :ref:`gmx tune_pme` can help automate this search.
 * For massively parallel runs (also ``gmx mdrun -multidir``), or with a slow
   network, global communication can become a bottleneck and you can reduce it
-  with ``gmx mdrun -gcom`` (note that this does affect the frequency of
-  temperature and pressure coupling).
+  by choosing larger periods for algorithms such as temperature and
+  pressure coupling).
 
 Checking and improving performance
 ^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^^
@@ -1345,7 +1336,7 @@ Checking and improving performance
   list to do more non-bonded computation to keep energy drift constant).
 
   * If ``Comm. energies`` takes a lot of time (a note will be printed in the log
-    file), increase nstcalcenergy or use ``mdrun -gcom``.
+    file), increase nstcalcenergy.
   * If all communication takes a lot of time, you might be running on too many
     cores, or you could try running combined MPI/OpenMP parallelization with 2
     or 4 OpenMP threads per MPI process.
index f5754a7f2fd72e8293c384a0ca962821999e4814..c2d9aef4c722a0114a5022d543e0af6310d1990d 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014,2015,2016,2017,2018, by the GROMACS development team, led by
+ * Copyright (c) 2014,2015,2016,2017,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.
@@ -453,8 +453,8 @@ struct gmx_domdec_comm_t // NOLINT (clang-analyzer-optin.performance.Padding)
     real PMELoadBal_max_cutoff;
 
     ivec tric_dir;                /**< tric_dir from \p gmx_ddbox_t is only stored here because dd_get_ns_ranges needs it */
-    rvec box0;                    /**< box lower corner, required with dim's without pbc and -gcom */
-    rvec box_size;                /**< box size, required with dim's without pbc and -gcom */
+    rvec box0;                    /**< box lower corner, required with dim's without pbc when avoiding communication */
+    rvec box_size;                /**< box size, required with dim's without pbc when avoiding communication */
 
     rvec cell_x0;                 /**< The DD cell lower corner, in triclinic space */
     rvec cell_x1;                 /**< The DD cell upper corner, in triclinic space */
index 28b3008347154979530f0b80963f26d70ffd6a43..789804bf756ad458f04d60cd21bcab3229decaa1 100644 (file)
@@ -3307,7 +3307,10 @@ void dd_partition_system(FILE                    *fplog,
         clearDDStateIndices(dd, dd->ncg_home, comm->atomRanges.numHomeAtoms());
         ncgindex_set = 0;
 
-        /* Avoid global communication for dim's without pbc and -gcom */
+        /* To avoid global communication, we do not recompute the extent
+         * of the system for dims without pbc. Therefore we need to copy
+         * the previously computed values when we do not communicate.
+         */
         if (!bNStGlobalComm)
         {
             copy_rvec(comm->box0, ddbox.box0    );
@@ -3319,7 +3322,7 @@ void dd_partition_system(FILE                    *fplog,
         bBoxChanged = TRUE;
         bRedist     = TRUE;
     }
-    /* For dim's without pbc and -gcom */
+    /* Copy needed for dim's without pbc when avoiding communication */
     copy_rvec(ddbox.box0, comm->box0    );
     copy_rvec(ddbox.box_size, comm->box_size);
 
index c4ef18259bfade2768fa6d75b0b02b40986af817..7774199656458608576e3a47a98bd00d7b191666 100644 (file)
@@ -880,12 +880,18 @@ void check_ir(const char *mdparin, t_inputrec *ir, t_gromppopts *opts,
     /* COMM STUFF */
     if (ir->nstcomm == 0)
     {
+        // TODO Change this behaviour. There should be exactly one way
+        // to turn off an algorithm.
         ir->comm_mode = ecmNO;
     }
     if (ir->comm_mode != ecmNO)
     {
         if (ir->nstcomm < 0)
         {
+            // TODO Such input was once valid. Now that we've been
+            // helpful for a few years, we should reject such input,
+            // lest we have to support every historical decision
+            // forever.
             warning(wi, "If you want to remove the rotation around the center of mass, you should set comm_mode = Angular instead of setting nstcomm < 0. nstcomm is modified to its absolute value");
             ir->nstcomm = abs(ir->nstcomm);
         }
index 0a77c4123ef4e80ddacc1daf2c587f79658e6896..51f09b661825a9d698bf1570cb77f60417494779 100644 (file)
@@ -319,21 +319,6 @@ void compute_globals(FILE *fplog, gmx_global_stat *gstat, t_commrec *cr, t_input
     }
 }
 
-/* check whether an 'nst'-style parameter p is a multiple of nst, and
-   set it to be one if not, with a warning. */
-static void check_nst_param(const gmx::MDLogger &mdlog,
-                            const char *desc_nst, int nst,
-                            const char *desc_p, int *p)
-{
-    if (*p > 0 && *p % nst != 0)
-    {
-        /* Round up to the next multiple of nst */
-        *p = ((*p)/nst + 1)*nst;
-        GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
-                "NOTE: %s changes %s to %d", desc_nst, desc_p, *p);
-    }
-}
-
 void setCurrentLambdasRerun(int64_t step, const t_lambda *fepvals,
                             const t_trxframe *rerun_fr, const double *lam0,
                             t_state *globalState)
@@ -446,14 +431,11 @@ static int lcd4(int i1, int i2, int i3, int i4)
     return nst;
 }
 
-int check_nstglobalcomm(const gmx::MDLogger &mdlog, int nstglobalcomm, t_inputrec *ir, const t_commrec * cr)
+int computeGlobalCommunicationPeriod(const gmx::MDLogger &mdlog,
+                                     t_inputrec          *ir,
+                                     const t_commrec     *cr)
 {
-    if (!EI_DYNAMICS(ir->eI))
-    {
-        nstglobalcomm = 1;
-    }
-
-    if (nstglobalcomm == -1)
+    int nstglobalcomm;
     {
         // Set up the default behaviour
         if (!(ir->nstcalcenergy > 0 ||
@@ -487,40 +469,9 @@ int check_nstglobalcomm(const gmx::MDLogger &mdlog, int nstglobalcomm, t_inputre
                                  ir->epc != epcNO ? ir->nstpcouple : 0);
         }
     }
-    else
-    {
-        // Check that the user's choice of mdrun -gcom will work
-        if (ir->nstlist > 0 &&
-            nstglobalcomm > ir->nstlist && nstglobalcomm % ir->nstlist != 0)
-        {
-            nstglobalcomm = (nstglobalcomm / ir->nstlist)*ir->nstlist;
-            GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
-                    "WARNING: nstglobalcomm is larger than nstlist, but not a multiple, setting it to %d",
-                    nstglobalcomm);
-        }
-        if (ir->nstcalcenergy > 0)
-        {
-            check_nst_param(mdlog, "-gcom", nstglobalcomm,
-                            "nstcalcenergy", &ir->nstcalcenergy);
-        }
-        if (ir->etc != etcNO && ir->nsttcouple > 0)
-        {
-            check_nst_param(mdlog, "-gcom", nstglobalcomm,
-                            "nsttcouple", &ir->nsttcouple);
-        }
-        if (ir->epc != epcNO && ir->nstpcouple > 0)
-        {
-            check_nst_param(mdlog, "-gcom", nstglobalcomm,
-                            "nstpcouple", &ir->nstpcouple);
-        }
-
-        check_nst_param(mdlog, "-gcom", nstglobalcomm,
-                        "nstenergy", &ir->nstenergy);
-
-        check_nst_param(mdlog, "-gcom", nstglobalcomm,
-                        "nstlog", &ir->nstlog);
-    }
 
+    // TODO change this behaviour. Instead grompp should print
+    // a (performance) note and mdrun should not change ir.
     if (ir->comm_mode != ecmNO && ir->nstcomm < nstglobalcomm)
     {
         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
index e218d383c0162f526e8c60f36940c99e385e6c12..87029dcd3cebe131b732e4e33849c8688f8a52ba 100644 (file)
@@ -93,11 +93,10 @@ class SimulationSignaller;
 
 /*! \brief Return the number of steps that will take place between
  * intra-simulation communications, given the constraints of the
- * inputrec and the value of mdrun -gcom. */
-int check_nstglobalcomm(const gmx::MDLogger &mdlog,
-                        int                  nstglobalcomm,
-                        t_inputrec          *ir,
-                        const t_commrec    * cr);
+ * inputrec. */
+int computeGlobalCommunicationPeriod(const gmx::MDLogger &mdlog,
+                                     t_inputrec          *ir,
+                                     const t_commrec     *cr);
 
 /*! \brief Return true if the \p value is equal across the set of multi-simulations
  *
index 8839d3670007c03c20b9d56552c829ab1a92bb0c..dca811f926e45419cab797abea4e3bce714b4bb1 100644 (file)
@@ -161,7 +161,7 @@ class LegacyMdrunOptions
 
         ImdOptions       &imdOptions = mdrunOptions.imdOptions;
 
-        t_pargs           pa[48] = {
+        t_pargs           pa[47] = {
 
             { "-dd",      FALSE, etRVEC, {&realddxyz},
               "Domain decomposition grid, 0 is optimize" },
@@ -212,8 +212,6 @@ class LegacyMdrunOptions
               "HIDDENA string containing a vector of the relative sizes in the z "
               "direction of the corresponding DD cells. Only effective with static "
               "load balancing." },
-            { "-gcom",    FALSE, etINT, {&mdrunOptions.globalCommunicationInterval},
-              "Global communication frequency" },
             { "-nb",      FALSE, etENUM, {nbpu_opt_choices},
               "Calculate non-bonded interactions on" },
             { "-nstlist", FALSE, etINT, {&nstlist_cmdline},
index 307e3182f040b87d889a69147d835697566069cb..16e7e890308eb4ed7aa0804225578060beebe6b5 100644 (file)
@@ -224,9 +224,8 @@ void gmx::Integrator::do_md()
     bTrotter = (EI_VV(ir->eI) && (inputrecNptTrotter(ir) || inputrecNphTrotter(ir) || inputrecNvtTrotter(ir)));
 
     const bool bRerunMD      = false;
-    int        nstglobalcomm = mdrunOptions.globalCommunicationInterval;
 
-    nstglobalcomm   = check_nstglobalcomm(mdlog, nstglobalcomm, ir, cr);
+    int        nstglobalcomm = computeGlobalCommunicationPeriod(mdlog, ir, cr);
     bGStatEveryStep = (nstglobalcomm == 1);
 
     SimulationGroups                  *groups = &top_global->groups;
index 8e7f626310950117209bafa8cdfbbdb7760091fe..4bc86b68d8aa165d7168fc56a25771e34bb20040 100644 (file)
@@ -106,8 +106,6 @@ struct MdrunOptions
     gmx_bool            rerun = FALSE;
     //! Re-construct virual sites durin a rerun simulation
     gmx_bool            rerunConstructVsites = FALSE;
-    //! Request to do global communication at this interval in steps, -1 is determine from inputrec (deprecated).
-    int                 globalCommunicationInterval = -1;
     //! Try to make the simulation binary reproducible
     gmx_bool            reproducible = FALSE;
     //! Write confout.gro at the end of the run
index 7aeb422a1285a1261e72dee278440d6eb0435a4a..34a00a67ebcb02ee603bbe3a6516e5bed4ada7b8 100644 (file)
@@ -1073,7 +1073,7 @@ void wallcycle_print(FILE *fplog, const gmx::MDLogger &mdlog, int nnodes, int np
     {
         GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
                 "NOTE: %d %% of the run time was spent communicating energies,\n"
-                "      you might want to use the -gcom option of mdrun\n",
+                "      you might want to increase some nst* mdp options\n",
                 gmx::roundToInt(100*cyc_sum[ewcMoveE]/tot));
     }
 }