Enable TPI with the Verlet cut-off scheme
[gromacs.git] / src / gromacs / taskassignment / resourcedivision.h
blob2bea6245a03d20a296b95f618cec1fb47c5b6c25
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,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
36 * \brief Declares utility functionality for dividing resources and
37 * checking for consistency and usefulness.
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \ingroup module_taskassignment
41 * \inlibraryapi
44 #ifndef GMX_TASKASSIGNMENT_RESOURCEDIVISION_H
45 #define GMX_TASKASSIGNMENT_RESOURCEDIVISION_H
47 #include <cstdio>
49 #include <vector>
51 #include "gromacs/ewald/pme.h"
52 #include "gromacs/utility/basedefinitions.h"
54 struct gmx_hw_info_t;
55 struct gmx_hw_opt_t;
56 struct gmx_mtop_t;
57 struct gmx_multisim_t;
58 struct t_commrec;
59 struct t_inputrec;
61 namespace gmx
63 class HardwareTopology;
64 class MDLogger;
65 class PhysicalNodeCommunicator;
68 /*! \brief Return the number of threads to use for thread-MPI based on how many
69 * were requested, which algorithms we're using,
70 * and how many particles there are.
71 * At the point we have already called check_and_update_hw_opt.
72 * Thus all options should be internally consistent and consistent
73 * with the hardware, except that ntmpi could be larger than number of GPUs.
74 * If necessary, this function will modify hw_opt->nthreads_omp.
76 int get_nthreads_mpi(const gmx_hw_info_t *hwinfo,
77 gmx_hw_opt_t *hw_opt,
78 const std::vector<int> &gpuIdsToUse,
79 bool nonbondedOnGpu,
80 bool pmeOnGpu,
81 const t_inputrec *inputrec,
82 const gmx_mtop_t *mtop,
83 const gmx::MDLogger &mdlog,
84 bool doMembed);
86 /*! \brief Check if the number of OpenMP threads is within reasonable range
87 * considering the hardware used. This is a crude check, but mainly
88 * intended to catch cases where the user starts 1 MPI rank per hardware
89 * thread or 1 rank per physical node.
90 * With a sub-optimal setup a note is printed to fplog and stderr when
91 * bNtOmpSet==TRUE; with bNtOptOptionSet==FALSE a fatal error is issued.
92 * This function should be called after thread-MPI and OpenMP are set up.
94 void check_resource_division_efficiency(const gmx_hw_info_t *hwinfo,
95 bool willUsePhysicalGpu,
96 gmx_bool bNtOmpOptionSet,
97 t_commrec *cr,
98 const gmx::MDLogger &mdlog);
100 /*! \brief Checks what our hardware options are based on how Gromacs was compiled
101 * and user-set options
103 * \param[in] mdlog Logger
104 * \param[in, out] hw_opt Hardware-related and threading options
105 * \param[in] isSimulationMasterRank
106 * \param[in] nPmeRanks Number of PME ranks
107 * \param[in] inputrec The input record, should point to a valid object when \p isSimulationMasterRank = true
108 * */
109 void checkAndUpdateHardwareOptions(const gmx::MDLogger &mdlog,
110 gmx_hw_opt_t *hw_opt,
111 bool isSimulationMasterRank,
112 int nPmeRanks,
113 const t_inputrec *inputrec);
115 /*! \brief Check, and if necessary update, the number of OpenMP threads requested
117 * Should be called when we know the MPI rank count and PME run mode.
119 void checkAndUpdateRequestedNumOpenmpThreads(gmx_hw_opt_t *hw_opt,
120 const gmx_hw_info_t &hwinfo,
121 const t_commrec *cr,
122 const gmx_multisim_t *ms,
123 int numRanksOnThisNode,
124 PmeRunMode pmeRunMode,
125 const gmx_mtop_t &mtop,
126 const t_inputrec &inputrec);
128 namespace gmx
131 /*! \brief Warns for oversubscribing the hardware threads, when that is the case
133 void checkHardwareOversubscription(int numThreadsOnThisRank,
134 int rank,
135 const HardwareTopology &hwTop,
136 const PhysicalNodeCommunicator &comm,
137 const MDLogger &mdlog);
139 } // namespace gmx
141 #endif