Move DomdecOptions into its own header
authorMark Abraham <mark.j.abraham@gmail.com>
Wed, 13 Mar 2019 09:26:21 +0000 (10:26 +0100)
committerBerk Hess <hess@kth.se>
Thu, 28 Mar 2019 08:56:30 +0000 (09:56 +0100)
This reduces compilation coupling and prepares for implementing proper
Options support into mdrun. Moved content into gmx namespace and
replaced gmx_bool. Renamed nr enum members to Count for future
compatibility with the new enumeration helpers.

Refs #2877

Change-Id: Ib5f5889111f858579f594f4c76f3aa8c222403b2

src/gromacs/domdec/domdec.cpp
src/gromacs/domdec/domdec.h
src/gromacs/domdec/options.h [new file with mode: 0644]
src/gromacs/mdrun/legacymdrunoptions.h
src/gromacs/mdrun/runner.h
src/programs/mdrun/mdrun.cpp

index 2090b4697ea4dd8358a460b5073a9eeec343c1e7..2191e995877b916d6b0d6bbe3b41edc8826c01ac 100644 (file)
@@ -55,6 +55,7 @@
 #include "gromacs/domdec/dlbtiming.h"
 #include "gromacs/domdec/domdec_network.h"
 #include "gromacs/domdec/ga2la.h"
+#include "gromacs/domdec/options.h"
 #include "gromacs/domdec/partition.h"
 #include "gromacs/gmxlib/network.h"
 #include "gromacs/gmxlib/nrnb.h"
 #include "redistribute.h"
 #include "utility.h"
 
+// TODO remove this when moving domdec into gmx namespace
+using gmx::DdRankOrder;
+using gmx::DlbOption;
+using gmx::DomdecOptions;
+
 static const char *edlbs_names[int(DlbState::nr)] = { "off", "auto", "locked", "on", "on" };
 
 /* The size per atom group of the cggl_flag buffer in gmx_domdec_comm_t */
index 5ea2ef1f8ee594338c2007438559e0131fa260b7..2ed78d2cbeafc93cf1a7abd8d81c6a91f8e2ea54 100644 (file)
@@ -86,8 +86,9 @@ class t_state;
 namespace gmx
 {
 class MDLogger;
-struct MdrunOptions;
 class LocalAtomSetManager;
+struct DomdecOptions;
+struct MdrunOptions;
 } // namespace
 
 /*! \brief Returns the global topology atom number belonging to local atom index i.
@@ -146,63 +147,11 @@ 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 The options for the domain decomposition MPI task ordering. */
-enum class DdRankOrder
-{
-    select,     //!< First value (needed to cope with command-line parsing)
-    interleave, //!< Interleave the PP and PME ranks
-    pp_pme,     //!< First all PP ranks, all PME rank at the end
-    cartesian,  //!< Use Cartesian communicators for PP, PME and PP-PME
-    nr          //!< The number of options
-};
-
-/*! \brief The options for the dynamic load balancing. */
-enum class DlbOption
-{
-    select,           //!< First value (needed to cope with command-line parsing)
-    turnOnWhenUseful, //!< Turn on DLB when we think it would improve performance
-    no,               //!< Never turn on DLB
-    yes,              //!< Turn on DLB from the start and keep it on
-    nr                //!< The number of options
-};
-
-/*! \libinternal \brief Structure containing all (command line) options for the domain decomposition */
-struct DomdecOptions
-{
-    //! If true, check that all bonded interactions have been assigned to exactly one domain/rank.
-    gmx_bool          checkBondedInteractions = TRUE;
-    //! If true, don't communicate all atoms between the non-bonded cut-off and the larger bonded cut-off, but only those that have non-local bonded interactions. This significantly reduces the communication volume.
-    gmx_bool          useBondedCommunication = TRUE;
-    //! The domain decomposition grid cell count, 0 means let domdec choose based on the number of ranks.
-    ivec              numCells = {0};
-    //! The number of separate PME ranks requested, -1 = auto.
-    int               numPmeRanks = -1;
-    //! Ordering of the PP and PME ranks, values from enum above.
-    DdRankOrder       rankOrder = DdRankOrder::interleave;
-    //! The minimum communication range, used for extended the communication range for bonded interactions (nm).
-    real              minimumCommunicationRange = 0;
-    //! Communication range for atom involved in constraints (P-LINCS) (nm).
-    real              constraintCommunicationRange = 0;
-    //! Dynamic load balancing option, values from enum above.
-    DlbOption         dlbOption = DlbOption::turnOnWhenUseful;
-    /*! \brief Fraction in (0,1) by whose reciprocal the initial
-     * DD cell size will be increased in order to provide a margin
-     * in which dynamic load balancing can act, while preserving
-     * the minimum cell size. */
-    real              dlbScaling = 0.8;
-    //! String containing a vector of the relative sizes in the x direction of the corresponding DD cells.
-    const char       *cellSizeX = nullptr;
-    //! String containing a vector of the relative sizes in the y direction of the corresponding DD cells.
-    const char       *cellSizeY = nullptr;
-    //! String containing a vector of the relative sizes in the z direction of the corresponding DD cells.
-    const char       *cellSizeZ = nullptr;
-};
-
 /*! \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,
                           t_commrec                      *cr,
-                          const DomdecOptions            &options,
+                          const gmx::DomdecOptions       &options,
                           const gmx::MdrunOptions        &mdrunOptions,
                           const gmx_mtop_t               *mtop,
                           const t_inputrec               *ir,
diff --git a/src/gromacs/domdec/options.h b/src/gromacs/domdec/options.h
new file mode 100644 (file)
index 0000000..d5e06b7
--- /dev/null
@@ -0,0 +1,110 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \libinternal \file
+ *
+ * \brief This file declares command-line options for mdrun related to
+ * domain decomposition.
+ *
+ * \author Berk Hess <hess@kth.se>
+ * \inlibraryapi
+ * \ingroup module_domdec
+ */
+
+#ifndef GMX_DOMDEC_OPTIONS_H
+#define GMX_DOMDEC_OPTIONS_H
+
+#include <vector>
+
+#include "gromacs/math/vectypes.h"
+#include "gromacs/utility/real.h"
+
+namespace gmx
+{
+
+/*! \brief The options for the domain decomposition MPI task ordering. */
+enum class DdRankOrder
+{
+    select,     //!< First value (needed to cope with command-line parsing)
+    interleave, //!< Interleave the PP and PME ranks
+    pp_pme,     //!< First all PP ranks, all PME rank at the end
+    cartesian,  //!< Use Cartesian communicators for PP, PME and PP-PME
+    Count       //!< The number of options
+};
+
+/*! \brief The options for the dynamic load balancing. */
+enum class DlbOption
+{
+    select,           //!< First value (needed to cope with command-line parsing)
+    turnOnWhenUseful, //!< Turn on DLB when we think it would improve performance
+    no,               //!< Never turn on DLB
+    yes,              //!< Turn on DLB from the start and keep it on
+    Count             //!< The number of options
+};
+
+/*! \libinternal \brief Structure containing all (command line) options for the domain decomposition */
+struct DomdecOptions
+{
+    //! If true, check that all bonded interactions have been assigned to exactly one domain/rank.
+    bool              checkBondedInteractions = true;
+    //! If true, don't communicate all atoms between the non-bonded cut-off and the larger bonded cut-off, but only those that have non-local bonded interactions. This significantly reduces the communication volume.
+    bool              useBondedCommunication = true;
+    //! The domain decomposition grid cell count, 0 means let domdec choose based on the number of ranks.
+    ivec              numCells = {0};
+    //! The number of separate PME ranks requested, -1 = auto.
+    int               numPmeRanks = -1;
+    //! Ordering of the PP and PME ranks, values from enum above.
+    DdRankOrder       rankOrder = DdRankOrder::interleave;
+    //! The minimum communication range, used for extended the communication range for bonded interactions (nm).
+    real              minimumCommunicationRange = 0;
+    //! Communication range for atom involved in constraints (P-LINCS) (nm).
+    real              constraintCommunicationRange = 0;
+    //! Dynamic load balancing option, values from enum above.
+    DlbOption         dlbOption = DlbOption::turnOnWhenUseful;
+    /*! \brief Fraction in (0,1) by whose reciprocal the initial
+     * DD cell size will be increased in order to provide a margin
+     * in which dynamic load balancing can act, while preserving
+     * the minimum cell size. */
+    real              dlbScaling = 0.8;
+    //! String containing a vector of the relative sizes in the x direction of the corresponding DD cells.
+    const char       *cellSizeX = nullptr;
+    //! String containing a vector of the relative sizes in the y direction of the corresponding DD cells.
+    const char       *cellSizeY = nullptr;
+    //! String containing a vector of the relative sizes in the z direction of the corresponding DD cells.
+    const char       *cellSizeZ = nullptr;
+};
+
+} // namespace gmx
+
+#endif
index 69602710e85733e56145047f4556a07e2bf8de53..8839d3670007c03c20b9d56552c829ab1a92bb0c 100644 (file)
@@ -51,7 +51,7 @@
 
 #include "gromacs/commandline/filenm.h"
 #include "gromacs/commandline/pargs.h"
