2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012, by the GROMACS development team, led by
5 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
6 * others, as listed in the AUTHORS file in the top-level source
7 * 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.
38 #if defined(HAVE_SCHED_H) && defined(HAVE_SCHED_GETAFFINITY)
41 #include <sys/syscall.h>
46 #include "types/commrec.h"
47 #include "types/hw_info.h"
48 #include "gmx_cpuid.h"
50 #include "gmx_omp_nthreads.h"
52 #include "md_logging.h"
54 #include "gmx_thread_affinity.h"
56 #include "thread_mpi/threads.h"
60 get_thread_affinity_layout(FILE *fplog,
62 const gmx_hw_info_t * hwinfo,
64 int pin_offset, int * pin_stride,
65 const int **locality_order)
67 int nhwthreads, npkg, ncores, nhwthreads_per_core, rc;
70 const int * hwthread_id;
74 gmx_fatal(FARGS, "Negative thread pinning offset requested");
78 gmx_fatal(FARGS, "Negative thread pinning stride requested");
81 rc = gmx_cpuid_topology(hwinfo->cpuid_info, &nhwthreads, &npkg, &ncores,
83 &pkg_id, &core_id, &hwthread_id, locality_order);
87 nhwthreads = hwinfo->nthreads_hw_avail;
88 *locality_order = NULL;
92 /* We don't know anything about the hardware, don't pin */
93 md_print_warn(cr, fplog,
94 "We don't know how many logical cores we have, will not pin threads");
100 if (pin_offset + nthreads > nhwthreads)
102 /* We are oversubscribing, don't pin */
103 md_print_warn(NULL, fplog,
104 "More threads requested than available logical cores, will not pin threads");
109 /* Check if we need to choose the pinning stride */
110 if (*pin_stride == 0)
112 if (rc == 0 && pin_offset + nthreads*nhwthreads_per_core <= nhwthreads)
114 /* Put one thread on each physical core */
115 *pin_stride = nhwthreads_per_core;
119 /* We don't know if we have SMT, and if we do, we don't know
120 * if hw threads in the same physical core are consecutive.
121 * Without SMT the pinning layout should not matter too much.
122 * so we assume a consecutive layout and maximally spread out"
123 * the threads at equal threads per core.
124 * Note that IBM is the major non-x86 case with cpuid support
125 * and probably threads are already pinned by the queuing system,
126 * so we wouldn't end up here in the first place.
128 *pin_stride = (nhwthreads - pin_offset)/nthreads;
133 fprintf(fplog, "Pinning threads with a logical core stride of %d\n",
139 if (pin_offset + nthreads*(*pin_stride) > nhwthreads)
141 /* We are oversubscribing, don't pin */
142 md_print_warn(NULL, fplog,
143 "The requested pinning stride is too large for the available logical cores, will not pin threads");
152 /* Set CPU affinity. Can be important for performance.
153 On some systems (e.g. Cray) CPU Affinity is set by default.
154 But default assigning doesn't work (well) with only some ranks
155 having threads. This causes very low performance.
156 External tools have cumbersome syntax for setting affinity
157 in the case that only some ranks have threads.
158 Thus it is important that GROMACS sets the affinity internally
159 if only PME is using threads.
162 gmx_set_thread_affinity(FILE *fplog,
164 gmx_hw_opt_t *hw_opt,
166 const gmx_hw_info_t *hwinfo,
167 const t_inputrec *inputrec)
169 int nth_affinity_set, thread_id_node, thread_id,
170 nthread_local, nthread_node, nthread_hw_max, nphyscore;
172 const int *locality_order;
175 if (hw_opt->thread_affinity == threadaffOFF)
182 /* If the tMPI thread affinity setting is not supported encourage the user
183 * to report it as it's either a bug or an exotic platform which we might
184 * want to support. */
185 if (tMPI_Thread_setaffinity_support() != TMPI_SETAFFINITY_SUPPORT_YES)
187 md_print_warn(NULL, fplog,
188 "Can not set thread affinities on the current platform. On NUMA systems this\n"
189 "can cause performance degradation. If you think your platform should support\n"
190 "setting affinities, contact the GROMACS developers.");
193 #endif /* __APPLE__ */
195 /* threads on this MPI process or TMPI thread */
196 if (cr->duty & DUTY_PP)
198 nthread_local = gmx_omp_nthreads_get(emntNonbonded);
202 nthread_local = gmx_omp_nthreads_get(emntPME);
205 /* map the current process to cores */
207 nthread_node = nthread_local;
209 if (PAR(cr) || MULTISIM(cr))
211 /* We need to determine a scan of the thread counts in this
216 MPI_Comm_split(MPI_COMM_WORLD, gmx_hostname_num(), cr->rank_intranode,
218 MPI_Scan(&nthread_local, &thread_id_node, 1, MPI_INT, MPI_SUM, comm_intra);
219 /* MPI_Scan is inclusive, but here we need exclusive */
220 thread_id_node -= nthread_local;
221 /* Get the total number of threads on this physical node */
222 MPI_Allreduce(&nthread_local, &nthread_node, 1, MPI_INT, MPI_SUM, comm_intra);
223 MPI_Comm_free(&comm_intra);
227 if (hw_opt->thread_affinity == threadaffAUTO &&
228 nthread_node != hwinfo->nthreads_hw_avail)
230 if (nthread_node > 1 && nthread_node < hwinfo->nthreads_hw_avail)
232 md_print_warn(cr, fplog,
233 "NOTE: The number of threads is not equal to the number of (logical) cores\n"
234 " and the -pin option is set to auto: will not pin thread to cores.\n"
235 " This can lead to significant performance degradation.\n"
236 " Consider using -pin on (and -pinoffset in case you run multiple jobs).\n");
243 if (hw_opt->core_pinning_offset != 0)
245 offset = hw_opt->core_pinning_offset;
246 md_print_info(cr, fplog, "Applying core pinning offset %d\n", offset);
249 rc = get_thread_affinity_layout(fplog, cr, hwinfo,
251 offset, &hw_opt->core_pinning_stride,
255 /* Incompatible layout, don't pin, warning was already issued */
259 /* Set the per-thread affinity. In order to be able to check the success
260 * of affinity settings, we will set nth_affinity_set to 1 on threads
261 * where the affinity setting succeded and to 0 where it failed.
262 * Reducing these 0/1 values over the threads will give the total number
263 * of threads on which we succeeded.
265 nth_affinity_set = 0;
266 #pragma omp parallel firstprivate(thread_id_node) num_threads(nthread_local) \
267 reduction(+:nth_affinity_set)
270 gmx_bool setaffinity_ret;
272 thread_id = gmx_omp_get_thread_num();
273 thread_id_node += thread_id;
274 index = offset + thread_id_node*hw_opt->core_pinning_stride;
275 if (locality_order != NULL)
277 core = locality_order[index];
284 setaffinity_ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
286 /* store the per-thread success-values of the setaffinity */
287 nth_affinity_set = (setaffinity_ret == 0);
291 fprintf(debug, "On rank %2d, thread %2d, core %2d the affinity setting returned %d\n",
292 cr->nodeid, gmx_omp_get_thread_num(), core, setaffinity_ret);
296 if (nth_affinity_set > nthread_local)
300 sprintf(msg, "Looks like we have set affinity for more threads than "
301 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
306 /* check & warn if some threads failed to set their affinities */
307 if (nth_affinity_set != nthread_local)
309 char sbuf1[STRLEN], sbuf2[STRLEN];
311 /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
312 sbuf1[0] = sbuf2[0] = '\0';
314 #ifdef GMX_THREAD_MPI
315 sprintf(sbuf1, "In thread-MPI thread #%d: ", cr->nodeid);
316 #else /* GMX_LIB_MPI */
317 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
321 if (nthread_local > 1)
323 sprintf(sbuf2, "of %d/%d thread%s ",
324 nthread_local - nth_affinity_set, nthread_local,
325 (nthread_local - nth_affinity_set) > 1 ? "s" : "");
328 md_print_warn(NULL, fplog,
329 "NOTE: %sAffinity setting %sfailed.\n"
330 " This can cause performance degradation!",
336 /* Check the process affinity mask and if it is found to be non-zero,
337 * will honor it and disable mdrun internal affinity setting.
338 * Note that this will only work on Linux as we use a GNU feature.
341 gmx_check_thread_affinity_set(FILE *fplog, const t_commrec *cr,
342 gmx_hw_opt_t *hw_opt, int ncpus,
343 gmx_bool bAfterOpenmpInit)
345 #ifdef HAVE_SCHED_GETAFFINITY
346 cpu_set_t mask_current;
347 int i, ret, cpu_count, cpu_set;
351 if (hw_opt->thread_affinity == threadaffOFF)
353 /* internal affinity setting is off, don't bother checking process affinity */
357 CPU_ZERO(&mask_current);
358 if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
360 /* failed to query affinity mask, will just return */
363 fprintf(debug, "Failed to query affinity mask (error %d)", ret);
368 /* Before proceeding with the actual check, make sure that the number of
369 * detected CPUs is >= the CPUs in the current set.
370 * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
372 if (ncpus < CPU_COUNT(&mask_current))
376 fprintf(debug, "%d CPUs detected, but %d was returned by CPU_COUNT",
377 ncpus, CPU_COUNT(&mask_current));
381 #endif /* CPU_COUNT */
384 for (i = 0; (i < ncpus && i < CPU_SETSIZE); i++)
386 bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
391 if (hw_opt->thread_affinity == threadaffAUTO)
393 if (!bAfterOpenmpInit)
395 md_print_warn(cr, fplog,
396 "Non-default thread affinity set, disabling internal thread affinity");
400 md_print_warn(cr, fplog,
401 "Non-default thread affinity set probably by the OpenMP library,\n"
402 "disabling internal thread affinity");
404 hw_opt->thread_affinity = threadaffOFF;
408 /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
409 if (bAfterOpenmpInit)
411 md_print_warn(cr, fplog,
412 "Overriding thread affinity set outside %s\n",
419 fprintf(debug, "Non-default affinity mask found\n");
426 fprintf(debug, "Default affinity mask found\n");
429 #endif /* HAVE_SCHED_GETAFFINITY */