Stop duplicate printing of detected GPUs
[gromacs.git] / src / gromacs / hardware / hardwareassign.cpp
blob25d664fce4d253b90e4a68711cb909c6cfc63581
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, 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 #include "gmxpre.h"
37 #include "hardwareassign.h"
39 #include "config.h"
41 #include <cstring>
43 #include <algorithm>
44 #include <string>
46 #include "gromacs/gmxlib/network.h"
47 #include "gromacs/gpu_utils/gpu_utils.h"
48 #include "gromacs/hardware/detecthardware.h"
49 #include "gromacs/hardware/gpu_hw_info.h"
50 #include "gromacs/mdtypes/commrec.h"
51 #include "gromacs/utility/basenetwork.h"
52 #include "gromacs/utility/exceptions.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/gmxassert.h"
55 #include "gromacs/utility/logger.h"
56 #include "gromacs/utility/smalloc.h"
57 #include "gromacs/utility/stringutil.h"
58 #include "gromacs/utility/sysinfo.h"
60 #define HOSTNAMELEN 80
62 /*! \internal \brief
63 * This function is responsible for mapping the GPUs to the processes on a single node
64 * (filling the gpu_opt->dev_use array).
66 * \param[in] compatibleGpus Vector of GPUs that are compatible
67 * \param[in,out] gpu_opt Input/output GPU assignment data.
68 * \param[in] nrank Number of PP GPU ranks on the node.
69 * \param[in] rank Index of PP GPU rank on the node.
71 * This selects the GPUs we will use. This is an operation local to each physical node.
72 * If we have less MPI ranks than GPUs, we will waste some GPUs.
74 static void assign_rank_gpu_ids(const std::vector<int> &compatibleGpus,
75 gmx_gpu_opt_t *gpu_opt, int nrank, int rank)
77 int numCompatibleGpus = static_cast<int>(compatibleGpus.size());
78 GMX_RELEASE_ASSERT(gpu_opt, "Invalid gpu_opt pointer passed");
79 GMX_RELEASE_ASSERT(nrank >= 1,
80 gmx::formatString("Invalid limit (%d) for the number of GPUs (detected %d compatible GPUs)",
81 rank, numCompatibleGpus).c_str());
83 if (numCompatibleGpus == 0)
85 char host[HOSTNAMELEN];
87 gmx_gethostname(host, HOSTNAMELEN);
88 gmx_fatal(FARGS, "A GPU was requested on host %s, but no compatible GPUs were detected. All nodes with PP ranks need to have GPUs. If you intended to use GPU acceleration in a parallel run, you can either avoid using the nodes that don't have GPUs or place PME ranks on these nodes.", host);
91 int nshare;
93 nshare = 1;
94 if (nrank > numCompatibleGpus)
96 if (nrank % numCompatibleGpus == 0)
98 nshare = gmx_gpu_sharing_supported() ? nrank/numCompatibleGpus : 1;
100 else
102 if (rank == 0)
104 gmx_fatal(FARGS, "The number of MPI ranks (%d) in a physical node is not a multiple of the number of GPUs (%d). Select a different number of MPI ranks or use the -gpu_id option to manually specify the GPU to be used.",
105 nrank, numCompatibleGpus);
108 #if GMX_MPI
109 /* We use a global barrier to prevent ranks from continuing with
110 * an invalid setup.
112 MPI_Barrier(MPI_COMM_WORLD);
113 #endif
117 /* Here we will waste GPUs when nrank < numCompatibleGpus */
118 gpu_opt->n_dev_use = std::min(numCompatibleGpus*nshare, nrank);
119 if (!gmx_multiple_gpu_per_node_supported())
121 gpu_opt->n_dev_use = std::min(gpu_opt->n_dev_use, 1);
123 snew(gpu_opt->dev_use, gpu_opt->n_dev_use);
124 for (int i = 0; i != gpu_opt->n_dev_use; ++i)
126 /* TODO: improve this implementation: either sort GPUs or remove the weakest here */
127 gpu_opt->dev_use[i] = compatibleGpus[i/nshare];
131 /*! \brief Return whether all selected GPUs are compatible.
133 * Given the list of selected GPU device IDs in \c gpu_opt and
134 * detected GPUs in \c gpu_info, return whether all selected GPUs are
135 * compatible. If not, place a suitable string in \c errorMessage.
137 * \param[in] gpu_info Information about detected GPUs
138 * \param[in] gpu_opt pointer to structure holding GPU options
139 * \param[out] errorMessage pointer to string to hold a possible error message (is not updated when returning true)
140 * \returns true if every requested GPU is compatible
142 static bool checkGpuSelection(const gmx_gpu_info_t &gpu_info,
143 const gmx_gpu_opt_t *gpu_opt,
144 std::string *errorMessage)
146 bool allOK = true;
147 std::string message = "Some of the requested GPUs do not exist, behave strangely, or are not compatible:\n";
148 for (int i = 0; i < gpu_opt->n_dev_use; i++)
150 GMX_ASSERT(gpu_opt, "Invalid gpu_opt");
151 GMX_ASSERT(gpu_opt->dev_use, "Invalid gpu_opt->dev_use");
153 int id = gpu_opt->dev_use[i];
154 if (!isGpuCompatible(gpu_info, id))
156 allOK = false;
157 message += gmx::formatString(" GPU #%d: %s\n",
159 getGpuCompatibilityDescription(gpu_info, id));
162 if (!allOK && errorMessage)
164 *errorMessage = message;
166 return allOK;
169 std::vector<int> getCompatibleGpus(const gmx_gpu_info_t &gpu_info)
171 // Possible minor over-allocation here, but not important for anything
172 std::vector<int> compatibleGpus;
173 compatibleGpus.reserve(gpu_info.n_dev);
174 for (int i = 0; i < gpu_info.n_dev; i++)
176 GMX_ASSERT(gpu_info.gpu_dev, "Invalid gpu_info.gpu_dev");
177 if (isGpuCompatible(gpu_info, i))
179 compatibleGpus.push_back(i);
182 return compatibleGpus;
185 void gmx_select_rank_gpu_ids(const t_commrec *cr,
186 const gmx_gpu_info_t &gpu_info,
187 bool userSetGpuIds,
188 gmx_gpu_opt_t *gpu_opt)
190 if (!(cr->duty & DUTY_PP))
192 /* Our rank is not doing PP, we don't use a GPU */
193 return;
196 if (userSetGpuIds)
198 /* Check the GPU IDs passed by the user.
199 * (GPU IDs have been parsed by gmx_parse_gpu_ids before)
201 std::string errorMessage;
202 if (!checkGpuSelection(gpu_info, gpu_opt, &errorMessage))
204 gmx_fatal(FARGS, errorMessage.c_str());
207 else
209 auto compatibleGpus = getCompatibleGpus(gpu_info);
210 assign_rank_gpu_ids(compatibleGpus, gpu_opt, cr->nrank_pp_intranode, cr->rank_pp_intranode);