-#include "gromacs/domdec/domdec.h"
+#include "gromacs/domdec/options.h"
 #include "gromacs/hardware/hw_info.h"
 #include "gromacs/mdrun/logging.h"
 #include "gromacs/mdtypes/mdrunoptions.h"
@@ -140,10 +140,10 @@ class LegacyMdrunOptions
 
         /*! \brief Command line options, defaults, docs and storage for them to fill. */
         /*! \{ */
-        rvec              realddxyz                                               = {0, 0, 0};
-        const char       *ddrank_opt_choices[static_cast<int>(DdRankOrder::nr)+1] =
+        rvec              realddxyz                                                  = {0, 0, 0};
+        const char       *ddrank_opt_choices[static_cast<int>(DdRankOrder::Count)+1] =
         { nullptr, "interleave", "pp_pme", "cartesian", nullptr };
-        const char       *dddlb_opt_choices[static_cast<int>(DlbOption::nr)+1] =
+        const char       *dddlb_opt_choices[static_cast<int>(DlbOption::Count)+1] =
         { nullptr, "auto", "no", "yes", nullptr };
         const char       *thread_aff_opt_choices[threadaffNR+1] =
         { nullptr, "auto", "on", "off", nullptr };
index deb2fb1ca0e3aba811c96fa8e1fc17618328a3b9..a67dea3127c51dc5acbc57524628b27abf98077a 100644 (file)
@@ -49,7 +49,7 @@
 
 #include "gromacs/commandline/filenm.h"
 #include "gromacs/compat/pointers.h"
-#include "gromacs/domdec/domdec.h"
+#include "gromacs/domdec/options.h"
 #include "gromacs/hardware/hw_info.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/mdtypes/mdrunoptions.h"
index c8015149612cb417b914a50c6b85da97e6b750f2..d638bcf9600976b29c9cfb6ad9b833384a5080c5 100644 (file)
@@ -56,7 +56,7 @@
 
 #include "gromacs/commandline/pargs.h"
 #include "gromacs/compat/pointers.h"
-#include "gromacs/domdec/domdec.h"
+#include "gromacs/domdec/options.h"
 #include "gromacs/fileio/gmxfio.h"
 #include "gromacs/gmxlib/network.h"
 #include "gromacs/mdrun/legacymdrunoptions.h"