Update instructions in containers.rst
[gromacs.git] / src / gromacs / domdec / options.h
blob0ed617f677031468c7326e6ce822158c4b60d6ca
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \libinternal \file
37 * \brief This file declares command-line options for mdrun related to
38 * domain decomposition.
40 * \author Berk Hess <hess@kth.se>
41 * \inlibraryapi
42 * \ingroup module_domdec
45 #ifndef GMX_DOMDEC_OPTIONS_H
46 #define GMX_DOMDEC_OPTIONS_H
48 #include <vector>
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/real.h"
53 namespace gmx
56 /*! \brief The options for the domain decomposition MPI task ordering. */
57 enum class DdRankOrder
59 select, //!< First value (needed to cope with command-line parsing)
60 interleave, //!< Interleave the PP and PME ranks
61 pp_pme, //!< First all PP ranks, all PME rank at the end
62 cartesian, //!< Use Cartesian communicators for PP, PME and PP-PME
63 Count //!< The number of options
66 /*! \brief The options for the dynamic load balancing. */
67 enum class DlbOption
69 select, //!< First value (needed to cope with command-line parsing)
70 turnOnWhenUseful, //!< Turn on DLB when we think it would improve performance
71 no, //!< Never turn on DLB
72 yes, //!< Turn on DLB from the start and keep it on
73 Count //!< The number of options
76 /*! \libinternal \brief Structure containing all (command line) options for the domain decomposition */
77 struct DomdecOptions
79 //! If true, check that all bonded interactions have been assigned to exactly one domain/rank.
80 bool checkBondedInteractions = true;
81 //! 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.
82 bool useBondedCommunication = true;
83 //! The domain decomposition grid cell count, 0 means let domdec choose based on the number of ranks.
84 ivec numCells = { 0 };
85 //! The number of separate PME ranks requested, -1 = auto.
86 int numPmeRanks = -1;
87 //! Ordering of the PP and PME ranks, values from enum above.
88 DdRankOrder rankOrder = DdRankOrder::interleave;
89 //! The minimum communication range, used for extended the communication range for bonded interactions (nm).
90 real minimumCommunicationRange = 0;
91 //! Communication range for atom involved in constraints (P-LINCS) (nm).
92 real constraintCommunicationRange = 0;
93 //! Dynamic load balancing option, values from enum above.
94 DlbOption dlbOption = DlbOption::turnOnWhenUseful;
95 /*! \brief Fraction in (0,1) by whose reciprocal the initial
96 * DD cell size will be increased in order to provide a margin
97 * in which dynamic load balancing can act, while preserving
98 * the minimum cell size. */
99 real dlbScaling = 0.8;
100 //! String containing a vector of the relative sizes in the x direction of the corresponding DD cells.
101 const char* cellSizeX = nullptr;
102 //! String containing a vector of the relative sizes in the y direction of the corresponding DD cells.
103 const char* cellSizeY = nullptr;
104 //! String containing a vector of the relative sizes in the z direction of the corresponding DD cells.
105 const char* cellSizeZ = nullptr;
108 } // namespace gmx
110 #endif