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