Merge release-5-0 into master
[gromacs/AngularHB.git] / src / gromacs / gmxlib / gmx_thread_affinity.c
blob8814ec410f72b9c3bb89b871ba7bb76d995063a5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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 #ifdef HAVE_CONFIG_H
36 #include <config.h>
37 #endif
38 #if defined(HAVE_SCHED_H) && defined(HAVE_SCHED_GETAFFINITY)
39 #define _GNU_SOURCE
40 #include <sched.h>
41 #include <sys/syscall.h>
42 #endif
43 #include <string.h>
44 #include <errno.h>
45 #include <assert.h>
46 #include <stdio.h>
48 #include "thread_mpi/threads.h"
50 #include "typedefs.h"
51 #include "types/commrec.h"
52 #include "types/hw_info.h"
53 #include "copyrite.h"
54 #include "gmx_cpuid.h"
55 #include "gmx_omp_nthreads.h"
56 #include "md_logging.h"
57 #include "gmx_thread_affinity.h"
59 #include "gromacs/utility/basenetwork.h"
60 #include "gromacs/utility/cstringutil.h"
61 #include "gromacs/utility/fatalerror.h"
62 #include "gromacs/utility/gmxomp.h"
63 #include "gromacs/utility/smalloc.h"
65 static int
66 get_thread_affinity_layout(FILE *fplog,
67 const t_commrec *cr,
68 const gmx_hw_info_t * hwinfo,
69 int nthreads,
70 int pin_offset, int * pin_stride,
71 const int **locality_order)
73 int nhwthreads, npkg, ncores, nhwthreads_per_core, rc;
74 const int * pkg_id;
75 const int * core_id;
76 const int * hwthread_id;
77 gmx_bool bPickPinStride;
79 if (pin_offset < 0)
81 gmx_fatal(FARGS, "Negative thread pinning offset requested");
83 if (*pin_stride < 0)
85 gmx_fatal(FARGS, "Negative thread pinning stride requested");
88 rc = gmx_cpuid_topology(hwinfo->cpuid_info, &nhwthreads, &npkg, &ncores,
89 &nhwthreads_per_core,
90 &pkg_id, &core_id, &hwthread_id, locality_order);
92 if (rc != 0)
94 /* topology information not available or invalid, ignore it */
95 nhwthreads = hwinfo->nthreads_hw_avail;
96 *locality_order = NULL;
98 if (nhwthreads <= 0)
100 /* We don't know anything about the hardware, don't pin */
101 md_print_warn(cr, fplog,
102 "NOTE: We don't know how many logical cores we have, will not pin threads");
104 return -1;
108 if (nthreads > nhwthreads)
110 /* We are oversubscribing, don't pin */
111 md_print_warn(NULL, fplog,
112 "WARNING: Oversubscribing the CPU, will not pin threads");
114 return -1;
117 if (pin_offset + nthreads > nhwthreads)
119 /* We are oversubscribing, don't pin */
120 md_print_warn(NULL, fplog,
121 "WARNING: The requested pin offset is too large for the available logical cores,\n"
122 " will not pin threads");
124 return -1;
128 /* do we need to choose the pinning stride? */
129 bPickPinStride = (*pin_stride == 0);
131 if (bPickPinStride)
133 if (rc == 0 && pin_offset + nthreads*nhwthreads_per_core <= nhwthreads)
135 /* Put one thread on each physical core */
136 *pin_stride = nhwthreads_per_core;
138 else
140 /* We don't know if we have SMT, and if we do, we don't know
141 * if hw threads in the same physical core are consecutive.
142 * Without SMT the pinning layout should not matter too much.
143 * so we assume a consecutive layout and maximally spread out"
144 * the threads at equal threads per core.
145 * Note that IBM is the major non-x86 case with cpuid support
146 * and probably threads are already pinned by the queuing system,
147 * so we wouldn't end up here in the first place.
149 *pin_stride = (nhwthreads - pin_offset)/nthreads;
152 else
154 /* Check the placement of the thread with the largest index to make sure
155 * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
156 if (pin_offset + (nthreads-1)*(*pin_stride) >= nhwthreads)
158 /* We are oversubscribing, don't pin */
159 md_print_warn(NULL, fplog,
160 "WARNING: The requested pinning stride is too large for the available logical cores,\n"
161 " will not pin threads");
163 return -1;
167 if (fplog != NULL)
169 fprintf(fplog, "Pinning threads with a%s logical core stride of %d\n",
170 bPickPinStride ? "n auto-selected" : " user-specified",
171 *pin_stride);
174 return 0;
177 /* Set CPU affinity. Can be important for performance.
178 On some systems (e.g. Cray) CPU Affinity is set by default.
179 But default assigning doesn't work (well) with only some ranks
180 having threads. This causes very low performance.
181 External tools have cumbersome syntax for setting affinity
182 in the case that only some ranks have threads.
183 Thus it is important that GROMACS sets the affinity internally
184 if only PME is using threads.
186 void
187 gmx_set_thread_affinity(FILE *fplog,
188 const t_commrec *cr,
189 gmx_hw_opt_t *hw_opt,
190 const gmx_hw_info_t *hwinfo)
192 int nth_affinity_set, thread0_id_node,
193 nthread_local, nthread_node, nthread_hw_max, nphyscore;
194 int offset;
195 const int *locality_order;
196 int rc;
198 if (hw_opt->thread_affinity == threadaffOFF)
200 /* Nothing to do */
201 return;
204 /* If the tMPI thread affinity setting is not supported encourage the user
205 * to report it as it's either a bug or an exotic platform which we might
206 * want to support. */
207 if (tMPI_Thread_setaffinity_support() != TMPI_SETAFFINITY_SUPPORT_YES)
209 /* we know Mac OS doesn't support setting thread affinity, so there's
210 no point in warning the user in that case. In any other case
211 the user might be able to do something about it. */
212 #ifndef __APPLE__
213 md_print_warn(NULL, fplog,
214 "Can not set thread affinities on the current platform. On NUMA systems this\n"
215 "can cause performance degradation. If you think your platform should support\n"
216 "setting affinities, contact the GROMACS developers.");
217 #endif /* __APPLE__ */
218 return;
221 /* threads on this MPI process or TMPI thread */
222 if (cr->duty & DUTY_PP)
224 nthread_local = gmx_omp_nthreads_get(emntNonbonded);
226 else
228 nthread_local = gmx_omp_nthreads_get(emntPME);
231 /* map the current process to cores */
232 thread0_id_node = 0;
233 nthread_node = nthread_local;
234 #ifdef GMX_MPI
235 if (PAR(cr) || MULTISIM(cr))
237 /* We need to determine a scan of the thread counts in this
238 * compute node.
240 MPI_Comm comm_intra;
242 MPI_Comm_split(MPI_COMM_WORLD,
243 gmx_physicalnode_id_hash(), cr->rank_intranode,
244 &comm_intra);
245 MPI_Scan(&nthread_local, &thread0_id_node, 1, MPI_INT, MPI_SUM, comm_intra);
246 /* MPI_Scan is inclusive, but here we need exclusive */
247 thread0_id_node -= nthread_local;
248 /* Get the total number of threads on this physical node */
249 MPI_Allreduce(&nthread_local, &nthread_node, 1, MPI_INT, MPI_SUM, comm_intra);
250 MPI_Comm_free(&comm_intra);
252 #endif
254 if (hw_opt->thread_affinity == threadaffAUTO &&
255 nthread_node != hwinfo->nthreads_hw_avail)
257 if (nthread_node > 1 && nthread_node < hwinfo->nthreads_hw_avail)
259 md_print_warn(cr, fplog,
260 "NOTE: The number of threads is not equal to the number of (logical) cores\n"
261 " and the -pin option is set to auto: will not pin thread to cores.\n"
262 " This can lead to significant performance degradation.\n"
263 " Consider using -pin on (and -pinoffset in case you run multiple jobs).\n");
266 return;
269 offset = 0;
270 if (hw_opt->core_pinning_offset != 0)
272 offset = hw_opt->core_pinning_offset;
273 md_print_info(cr, fplog, "Applying core pinning offset %d\n", offset);
276 rc = get_thread_affinity_layout(fplog, cr, hwinfo,
277 nthread_node,
278 offset, &hw_opt->core_pinning_stride,
279 &locality_order);
281 if (rc != 0)
283 /* Incompatible layout, don't pin, warning was already issued */
284 return;
287 /* Set the per-thread affinity. In order to be able to check the success
288 * of affinity settings, we will set nth_affinity_set to 1 on threads
289 * where the affinity setting succeded and to 0 where it failed.
290 * Reducing these 0/1 values over the threads will give the total number
291 * of threads on which we succeeded.
293 nth_affinity_set = 0;
294 #pragma omp parallel num_threads(nthread_local) reduction(+:nth_affinity_set)
296 int thread_id, thread_id_node;
297 int index, core;
298 gmx_bool setaffinity_ret;
300 thread_id = gmx_omp_get_thread_num();
301 thread_id_node = thread0_id_node + thread_id;
302 index = offset + thread_id_node*hw_opt->core_pinning_stride;
303 if (locality_order != NULL)
305 core = locality_order[index];
307 else
309 core = index;
312 setaffinity_ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
314 /* store the per-thread success-values of the setaffinity */
315 nth_affinity_set = (setaffinity_ret == 0);
317 if (debug)
319 fprintf(debug, "On rank %2d, thread %2d, index %2d, core %2d the affinity setting returned %d\n",
320 cr->nodeid, gmx_omp_get_thread_num(), index, core, setaffinity_ret);
324 if (nth_affinity_set > nthread_local)
326 char msg[STRLEN];
328 sprintf(msg, "Looks like we have set affinity for more threads than "
329 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
330 gmx_incons(msg);
332 else
334 /* check & warn if some threads failed to set their affinities */
335 if (nth_affinity_set != nthread_local)
337 char sbuf1[STRLEN], sbuf2[STRLEN];
339 /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
340 sbuf1[0] = sbuf2[0] = '\0';
341 /* Only add rank info if we have more than one rank. */
342 if (cr->nnodes > 1)
344 #ifdef GMX_MPI
345 #ifdef GMX_THREAD_MPI
346 sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
347 #else /* GMX_LIB_MPI */
348 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
349 #endif
350 #endif /* GMX_MPI */
353 if (nthread_local > 1)
355 sprintf(sbuf2, "for %d/%d thread%s ",
356 nthread_local - nth_affinity_set, nthread_local,
357 nthread_local > 1 ? "s" : "");
360 md_print_warn(NULL, fplog,
361 "WARNING: %sAffinity setting %sfailed.\n"
362 " This can cause performance degradation! If you think your setting are\n"
363 " correct, contact the GROMACS developers.",
364 sbuf1, sbuf2);
367 return;
370 /* Check the process affinity mask and if it is found to be non-zero,
371 * will honor it and disable mdrun internal affinity setting.
372 * Note that this will only work on Linux as we use a GNU feature.
374 void
375 gmx_check_thread_affinity_set(FILE *fplog,
376 const t_commrec *cr,
377 gmx_hw_opt_t *hw_opt,
378 int gmx_unused ncpus,
379 gmx_bool bAfterOpenmpInit)
381 #ifdef HAVE_SCHED_GETAFFINITY
382 cpu_set_t mask_current;
383 int i, ret, cpu_count, cpu_set;
384 gmx_bool bAllSet;
385 #endif
387 assert(hw_opt);
388 if (!bAfterOpenmpInit)
390 /* Check for externally set OpenMP affinity and turn off internal
391 * pinning if any is found. We need to do this check early to tell
392 * thread-MPI whether it should do pinning when spawning threads.
393 * TODO: the above no longer holds, we should move these checks later
395 if (hw_opt->thread_affinity != threadaffOFF)
397 char *message;
398 if (!gmx_omp_check_thread_affinity(&message))
400 /* TODO: with -pin auto we should only warn when using all cores */
401 md_print_warn(cr, fplog, "%s", message);
402 sfree(message);
403 hw_opt->thread_affinity = threadaffOFF;
407 /* With thread-MPI this is needed as pinning might get turned off,
408 * which needs to be known before starting thread-MPI.
409 * With thread-MPI hw_opt is processed here on the master rank
410 * and passed to the other ranks later, so we only do this on master.
412 if (!SIMMASTER(cr))
414 return;
416 #ifndef GMX_THREAD_MPI
417 return;
418 #endif
421 #ifdef HAVE_SCHED_GETAFFINITY
422 if (hw_opt->thread_affinity == threadaffOFF)
424 /* internal affinity setting is off, don't bother checking process affinity */
425 return;
428 CPU_ZERO(&mask_current);
429 if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
431 /* failed to query affinity mask, will just return */
432 if (debug)
434 fprintf(debug, "Failed to query affinity mask (error %d)", ret);
436 return;
439 /* Before proceeding with the actual check, make sure that the number of
440 * detected CPUs is >= the CPUs in the current set.
441 * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
442 #ifdef CPU_COUNT
443 if (ncpus < CPU_COUNT(&mask_current))
445 if (debug)
447 fprintf(debug, "%d CPUs detected, but %d was returned by CPU_COUNT",
448 ncpus, CPU_COUNT(&mask_current));
450 return;
452 #endif /* CPU_COUNT */
454 bAllSet = TRUE;
455 for (i = 0; (i < ncpus && i < CPU_SETSIZE); i++)
457 bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
460 if (!bAllSet)
462 if (hw_opt->thread_affinity == threadaffAUTO)
464 if (!bAfterOpenmpInit)
466 md_print_warn(cr, fplog,
467 "Non-default thread affinity set, disabling internal thread affinity");
469 else
471 md_print_warn(cr, fplog,
472 "Non-default thread affinity set probably by the OpenMP library,\n"
473 "disabling internal thread affinity");
475 hw_opt->thread_affinity = threadaffOFF;
477 else
479 /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
480 if (bAfterOpenmpInit)
482 md_print_warn(cr, fplog,
483 "Overriding thread affinity set outside %s\n",
484 ShortProgram());
488 if (debug)
490 fprintf(debug, "Non-default affinity mask found\n");
493 else
495 if (debug)
497 fprintf(debug, "Default affinity mask found\n");
500 #endif /* HAVE_SCHED_GETAFFINITY */