2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2019,2020,2021, 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.
37 * \brief Declares functions for choosing the DD grid setup
39 * \author Berk Hess <hess@kth.se>
40 * \ingroup module_domdec
42 #ifndef GMX_DOMDEC_DOMDEC_SETUP_H
43 #define GMX_DOMDEC_DOMDEC_SETUP_H
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/gmxmpi.h"
59 struct MDModulesNotifiers;
61 class SeparatePmeRanksPermitted;
66 /*! \brief Returns the volume fraction of the system that is communicated */
67 real comm_box_frac(const gmx::IVec& dd_nc, real cutoff, const gmx_ddbox_t& ddbox);
70 * \brief Describes the DD grid setup
72 * This struct is for temporary use when choosing and initializing
73 * the domain decomposition grid.
77 //! The number of separate PME ranks, 0 if none or all ranks do PME
78 int numPmeOnlyRanks = 0;
79 //! The number of domains along each dimension
80 ivec numDomains = { 0, 0, 0 };
81 //! The number of dimensions which we decompose in domains
82 int numDDDimensions = 0;
83 //! The domain decomposition dimensions, the first numDDDimensions entries are used
84 ivec ddDimensions = { -1, -1, -1 };
87 /*! \brief Checks for ability to use separate PME ranks
89 * Disables automatic usage if:
90 * some MDModule could not use separate PME ranks,
91 * GPU setup is not compatible with separate PME ranks,
92 * user provided explicit DD grid
93 * or total number of ranks is not large enough to use PME ranks
95 gmx::SeparatePmeRanksPermitted checkForSeparatePmeRanks(const gmx::MDModulesNotifiers& notifiers,
96 const gmx::DomdecOptions& options,
97 int numRanksRequested,
98 bool useGpuForNonbonded,
101 /*! \brief Checks that requests for PP and PME ranks honor basic expectations
103 * Issues a fatal error if there are more PME ranks than PP, if the
104 * count of PP ranks has a prime factor that is too large to be likely
105 * to have good performance or PME-only ranks could not be used,
106 * but requested with -npme > 0 */
107 void checkForValidRankCountRequests(int numRanksRequested,
109 int numPmeRanksRequested,
110 const gmx::SeparatePmeRanksPermitted& separatePmeRanksPermitted,
111 bool checkForLargePrimeFactors);
113 /*! \brief Return the minimum cell size (in nm) required for DD */
114 real getDDGridSetupCellSizeLimit(const gmx::MDLogger& mdlog,
117 const t_inputrec& ir,
118 real systemInfoCellSizeLimit);
120 /*! \brief Determines the DD grid setup
122 * Either implements settings required by the user, or otherwise
123 * chooses estimated optimal number of separate PME ranks and DD grid
124 * cell setup, DD cell size limits, and the initial ddbox.
126 DDGridSetup getDDGridSetup(const gmx::MDLogger& mdlog,
128 MPI_Comm communicator,
129 int numRanksRequested,
130 const gmx::DomdecOptions& options,
131 const DDSettings& ddSettings,
132 const DDSystemInfo& systemInfo,
134 const gmx_mtop_t& mtop,
135 const t_inputrec& ir,
136 const gmx::SeparatePmeRanksPermitted& separatePmeRanksPermitted,
138 gmx::ArrayRef<const gmx::RVec> xGlobal,