Move some mdrun code out of gmxlib/
[gromacs.git] / src / gromacs / mdrunutility / threadaffinity.cpp
blobaabe94144d873af75a7a0676b62d47263e9cee4a
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, 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 "threadaffinity.h"
39 #include "config.h"
41 #include <cerrno>
42 #include <cstdio>
43 #include <cstring>
45 #ifdef HAVE_SCHED_AFFINITY
46 # include <sched.h>
47 # include <sys/syscall.h>
48 #endif
50 #include "thread_mpi/threads.h"
52 #include "gromacs/gmxlib/gmx_omp_nthreads.h"
53 #include "gromacs/gmxlib/md_logging.h"
54 #include "gromacs/hardware/hardwaretopology.h"
55 #include "gromacs/hardware/hw_info.h"
56 #include "gromacs/mdtypes/commrec.h"
57 #include "gromacs/utility/basenetwork.h"
58 #include "gromacs/utility/cstringutil.h"
59 #include "gromacs/utility/exceptions.h"
60 #include "gromacs/utility/fatalerror.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/gmxomp.h"
63 #include "gromacs/utility/programcontext.h"
64 #include "gromacs/utility/scoped_cptr.h"
65 #include "gromacs/utility/smalloc.h"
68 static int
69 get_thread_affinity_layout(FILE *fplog,
70 const t_commrec *cr,
71 const gmx_hw_info_t * hwinfo,
72 int threads,
73 int pin_offset, int * pin_stride,
74 int **localityOrder)
76 int hwThreads;
77 int hwThreadsPerCore = 0;
78 bool bPickPinStride;
79 bool haveTopology;
81 const gmx::HardwareTopology &hwTop = *hwinfo->hardwareTopology;
83 haveTopology = (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic);
85 if (pin_offset < 0)
87 gmx_fatal(FARGS, "Negative thread pinning offset requested");
89 if (*pin_stride < 0)
91 gmx_fatal(FARGS, "Negative thread pinning stride requested");
94 if (haveTopology)
96 hwThreads = hwTop.machine().logicalProcessorCount;
97 // Just use the value for the first core
98 hwThreadsPerCore = hwTop.machine().sockets[0].cores[0].hwThreads.size();
99 snew(*localityOrder, hwThreads);
100 int i = 0;
101 for (auto &s : hwTop.machine().sockets)
103 for (auto &c : s.cores)
105 for (auto &t : c.hwThreads)
107 (*localityOrder)[i++] = t.logicalProcessorId;
112 else
114 /* topology information not available or invalid, ignore it */
115 hwThreads = hwinfo->nthreads_hw_avail;
116 *localityOrder = NULL;
118 if (hwThreads <= 0)
120 /* We don't know anything about the hardware, don't pin */
121 md_print_warn(cr, fplog,
122 "NOTE: No information on available cores, thread pinning disabled.");
124 return -1;
128 if (threads > hwThreads)
130 /* We are oversubscribing, don't pin */
131 md_print_warn(NULL, fplog,
132 "NOTE: Oversubscribing the CPU, will not pin threads");
134 return -1;
137 if (pin_offset + threads > hwThreads)
139 /* We are oversubscribing, don't pin */
140 md_print_warn(NULL, fplog,
141 "WARNING: Requested offset too large for available cores, thread pinning disabled.");
143 return -1;
147 /* do we need to choose the pinning stride? */
148 bPickPinStride = (*pin_stride == 0);
150 if (bPickPinStride)
152 if (haveTopology && pin_offset + threads*hwThreadsPerCore <= hwThreads)
154 /* Put one thread on each physical core */
155 *pin_stride = hwThreadsPerCore;
157 else
159 /* We don't know if we have SMT, and if we do, we don't know
160 * if hw threads in the same physical core are consecutive.
161 * Without SMT the pinning layout should not matter too much.
162 * so we assume a consecutive layout and maximally spread out"
163 * the threads at equal threads per core.
164 * Note that IBM is the major non-x86 case with cpuid support
165 * and probably threads are already pinned by the queuing system,
166 * so we wouldn't end up here in the first place.
168 *pin_stride = (hwThreads - pin_offset)/threads;
171 else
173 /* Check the placement of the thread with the largest index to make sure
174 * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
175 if (pin_offset + (threads-1)*(*pin_stride) >= hwThreads)
177 /* We are oversubscribing, don't pin */
178 md_print_warn(NULL, fplog,
179 "WARNING: Requested stride too large for available cores, thread pinning disabled.");
181 return -1;
185 if (fplog != NULL)
187 fprintf(fplog, "Pinning threads with a%s logical core stride of %d\n",
188 bPickPinStride ? "n auto-selected" : " user-specified",
189 *pin_stride);
192 return 0;
195 /* Set CPU affinity. Can be important for performance.
196 On some systems (e.g. Cray) CPU Affinity is set by default.
197 But default assigning doesn't work (well) with only some ranks
198 having threads. This causes very low performance.
199 External tools have cumbersome syntax for setting affinity
200 in the case that only some ranks have threads.
201 Thus it is important that GROMACS sets the affinity internally
202 if only PME is using threads.
204 void
205 gmx_set_thread_affinity(FILE *fplog,
206 const t_commrec *cr,
207 const gmx_hw_opt_t *hw_opt,
208 const gmx_hw_info_t *hwinfo)
210 int nth_affinity_set, thread0_id_node,
211 nthread_local, nthread_node;
212 int offset;
213 int * localityOrder = nullptr;
214 int rc;
216 if (hw_opt->thread_affinity == threadaffOFF)
218 /* Nothing to do */
219 return;
222 /* If the tMPI thread affinity setting is not supported encourage the user
223 * to report it as it's either a bug or an exotic platform which we might
224 * want to support. */
225 if (tMPI_Thread_setaffinity_support() != TMPI_SETAFFINITY_SUPPORT_YES)
227 /* we know Mac OS & BlueGene do not support setting thread affinity, so there's
228 no point in warning the user in that case. In any other case
229 the user might be able to do something about it. */
230 #if !defined(__APPLE__) && !defined(__bg__)
231 md_print_warn(NULL, fplog,
232 "NOTE: Cannot set thread affinities on the current platform.");
233 #endif /* __APPLE__ */
234 return;
237 /* threads on this MPI process or TMPI thread */
238 if (cr->duty & DUTY_PP)
240 nthread_local = gmx_omp_nthreads_get(emntNonbonded);
242 else
244 nthread_local = gmx_omp_nthreads_get(emntPME);
247 /* map the current process to cores */
248 thread0_id_node = 0;
249 nthread_node = nthread_local;
250 #ifdef GMX_MPI
251 if (PAR(cr) || MULTISIM(cr))
253 /* We need to determine a scan of the thread counts in this
254 * compute node.
256 MPI_Comm comm_intra;
258 MPI_Comm_split(MPI_COMM_WORLD,
259 gmx_physicalnode_id_hash(), cr->rank_intranode,
260 &comm_intra);
261 MPI_Scan(&nthread_local, &thread0_id_node, 1, MPI_INT, MPI_SUM, comm_intra);
262 /* MPI_Scan is inclusive, but here we need exclusive */
263 thread0_id_node -= nthread_local;
264 /* Get the total number of threads on this physical node */
265 MPI_Allreduce(&nthread_local, &nthread_node, 1, MPI_INT, MPI_SUM, comm_intra);
266 MPI_Comm_free(&comm_intra);
268 #endif
270 if (hw_opt->thread_affinity == threadaffAUTO &&
271 nthread_node != hwinfo->nthreads_hw_avail)
273 if (nthread_node > 1 && nthread_node < hwinfo->nthreads_hw_avail)
275 md_print_warn(cr, fplog,
276 "NOTE: The number of threads is not equal to the number of (logical) cores\n"
277 " and the -pin option is set to auto: will not pin thread to cores.\n"
278 " This can lead to significant performance degradation.\n"
279 " Consider using -pin on (and -pinoffset in case you run multiple jobs).\n");
282 return;
285 offset = 0;
286 if (hw_opt->core_pinning_offset != 0)
288 offset = hw_opt->core_pinning_offset;
289 md_print_info(cr, fplog, "Applying core pinning offset %d\n", offset);
292 int core_pinning_stride = hw_opt->core_pinning_stride;
293 rc = get_thread_affinity_layout(fplog, cr, hwinfo,
294 nthread_node,
295 offset, &core_pinning_stride,
296 &localityOrder);
297 gmx::scoped_guard_sfree localityOrderGuard(localityOrder);
299 if (rc != 0)
301 /* Incompatible layout, don't pin, warning was already issued */
302 return;
305 /* Set the per-thread affinity. In order to be able to check the success
306 * of affinity settings, we will set nth_affinity_set to 1 on threads
307 * where the affinity setting succeded and to 0 where it failed.
308 * Reducing these 0/1 values over the threads will give the total number
309 * of threads on which we succeeded.
312 // To avoid warnings from the static analyzer we initialize nth_affinity_set
313 // to zero outside the OpenMP block, and then add to it inside the block.
314 // The value will still always be 0 or 1 from each thread.
315 nth_affinity_set = 0;
316 #pragma omp parallel num_threads(nthread_local) reduction(+:nth_affinity_set)
320 int thread_id, thread_id_node;
321 int index, core;
322 gmx_bool setaffinity_ret;
324 thread_id = gmx_omp_get_thread_num();
325 thread_id_node = thread0_id_node + thread_id;
326 index = offset + thread_id_node*core_pinning_stride;
327 if (localityOrder != nullptr)
329 core = localityOrder[index];
331 else
333 core = index;
336 setaffinity_ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
338 /* store the per-thread success-values of the setaffinity */
339 nth_affinity_set += (setaffinity_ret == 0);
341 if (debug)
343 fprintf(debug, "On rank %2d, thread %2d, index %2d, core %2d the affinity setting returned %d\n",
344 cr->nodeid, gmx_omp_get_thread_num(), index, core, setaffinity_ret);
347 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
350 if (nth_affinity_set > nthread_local)
352 char msg[STRLEN];
354 sprintf(msg, "Looks like we have set affinity for more threads than "
355 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
356 gmx_incons(msg);
358 else
360 /* check & warn if some threads failed to set their affinities */
361 if (nth_affinity_set != nthread_local)
363 char sbuf1[STRLEN], sbuf2[STRLEN];
365 /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
366 sbuf1[0] = sbuf2[0] = '\0';
367 /* Only add rank info if we have more than one rank. */
368 if (cr->nnodes > 1)
370 #ifdef GMX_MPI
371 #ifdef GMX_THREAD_MPI
372 sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
373 #else /* GMX_LIB_MPI */
374 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
375 #endif
376 #endif /* GMX_MPI */
379 if (nthread_local > 1)
381 sprintf(sbuf2, "for %d/%d thread%s ",
382 nthread_local - nth_affinity_set, nthread_local,
383 nthread_local > 1 ? "s" : "");
386 md_print_warn(NULL, fplog,
387 "NOTE: %sAffinity setting %sfailed. This can cause performance degradation.\n"
388 " If you think your settings are correct, ask on the gmx-users list.",
389 sbuf1, sbuf2);
392 return;
395 /* Check the process affinity mask and if it is found to be non-zero,
396 * will honor it and disable mdrun internal affinity setting.
397 * Note that this will only work on Linux as we use a GNU feature.
399 void
400 gmx_check_thread_affinity_set(FILE *fplog,
401 const t_commrec *cr,
402 gmx_hw_opt_t *hw_opt,
403 int gmx_unused nthreads_hw_avail,
404 gmx_bool bAfterOpenmpInit)
406 GMX_RELEASE_ASSERT(hw_opt, "hw_opt must be a non-NULL pointer");
408 if (!bAfterOpenmpInit)
410 /* Check for externally set OpenMP affinity and turn off internal
411 * pinning if any is found. We need to do this check early to tell
412 * thread-MPI whether it should do pinning when spawning threads.
413 * TODO: the above no longer holds, we should move these checks later
415 if (hw_opt->thread_affinity != threadaffOFF)
417 char *message;
418 if (!gmx_omp_check_thread_affinity(&message))
420 /* TODO: with -pin auto we should only warn when using all cores */
421 md_print_warn(cr, fplog, "%s", message);
422 sfree(message);
423 hw_opt->thread_affinity = threadaffOFF;
427 /* With thread-MPI this is needed as pinning might get turned off,
428 * which needs to be known before starting thread-MPI.
429 * With thread-MPI hw_opt is processed here on the master rank
430 * and passed to the other ranks later, so we only do this on master.
432 if (!SIMMASTER(cr))
434 return;
436 #ifndef GMX_THREAD_MPI
437 return;
438 #endif
441 #ifdef HAVE_SCHED_AFFINITY
442 int ret;
443 cpu_set_t mask_current;
445 if (hw_opt->thread_affinity == threadaffOFF)
447 /* internal affinity setting is off, don't bother checking process affinity */
448 return;
451 CPU_ZERO(&mask_current);
452 if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
454 /* failed to query affinity mask, will just return */
455 if (debug)
457 fprintf(debug, "Failed to query affinity mask (error %d)", ret);
459 return;
462 /* Before proceeding with the actual check, make sure that the number of
463 * detected CPUs is >= the CPUs in the current set.
464 * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
465 #ifdef CPU_COUNT
466 if (nthreads_hw_avail < CPU_COUNT(&mask_current))
468 if (debug)
470 fprintf(debug, "%d hardware threads detected, but %d was returned by CPU_COUNT",
471 nthreads_hw_avail, CPU_COUNT(&mask_current));
473 return;
475 #endif /* CPU_COUNT */
477 gmx_bool bAllSet = TRUE;
478 for (int i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
480 bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
483 #ifdef GMX_LIB_MPI
484 gmx_bool bAllSet_All;
486 MPI_Allreduce(&bAllSet, &bAllSet_All, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
487 bAllSet = bAllSet_All;
488 #endif
490 if (!bAllSet)
492 if (hw_opt->thread_affinity == threadaffAUTO)
494 if (!bAfterOpenmpInit)
496 md_print_warn(cr, fplog,
497 "Non-default thread affinity set, disabling internal thread affinity");
499 else
501 md_print_warn(cr, fplog,
502 "Non-default thread affinity set probably by the OpenMP library,\n"
503 "disabling internal thread affinity");
505 hw_opt->thread_affinity = threadaffOFF;
507 else
509 /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
510 if (bAfterOpenmpInit)
512 md_print_warn(cr, fplog,
513 "Overriding thread affinity set outside %s\n",
514 gmx::getProgramContext().displayName());
518 if (debug)
520 fprintf(debug, "Non-default affinity mask found\n");
523 else
525 if (debug)
527 fprintf(debug, "Default affinity mask found\n");
530 #endif /* HAVE_SCHED_AFFINITY */