2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5 * Copyright (c) 2017,2018,2019,2020, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
38 #include "threadaffinity.h"
46 #if HAVE_SCHED_AFFINITY
50 #include "thread_mpi/threads.h"
52 #include "gromacs/hardware/hardwaretopology.h"
53 #include "gromacs/hardware/hw_info.h"
54 #include "gromacs/mdtypes/commrec.h"
55 #include "gromacs/utility/basenetwork.h"
56 #include "gromacs/utility/cstringutil.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/fatalerror.h"
59 #include "gromacs/utility/gmxassert.h"
60 #include "gromacs/utility/gmxomp.h"
61 #include "gromacs/utility/logger.h"
62 #include "gromacs/utility/physicalnodecommunicator.h"
63 #include "gromacs/utility/programcontext.h"
64 #include "gromacs/utility/smalloc.h"
65 #include "gromacs/utility/unique_cptr.h"
70 class DefaultThreadAffinityAccess : public gmx::IThreadAffinityAccess
73 bool isThreadAffinitySupported() const override
75 return tMPI_Thread_setaffinity_support() == TMPI_SETAFFINITY_SUPPORT_YES;
77 bool setCurrentThreadAffinityToCore(int core) override
79 const int ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
84 //! Global instance of DefaultThreadAffinityAccess
85 DefaultThreadAffinityAccess g_defaultAffinityAccess;
89 gmx::IThreadAffinityAccess::~IThreadAffinityAccess() {}
91 static bool invalidWithinSimulation(const t_commrec* cr, bool invalidLocally)
96 int value = invalidLocally ? 1 : 0;
98 MPI_Reduce(&value, &globalValue, 1, MPI_INT, MPI_LOR, MASTERRANK(cr), cr->mpi_comm_mysim);
99 return SIMMASTER(cr) ? (globalValue != 0) : invalidLocally;
102 GMX_UNUSED_VALUE(cr);
104 return invalidLocally;
107 static bool get_thread_affinity_layout(const gmx::MDLogger& mdlog,
109 const gmx::HardwareTopology& hwTop,
111 bool affinityIsAutoAndNumThreadsIsNotAuto,
118 int hwThreadsPerCore = 0;
123 haveTopology = (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic);
127 gmx_fatal(FARGS, "Negative thread pinning offset requested");
131 gmx_fatal(FARGS, "Negative thread pinning stride requested");
136 hwThreads = hwTop.machine().logicalProcessorCount;
137 // Just use the value for the first core
138 hwThreadsPerCore = hwTop.machine().sockets[0].cores[0].hwThreads.size();
139 snew(*localityOrder, hwThreads);
141 for (auto& s : hwTop.machine().sockets)
143 for (auto& c : s.cores)
145 for (auto& t : c.hwThreads)
147 (*localityOrder)[i++] = t.logicalProcessorId;
154 /* topology information not available or invalid, ignore it */
155 hwThreads = hwTop.machine().logicalProcessorCount;
156 *localityOrder = nullptr;
158 // Only warn about the first problem per node. Otherwise, the first test
159 // failing would essentially always cause also the other problems get
160 // reported, leading to bogus warnings. The order in the conditionals
161 // with this variable is important, since the MPI_Reduce() in
162 // invalidWithinSimulation() needs to always happen.
163 bool alreadyWarned = false;
164 invalidValue = (hwThreads <= 0);
165 if (invalidWithinSimulation(cr, invalidValue))
167 /* We don't know anything about the hardware, don't pin */
168 GMX_LOG(mdlog.warning)
170 .appendText("NOTE: No information on available cores, thread pinning disabled.");
171 alreadyWarned = true;
173 bool validLayout = !invalidValue;
175 if (affinityIsAutoAndNumThreadsIsNotAuto)
177 invalidValue = (threads != hwThreads);
178 bool warn = (invalidValue && threads > 1 && threads < hwThreads);
179 if (invalidWithinSimulation(cr, warn) && !alreadyWarned)
181 GMX_LOG(mdlog.warning)
184 "NOTE: The number of threads is not equal to the number of (logical) "
186 " and the -pin option is set to auto: will not pin threads to "
188 " This can lead to significant performance degradation.\n"
189 " Consider using -pin on (and -pinoffset in case you run multiple "
191 alreadyWarned = true;
193 validLayout = validLayout && !invalidValue;
196 invalidValue = (threads > hwThreads);
197 if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
199 GMX_LOG(mdlog.warning)
201 .appendText("NOTE: Oversubscribing the CPU, will not pin threads");
202 alreadyWarned = true;
204 validLayout = validLayout && !invalidValue;
206 invalidValue = (pin_offset + threads > hwThreads);
207 if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
209 GMX_LOG(mdlog.warning)
212 "WARNING: Requested offset too large for available cores, thread pinning "
214 alreadyWarned = true;
216 validLayout = validLayout && !invalidValue;
218 invalidValue = false;
219 /* do we need to choose the pinning stride? */
220 bPickPinStride = (*pin_stride == 0);
224 if (haveTopology && pin_offset + threads * hwThreadsPerCore <= hwThreads)
226 /* Put one thread on each physical core */
227 *pin_stride = hwThreadsPerCore;
231 /* We don't know if we have SMT, and if we do, we don't know
232 * if hw threads in the same physical core are consecutive.
233 * Without SMT the pinning layout should not matter too much.
234 * so we assume a consecutive layout and maximally spread out"
235 * the threads at equal threads per core.
236 * Note that IBM is the major non-x86 case with cpuid support
237 * and probably threads are already pinned by the queuing system,
238 * so we wouldn't end up here in the first place.
240 *pin_stride = (hwThreads - pin_offset) / threads;
245 /* Check the placement of the thread with the largest index to make sure
246 * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
247 invalidValue = (pin_offset + (threads - 1) * (*pin_stride) >= hwThreads);
249 if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
251 /* We are oversubscribing, don't pin */
252 GMX_LOG(mdlog.warning)
255 "WARNING: Requested stride too large for available cores, thread pinning "
257 alreadyWarned = true;
259 validLayout = validLayout && !invalidValue;
264 .appendTextFormatted("Pinning threads with a%s logical core stride of %d",
265 bPickPinStride ? "n auto-selected" : " user-specified",
269 *issuedWarning = alreadyWarned;
274 static bool set_affinity(const t_commrec* cr,
276 int intraNodeThreadOffset,
278 int core_pinning_stride,
279 const int* localityOrder,
280 gmx::IThreadAffinityAccess* affinityAccess)
282 // Set the per-thread affinity. In order to be able to check the success
283 // of affinity settings, we will set nth_affinity_set to 1 on threads
284 // where the affinity setting succeded and to 0 where it failed.
285 // Reducing these 0/1 values over the threads will give the total number
286 // of threads on which we succeeded.
288 // To avoid warnings from the static analyzer we initialize nth_affinity_set
289 // to zero outside the OpenMP block, and then add to it inside the block.
290 // The value will still always be 0 or 1 from each thread.
291 int nth_affinity_set = 0;
292 #pragma omp parallel num_threads(nthread_local) reduction(+ : nth_affinity_set)
296 int thread_id, thread_id_node;
299 thread_id = gmx_omp_get_thread_num();
300 thread_id_node = intraNodeThreadOffset + thread_id;
301 index = offset + thread_id_node * core_pinning_stride;
302 if (localityOrder != nullptr)
304 core = localityOrder[index];
311 const bool ret = affinityAccess->setCurrentThreadAffinityToCore(core);
313 /* store the per-thread success-values of the setaffinity */
314 nth_affinity_set += (ret ? 1 : 0);
319 "On rank %2d, thread %2d, index %2d, core %2d the affinity setting "
322 gmx_omp_get_thread_num(),
328 GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
331 if (nth_affinity_set > nthread_local)
336 "Looks like we have set affinity for more threads than "
337 "we have (%d > %d)!\n",
343 /* check & warn if some threads failed to set their affinities */
344 const bool allAffinitiesSet = (nth_affinity_set == nthread_local);
345 if (!allAffinitiesSet)
347 char sbuf1[STRLEN], sbuf2[STRLEN];
349 /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
350 sbuf1[0] = sbuf2[0] = '\0';
351 /* Only add rank info if we have more than one rank. */
356 sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
357 # else /* GMX_LIB_MPI */
358 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
363 if (nthread_local > 1)
366 "for %d/%d thread%s ",
367 nthread_local - nth_affinity_set,
369 nthread_local > 1 ? "s" : "");
372 // TODO: This output should also go through mdlog.
373 fprintf(stderr, "NOTE: %sAffinity setting %sfailed.\n", sbuf1, sbuf2);
375 return allAffinitiesSet;
378 void analyzeThreadsOnThisNode(const gmx::PhysicalNodeCommunicator& physicalNodeComm,
379 int numThreadsOnThisRank,
380 int* numThreadsOnThisNode,
381 int* intraNodeThreadOffset)
383 *intraNodeThreadOffset = 0;
384 *numThreadsOnThisNode = numThreadsOnThisRank;
386 if (physicalNodeComm.size_ > 1)
388 /* We need to determine a scan of the thread counts in this
390 MPI_Scan(&numThreadsOnThisRank, intraNodeThreadOffset, 1, MPI_INT, MPI_SUM, physicalNodeComm.comm_);
391 /* MPI_Scan is inclusive, but here we need exclusive */
392 *intraNodeThreadOffset -= numThreadsOnThisRank;
393 /* Get the total number of threads on this physical node */
395 &numThreadsOnThisRank, numThreadsOnThisNode, 1, MPI_INT, MPI_SUM, physicalNodeComm.comm_);
398 GMX_UNUSED_VALUE(physicalNodeComm);
402 /* Set CPU affinity. Can be important for performance.
403 On some systems (e.g. Cray) CPU Affinity is set by default.
404 But default assigning doesn't work (well) with only some ranks
405 having threads. This causes very low performance.
406 External tools have cumbersome syntax for setting affinity
407 in the case that only some ranks have threads.
408 Thus it is important that GROMACS sets the affinity internally
409 if only PME is using threads.
411 void gmx_set_thread_affinity(const gmx::MDLogger& mdlog,
413 const gmx_hw_opt_t* hw_opt,
414 const gmx::HardwareTopology& hwTop,
415 int numThreadsOnThisRank,
416 int numThreadsOnThisNode,
417 int intraNodeThreadOffset,
418 gmx::IThreadAffinityAccess* affinityAccess)
420 int* localityOrder = nullptr;
422 if (hw_opt->threadAffinity == ThreadAffinity::Off)
428 if (affinityAccess == nullptr)
430 affinityAccess = &g_defaultAffinityAccess;
433 /* If the tMPI thread affinity setting is not supported encourage the user
434 * to report it as it's either a bug or an exotic platform which we might
435 * want to support. */
436 if (!affinityAccess->isThreadAffinitySupported())
438 /* we know Mac OS does not support setting thread affinity, so there's
439 no point in warning the user in that case. In any other case
440 the user might be able to do something about it. */
441 #if !defined(__APPLE__)
442 GMX_LOG(mdlog.warning)
444 .appendText("NOTE: Cannot set thread affinities on the current platform.");
445 #endif /* __APPLE__ */
449 int offset = hw_opt->core_pinning_offset;
450 int core_pinning_stride = hw_opt->core_pinning_stride;
453 GMX_LOG(mdlog.warning).appendTextFormatted("Applying core pinning offset %d", offset);
456 bool affinityIsAutoAndNumThreadsIsNotAuto =
457 (hw_opt->threadAffinity == ThreadAffinity::Auto && !hw_opt->totNumThreadsIsAuto);
459 bool validLayout = get_thread_affinity_layout(mdlog,
462 numThreadsOnThisNode,
463 affinityIsAutoAndNumThreadsIsNotAuto,
465 &core_pinning_stride,
468 const gmx::sfree_guard localityOrderGuard(localityOrder);
470 bool allAffinitiesSet;
473 allAffinitiesSet = set_affinity(
474 cr, numThreadsOnThisRank, intraNodeThreadOffset, offset, core_pinning_stride, localityOrder, affinityAccess);
478 // Produce the warning if any rank fails.
479 allAffinitiesSet = false;
481 if (invalidWithinSimulation(cr, !allAffinitiesSet) && !issuedWarning)
483 GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: Thread affinity was not set.");
487 /* Detects and returns whether we have the default affinity mask
489 * Returns true when we can query thread affinities and CPU count is
490 * consistent and we have default affinity mask on all ranks.
492 * Should be called simultaneously by all MPI ranks.
494 static bool detectDefaultAffinityMask(const int nthreads_hw_avail)
496 bool detectedDefaultAffinityMask = true;
498 #if HAVE_SCHED_AFFINITY
499 cpu_set_t mask_current;
500 CPU_ZERO(&mask_current);
502 if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
504 /* failed to query affinity mask, will just return */
507 fprintf(debug, "Failed to query affinity mask (error %d)", ret);
509 detectedDefaultAffinityMask = false;
512 /* Before proceeding with the actual check, make sure that the number of
513 * detected CPUs is >= the CPUs in the current set.
514 * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
516 if (detectedDefaultAffinityMask && nthreads_hw_avail < CPU_COUNT(&mask_current))
521 "%d hardware threads detected, but %d was returned by CPU_COUNT",
523 CPU_COUNT(&mask_current));
525 detectedDefaultAffinityMask = false;
527 # endif /* CPU_COUNT */
529 if (detectedDefaultAffinityMask)
531 /* Here we check whether all bits of the affinity mask are set.
532 * Note that this mask can change while the program is executing,
533 * e.g. by the system dynamically enabling or disabling cores or
534 * by some scheduler changing the affinities. Thus we can not assume
535 * that the result is the same on all ranks within a node
536 * when running this code at different times.
538 bool allBitsAreSet = true;
539 for (int i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
541 allBitsAreSet = allBitsAreSet && (CPU_ISSET(i, &mask_current) != 0);
545 fprintf(debug, "%s affinity mask found\n", allBitsAreSet ? "Default" : "Non-default");
549 detectedDefaultAffinityMask = false;
553 GMX_UNUSED_VALUE(nthreads_hw_avail);
554 #endif /* HAVE_SCHED_AFFINITY */
557 int mpiIsInitialized;
558 MPI_Initialized(&mpiIsInitialized);
559 if (mpiIsInitialized)
561 bool maskToReduce = detectedDefaultAffinityMask;
562 MPI_Allreduce(&maskToReduce, &detectedDefaultAffinityMask, 1, MPI_C_BOOL, MPI_LAND, MPI_COMM_WORLD);
566 return detectedDefaultAffinityMask;
569 /* Check the process affinity mask and if it is found to be non-zero,
570 * will honor it and disable mdrun internal affinity setting.
571 * Note that this will only work on Linux as we use a GNU feature.
573 void gmx_check_thread_affinity_set(const gmx::MDLogger& mdlog,
574 gmx_hw_opt_t* hw_opt,
575 int gmx_unused nthreads_hw_avail,
576 gmx_bool bAfterOpenmpInit)
578 GMX_RELEASE_ASSERT(hw_opt, "hw_opt must be a non-NULL pointer");
580 if (!bAfterOpenmpInit)
582 /* Check for externally set OpenMP affinity and turn off internal
583 * pinning if any is found. We need to do this check early to tell
584 * thread-MPI whether it should do pinning when spawning threads.
585 * TODO: the above no longer holds, we should move these checks later
587 if (hw_opt->threadAffinity != ThreadAffinity::Off)
590 if (!gmx_omp_check_thread_affinity(&message))
592 /* We only pin automatically with totNumThreadsIsAuto=true */
593 if (hw_opt->threadAffinity == ThreadAffinity::On || hw_opt->totNumThreadsIsAuto)
595 GMX_LOG(mdlog.warning).asParagraph().appendText(message);
598 hw_opt->threadAffinity = ThreadAffinity::Off;
603 if (!detectDefaultAffinityMask(nthreads_hw_avail))
605 if (hw_opt->threadAffinity == ThreadAffinity::Auto)
607 if (!bAfterOpenmpInit)
609 GMX_LOG(mdlog.warning)
612 "Non-default thread affinity set, disabling internal thread "
617 GMX_LOG(mdlog.warning)
620 "Non-default thread affinity set probably by the OpenMP library,\n"
621 "disabling internal thread affinity");
623 hw_opt->threadAffinity = ThreadAffinity::Off;
627 /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
628 if (bAfterOpenmpInit)
630 GMX_LOG(mdlog.warning)
632 .appendTextFormatted("Overriding thread affinity set outside %s",
633 gmx::getProgramContext().displayName());