Apply clang-format-11
[alexxy/gromacs.git] / src / gromacs / mdrunutility / threadaffinity.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, 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.
9  *
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.
14  *
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.
19  *
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.
24  *
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.
32  *
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.
35  */
36 #include "gmxpre.h"
37
38 #include "threadaffinity.h"
39
40 #include "config.h"
41
42 #include <cerrno>
43 #include <cstdio>
44 #include <cstring>
45
46 #if HAVE_SCHED_AFFINITY
47 #    include <sched.h>
48 #endif
49
50 #include "thread_mpi/threads.h"
51
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"
66
67 namespace
68 {
69
70 class DefaultThreadAffinityAccess : public gmx::IThreadAffinityAccess
71 {
72 public:
73     bool isThreadAffinitySupported() const override
74     {
75         return tMPI_Thread_setaffinity_support() == TMPI_SETAFFINITY_SUPPORT_YES;
76     }
77     bool setCurrentThreadAffinityToCore(int core) override
78     {
79         const int ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
80         return ret == 0;
81     }
82 };
83
84 //! Global instance of DefaultThreadAffinityAccess
85 // NOLINTNEXTLINE(cppcoreguidelines-avoid-non-const-global-variables)
86 DefaultThreadAffinityAccess g_defaultAffinityAccess;
87
88 } // namespace
89
90 gmx::IThreadAffinityAccess::~IThreadAffinityAccess() {}
91
92 static bool invalidWithinSimulation(const t_commrec* cr, bool invalidLocally)
93 {
94 #if GMX_MPI
95     if (cr->nnodes > 1)
96     {
97         int value = invalidLocally ? 1 : 0;
98         int globalValue;
99         MPI_Reduce(&value, &globalValue, 1, MPI_INT, MPI_LOR, MASTERRANK(cr), cr->mpi_comm_mysim);
100         return SIMMASTER(cr) ? (globalValue != 0) : invalidLocally;
101     }
102 #else
103     GMX_UNUSED_VALUE(cr);
104 #endif
105     return invalidLocally;
106 }
107
108 static bool get_thread_affinity_layout(const gmx::MDLogger&         mdlog,
109                                        const t_commrec*             cr,
110                                        const gmx::HardwareTopology& hwTop,
111                                        int                          threads,
112                                        bool  affinityIsAutoAndNumThreadsIsNotAuto,
113                                        int   pin_offset,
114                                        int*  pin_stride,
115                                        int** localityOrder,
116                                        bool* issuedWarning)
117 {
118     int  hwThreads;
119     int  hwThreadsPerCore = 0;
120     bool bPickPinStride;
121     bool haveTopology;
122     bool invalidValue;
123
124     haveTopology = (hwTop.supportLevel() >= gmx::HardwareTopology::SupportLevel::Basic);
125
126     if (pin_offset < 0)
127     {
128         gmx_fatal(FARGS, "Negative thread pinning offset requested");
129     }
130     if (*pin_stride < 0)
131     {
132         gmx_fatal(FARGS, "Negative thread pinning stride requested");
133     }
134
135     if (haveTopology)
136     {
137         hwThreads = hwTop.machine().logicalProcessorCount;
138         // Just use the value for the first core
139         hwThreadsPerCore = hwTop.machine().sockets[0].cores[0].hwThreads.size();
140         snew(*localityOrder, hwThreads);
141         int i = 0;
142         for (const auto& s : hwTop.machine().sockets)
143         {
144             for (const auto& c : s.cores)
145             {
146                 for (const auto& t : c.hwThreads)
147                 {
148                     (*localityOrder)[i++] = t.logicalProcessorId;
149                 }
150             }
151         }
152     }
153     else
154     {
155         /* topology information not available or invalid, ignore it */
156         hwThreads      = hwTop.machine().logicalProcessorCount;
157         *localityOrder = nullptr;
158     }
159     // Only warn about the first problem per node.  Otherwise, the first test
160     // failing would essentially always cause also the other problems get
161     // reported, leading to bogus warnings.  The order in the conditionals
162     // with this variable is important, since the MPI_Reduce() in
163     // invalidWithinSimulation() needs to always happen.
164     bool alreadyWarned = false;
165     invalidValue       = (hwThreads <= 0);
166     if (invalidWithinSimulation(cr, invalidValue))
167     {
168         /* We don't know anything about the hardware, don't pin */
169         GMX_LOG(mdlog.warning)
170                 .asParagraph()
171                 .appendText("NOTE: No information on available cores, thread pinning disabled.");
172         alreadyWarned = true;
173     }
174     bool validLayout = !invalidValue;
175
176     if (affinityIsAutoAndNumThreadsIsNotAuto)
177     {
178         invalidValue = (threads != hwThreads);
179         bool warn    = (invalidValue && threads > 1 && threads < hwThreads);
180         if (invalidWithinSimulation(cr, warn) && !alreadyWarned)
181         {
182             GMX_LOG(mdlog.warning)
183                     .asParagraph()
184                     .appendText(
185                             "NOTE: The number of threads is not equal to the number of (logical) "
186                             "cores\n"
187                             "      and the -pin option is set to auto: will not pin threads to "
188                             "cores.\n"
189                             "      This can lead to significant performance degradation.\n"
190                             "      Consider using -pin on (and -pinoffset in case you run multiple "
191                             "jobs).");
192             alreadyWarned = true;
193         }
194         validLayout = validLayout && !invalidValue;
195     }
196
197     invalidValue = (threads > hwThreads);
198     if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
199     {
200         GMX_LOG(mdlog.warning)
201                 .asParagraph()
202                 .appendText("NOTE: Oversubscribing the CPU, will not pin threads");
203         alreadyWarned = true;
204     }
205     validLayout = validLayout && !invalidValue;
206
207     invalidValue = (pin_offset + threads > hwThreads);
208     if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
209     {
210         GMX_LOG(mdlog.warning)
211                 .asParagraph()
212                 .appendText(
213                         "WARNING: Requested offset too large for available cores, thread pinning "
214                         "disabled.");
215         alreadyWarned = true;
216     }
217     validLayout = validLayout && !invalidValue;
218
219     invalidValue = false;
220     /* do we need to choose the pinning stride? */
221     bPickPinStride = (*pin_stride == 0);
222
223     if (bPickPinStride)
224     {
225         if (haveTopology && pin_offset + threads * hwThreadsPerCore <= hwThreads)
226         {
227             /* Put one thread on each physical core */
228             *pin_stride = hwThreadsPerCore;
229         }
230         else
231         {
232             /* We don't know if we have SMT, and if we do, we don't know
233              * if hw threads in the same physical core are consecutive.
234              * Without SMT the pinning layout should not matter too much.
235              * so we assume a consecutive layout and maximally spread out"
236              * the threads at equal threads per core.
237              * Note that IBM is the major non-x86 case with cpuid support
238              * and probably threads are already pinned by the queuing system,
239              * so we wouldn't end up here in the first place.
240              */
241             *pin_stride = (hwThreads - pin_offset) / threads;
242         }
243     }
244     else
245     {
246         /* Check the placement of the thread with the largest index to make sure
247          * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
248         invalidValue = (pin_offset + (threads - 1) * (*pin_stride) >= hwThreads);
249     }
250     if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
251     {
252         /* We are oversubscribing, don't pin */
253         GMX_LOG(mdlog.warning)
254                 .asParagraph()
255                 .appendText(
256                         "WARNING: Requested stride too large for available cores, thread pinning "
257                         "disabled.");
258         alreadyWarned = true;
259     }
260     validLayout = validLayout && !invalidValue;
261
262     if (validLayout)
263     {
264         GMX_LOG(mdlog.info)
265                 .appendTextFormatted("Pinning threads with a%s logical core stride of %d",
266                                      bPickPinStride ? "n auto-selected" : " user-specified",
267                                      *pin_stride);
268     }
269
270     *issuedWarning = alreadyWarned;
271
272     return validLayout;
273 }
274
275 static bool set_affinity(const t_commrec*            cr,
276                          int                         nthread_local,
277                          int                         intraNodeThreadOffset,
278                          int                         offset,
279                          int                         core_pinning_stride,
280                          const int*                  localityOrder,
281                          gmx::IThreadAffinityAccess* affinityAccess)
282 {
283     // Set the per-thread affinity. In order to be able to check the success
284     // of affinity settings, we will set nth_affinity_set to 1 on threads
285     // where the affinity setting succeded and to 0 where it failed.
286     // Reducing these 0/1 values over the threads will give the total number
287     // of threads on which we succeeded.
288
289     // To avoid warnings from the static analyzer we initialize nth_affinity_set
290     // to zero outside the OpenMP block, and then add to it inside the block.
291     // The value will still always be 0 or 1 from each thread.
292     int nth_affinity_set = 0;
293 #pragma omp parallel num_threads(nthread_local) reduction(+ : nth_affinity_set)
294     {
295         try
296         {
297             int thread_id, thread_id_node;
298             int index, core;
299
300             thread_id      = gmx_omp_get_thread_num();
301             thread_id_node = intraNodeThreadOffset + thread_id;
302             index          = offset + thread_id_node * core_pinning_stride;
303             if (localityOrder != nullptr)
304             {
305                 core = localityOrder[index];
306             }
307             else
308             {
309                 core = index;
310             }
311
312             const bool ret = affinityAccess->setCurrentThreadAffinityToCore(core);
313
314             /* store the per-thread success-values of the setaffinity */
315             nth_affinity_set += (ret ? 1 : 0);
316
317             if (debug)
318             {
319                 fprintf(debug,
320                         "On rank %2d, thread %2d, index %2d, core %2d the affinity setting "
321                         "returned %d\n",
322                         cr->nodeid,
323                         gmx_omp_get_thread_num(),
324                         index,
325                         core,
326                         ret ? 1 : 0);
327             }
328         }
329         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR
330     }
331
332     if (nth_affinity_set > nthread_local)
333     {
334         char msg[STRLEN];
335
336         sprintf(msg,
337                 "Looks like we have set affinity for more threads than "
338                 "we have (%d > %d)!\n",
339                 nth_affinity_set,
340                 nthread_local);
341         gmx_incons(msg);
342     }
343
344     /* check & warn if some threads failed to set their affinities */
345     const bool allAffinitiesSet = (nth_affinity_set == nthread_local);
346     if (!allAffinitiesSet)
347     {
348         char sbuf1[STRLEN], sbuf2[STRLEN];
349
350         /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
351         sbuf1[0] = sbuf2[0] = '\0';
352         /* Only add rank info if we have more than one rank. */
353         if (cr->nnodes > 1)
354         {
355 #if GMX_MPI
356 #    if GMX_THREAD_MPI
357             sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
358 #    else /* GMX_LIB_MPI */
359             sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
360 #    endif
361 #endif /* GMX_MPI */
362         }
363
364         if (nthread_local > 1)
365         {
366             sprintf(sbuf2,
367                     "for %d/%d thread%s ",
368                     nthread_local - nth_affinity_set,
369                     nthread_local,
370                     nthread_local > 1 ? "s" : "");
371         }
372
373         // TODO: This output should also go through mdlog.
374         fprintf(stderr, "NOTE: %sAffinity setting %sfailed.\n", sbuf1, sbuf2);
375     }
376     return allAffinitiesSet;
377 }
378
379 void analyzeThreadsOnThisNode(const gmx::PhysicalNodeCommunicator& physicalNodeComm,
380                               int                                  numThreadsOnThisRank,
381                               int*                                 numThreadsOnThisNode,
382                               int*                                 intraNodeThreadOffset)
383 {
384     *intraNodeThreadOffset = 0;
385     *numThreadsOnThisNode  = numThreadsOnThisRank;
386 #if GMX_MPI
387     if (physicalNodeComm.size_ > 1)
388     {
389         /* We need to determine a scan of the thread counts in this
390          * compute node. */
391         MPI_Scan(&numThreadsOnThisRank, intraNodeThreadOffset, 1, MPI_INT, MPI_SUM, physicalNodeComm.comm_);
392         /* MPI_Scan is inclusive, but here we need exclusive */
393         *intraNodeThreadOffset -= numThreadsOnThisRank;
394         /* Get the total number of threads on this physical node */
395         MPI_Allreduce(
396                 &numThreadsOnThisRank, numThreadsOnThisNode, 1, MPI_INT, MPI_SUM, physicalNodeComm.comm_);
397     }
398 #else
399     GMX_UNUSED_VALUE(physicalNodeComm);
400 #endif
401 }
402
403 /* Set CPU affinity. Can be important for performance.
404    On some systems (e.g. Cray) CPU Affinity is set by default.
405    But default assigning doesn't work (well) with only some ranks
406    having threads. This causes very low performance.
407    External tools have cumbersome syntax for setting affinity
408    in the case that only some ranks have threads.
409    Thus it is important that GROMACS sets the affinity internally
410    if only PME is using threads.
411  */
412 void gmx_set_thread_affinity(const gmx::MDLogger&         mdlog,
413                              const t_commrec*             cr,
414                              const gmx_hw_opt_t*          hw_opt,
415                              const gmx::HardwareTopology& hwTop,
416                              int                          numThreadsOnThisRank,
417                              int                          numThreadsOnThisNode,
418                              int                          intraNodeThreadOffset,
419                              gmx::IThreadAffinityAccess*  affinityAccess)
420 {
421     int* localityOrder = nullptr;
422
423     if (hw_opt->threadAffinity == ThreadAffinity::Off)
424     {
425         /* Nothing to do */
426         return;
427     }
428
429     if (affinityAccess == nullptr)
430     {
431         affinityAccess = &g_defaultAffinityAccess;
432     }
433
434     /* If the tMPI thread affinity setting is not supported encourage the user
435      * to report it as it's either a bug or an exotic platform which we might
436      * want to support. */
437     if (!affinityAccess->isThreadAffinitySupported())
438     {
439         /* we know Mac OS does not support setting thread affinity, so there's
440            no point in warning the user in that case. In any other case
441            the user might be able to do something about it. */
442 #if !defined(__APPLE__)
443         GMX_LOG(mdlog.warning)
444                 .asParagraph()
445                 .appendText("NOTE: Cannot set thread affinities on the current platform.");
446 #endif /* __APPLE__ */
447         return;
448     }
449
450     int offset              = hw_opt->core_pinning_offset;
451     int core_pinning_stride = hw_opt->core_pinning_stride;
452     if (offset != 0)
453     {
454         GMX_LOG(mdlog.warning).appendTextFormatted("Applying core pinning offset %d", offset);
455     }
456
457     bool affinityIsAutoAndNumThreadsIsNotAuto =
458             (hw_opt->threadAffinity == ThreadAffinity::Auto && !hw_opt->totNumThreadsIsAuto);
459     bool                   issuedWarning;
460     bool                   validLayout = get_thread_affinity_layout(mdlog,
461                                                   cr,
462                                                   hwTop,
463                                                   numThreadsOnThisNode,
464                                                   affinityIsAutoAndNumThreadsIsNotAuto,
465                                                   offset,
466                                                   &core_pinning_stride,
467                                                   &localityOrder,
468                                                   &issuedWarning);
469     const gmx::sfree_guard localityOrderGuard(localityOrder);
470
471     bool allAffinitiesSet;
472     if (validLayout)
473     {
474         allAffinitiesSet = set_affinity(
475                 cr, numThreadsOnThisRank, intraNodeThreadOffset, offset, core_pinning_stride, localityOrder, affinityAccess);
476     }
477     else
478     {
479         // Produce the warning if any rank fails.
480         allAffinitiesSet = false;
481     }
482     if (invalidWithinSimulation(cr, !allAffinitiesSet) && !issuedWarning)
483     {
484         GMX_LOG(mdlog.warning).asParagraph().appendText("NOTE: Thread affinity was not set.");
485     }
486 }
487
488 /* Detects and returns whether we have the default affinity mask
489  *
490  * Returns true when we can query thread affinities and CPU count is
491  * consistent and we have default affinity mask on all ranks.
492  *
493  * Should be called simultaneously by all MPI ranks.
494  */
495 static bool detectDefaultAffinityMask(const int nthreads_hw_avail)
496 {
497     bool detectedDefaultAffinityMask = true;
498
499 #if HAVE_SCHED_AFFINITY
500     cpu_set_t mask_current;
501     CPU_ZERO(&mask_current);
502     int ret;
503     if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
504     {
505         /* failed to query affinity mask, will just return */
506         if (debug)
507         {
508             fprintf(debug, "Failed to query affinity mask (error %d)", ret);
509         }
510         detectedDefaultAffinityMask = false;
511     }
512
513     /* Before proceeding with the actual check, make sure that the number of
514      * detected CPUs is >= the CPUs in the current set.
515      * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
516 #    ifdef CPU_COUNT
517     if (detectedDefaultAffinityMask && nthreads_hw_avail < CPU_COUNT(&mask_current))
518     {
519         if (debug)
520         {
521             fprintf(debug,
522                     "%d hardware threads detected, but %d was returned by CPU_COUNT",
523                     nthreads_hw_avail,
524                     CPU_COUNT(&mask_current));
525         }
526         detectedDefaultAffinityMask = false;
527     }
528 #    endif /* CPU_COUNT */
529
530     if (detectedDefaultAffinityMask)
531     {
532         /* Here we check whether all bits of the affinity mask are set.
533          * Note that this mask can change while the program is executing,
534          * e.g. by the system dynamically enabling or disabling cores or
535          * by some scheduler changing the affinities. Thus we can not assume
536          * that the result is the same on all ranks within a node
537          * when running this code at different times.
538          */
539         bool allBitsAreSet = true;
540         for (int i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
541         {
542             allBitsAreSet = allBitsAreSet && (CPU_ISSET(i, &mask_current) != 0);
543         }
544         if (debug)
545         {
546             fprintf(debug, "%s affinity mask found\n", allBitsAreSet ? "Default" : "Non-default");
547         }
548         if (!allBitsAreSet)
549         {
550             detectedDefaultAffinityMask = false;
551         }
552     }
553 #else
554     GMX_UNUSED_VALUE(nthreads_hw_avail);
555 #endif /* HAVE_SCHED_AFFINITY */
556
557 #if GMX_MPI
558     int mpiIsInitialized;
559     MPI_Initialized(&mpiIsInitialized);
560     if (mpiIsInitialized)
561     {
562         bool maskToReduce = detectedDefaultAffinityMask;
563         MPI_Allreduce(&maskToReduce, &detectedDefaultAffinityMask, 1, MPI_C_BOOL, MPI_LAND, MPI_COMM_WORLD);
564     }
565 #endif
566
567     return detectedDefaultAffinityMask;
568 }
569
570 /* Check the process affinity mask and if it is found to be non-zero,
571  * will honor it and disable mdrun internal affinity setting.
572  * Note that this will only work on Linux as we use a GNU feature.
573  */
574 void gmx_check_thread_affinity_set(const gmx::MDLogger& mdlog,
575                                    gmx_hw_opt_t*        hw_opt,
576                                    int gmx_unused       nthreads_hw_avail,
577                                    gmx_bool             bAfterOpenmpInit)
578 {
579     GMX_RELEASE_ASSERT(hw_opt, "hw_opt must be a non-NULL pointer");
580
581     if (!bAfterOpenmpInit)
582     {
583         /* Check for externally set OpenMP affinity and turn off internal
584          * pinning if any is found. We need to do this check early to tell
585          * thread-MPI whether it should do pinning when spawning threads.
586          * TODO: the above no longer holds, we should move these checks later
587          */
588         if (hw_opt->threadAffinity != ThreadAffinity::Off)
589         {
590             char* message;
591             if (!gmx_omp_check_thread_affinity(&message))
592             {
593                 /* We only pin automatically with totNumThreadsIsAuto=true */
594                 if (hw_opt->threadAffinity == ThreadAffinity::On || hw_opt->totNumThreadsIsAuto)
595                 {
596                     GMX_LOG(mdlog.warning).asParagraph().appendText(message);
597                 }
598                 sfree(message);
599                 hw_opt->threadAffinity = ThreadAffinity::Off;
600             }
601         }
602     }
603
604     if (!detectDefaultAffinityMask(nthreads_hw_avail))
605     {
606         if (hw_opt->threadAffinity == ThreadAffinity::Auto)
607         {
608             if (!bAfterOpenmpInit)
609             {
610                 GMX_LOG(mdlog.warning)
611                         .asParagraph()
612                         .appendText(
613                                 "Non-default thread affinity set, disabling internal thread "
614                                 "affinity");
615             }
616             else
617             {
618                 GMX_LOG(mdlog.warning)
619                         .asParagraph()
620                         .appendText(
621                                 "Non-default thread affinity set probably by the OpenMP library,\n"
622                                 "disabling internal thread affinity");
623             }
624             hw_opt->threadAffinity = ThreadAffinity::Off;
625         }
626         else
627         {
628             /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
629             if (bAfterOpenmpInit)
630             {
631                 GMX_LOG(mdlog.warning)
632                         .asParagraph()
633                         .appendTextFormatted("Overriding thread affinity set outside %s",
634                                              gmx::getProgramContext().displayName());
635             }
636         }
637     }
638 }