d982ed328a5d20b7e98e92a7fa51faecbeaae8fc
[alexxy/gromacs.git] / src / gromacs / domdec / domdec_setup.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2019,2020, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*! \internal \file
36  *
37  * \brief Declares functions for choosing the DD grid setup
38  *
39  * \author Berk Hess <hess@kth.se>
40  * \ingroup module_domdec
41  */
42 #ifndef GMX_DOMDEC_DOMDEC_SETUP_H
43 #define GMX_DOMDEC_DOMDEC_SETUP_H
44
45 #include "gromacs/math/vec.h"
46 #include "gromacs/utility/gmxmpi.h"
47
48 enum class DDRole;
49 struct DDSettings;
50 struct DDSystemInfo;
51 struct gmx_ddbox_t;
52 struct gmx_mtop_t;
53 struct t_commrec;
54 struct t_inputrec;
55
56 namespace gmx
57 {
58 struct DomdecOptions;
59 class MDLogger;
60 template<typename T>
61 class ArrayRef;
62 } // namespace gmx
63
64 /*! \brief Returns the volume fraction of the system that is communicated */
65 real comm_box_frac(const gmx::IVec& dd_nc, real cutoff, const gmx_ddbox_t& ddbox);
66
67 /*! \internal
68  * \brief Describes the DD grid setup
69  *
70  * This struct is for temporary use when choosing and initializing
71  * the domain decomposition grid.
72  */
73 struct DDGridSetup
74 {
75     //! The number of separate PME ranks, 0 if none or all ranks do PME
76     int numPmeOnlyRanks = 0;
77     //! The number of domains along each dimension
78     ivec numDomains = { 0, 0, 0 };
79     //! The number of dimensions which we decompose in domains
80     int numDDDimensions = 0;
81     //! The domain decomposition dimensions, the first numDDDimensions entries are used
82     ivec ddDimensions = { -1, -1, -1 };
83 };
84
85 /*! \brief Checks that requests for PP and PME ranks honor basic expectations
86  *
87  * Issues a fatal error if there are more PME ranks than PP, or if the
88  * count of PP ranks has a prime factor that is too large to be likely
89  * to have good performance. */
90 void checkForValidRankCountRequests(int  numRanksRequested,
91                                     bool usingPme,
92                                     int  numPmeRanksRequested,
93                                     bool checkForLargePrimeFactors);
94
95 /*! \brief Return the minimum cell size (in nm) required for DD */
96 real getDDGridSetupCellSizeLimit(const gmx::MDLogger& mdlog,
97                                  bool                 bDynLoadBal,
98                                  real                 dlb_scale,
99                                  const t_inputrec&    ir,
100                                  real                 systemInfoCellSizeLimit);
101
102 /*! \brief Determines the DD grid setup
103  *
104  * Either implements settings required by the user, or otherwise
105  * chooses estimated optimal number of separate PME ranks and DD grid
106  * cell setup, DD cell size limits, and the initial ddbox.
107  */
108 DDGridSetup getDDGridSetup(const gmx::MDLogger&           mdlog,
109                            DDRole                         ddRole,
110                            MPI_Comm                       communicator,
111                            int                            numRanksRequested,
112                            const gmx::DomdecOptions&      options,
113                            const DDSettings&              ddSettings,
114                            const DDSystemInfo&            systemInfo,
115                            real                           cellSizeLimit,
116                            const gmx_mtop_t&              mtop,
117                            const t_inputrec&              ir,
118                            const matrix                   box,
119                            gmx::ArrayRef<const gmx::RVec> xGlobal,
120                            gmx_ddbox_t*                   ddbox);
121
122 #endif