Replace scoped_cptr with unique_cptr
[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, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "threadaffinity.h"
38
39 #include "config.h"
40
41 #include <cerrno>
42 #include <cstdio>
43 #include <cstring>
44
45 #if HAVE_SCHED_AFFINITY
46 #  include <sched.h>
47 #  include <sys/syscall.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/programcontext.h"
63 #include "gromacs/utility/smalloc.h"
64 #include "gromacs/utility/unique_cptr.h"
65
66 namespace
67 {
68
69 class DefaultThreadAffinityAccess : public gmx::IThreadAffinityAccess
70 {
71     public:
72         virtual bool isThreadAffinitySupported() const
73         {
74             return tMPI_Thread_setaffinity_support() == TMPI_SETAFFINITY_SUPPORT_YES;
75         }
76         virtual bool setCurrentThreadAffinityToCore(int core)
77         {
78             const int ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
79             return ret == 0;
80         }
81 };
82
83 //! Global instance of DefaultThreadAffinityAccess
84 DefaultThreadAffinityAccess g_defaultAffinityAccess;
85
86 } // namespace
87
88 gmx::IThreadAffinityAccess::~IThreadAffinityAccess()
89 {
90 }
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),
100                    cr->mpi_comm_mysim);
101         return SIMMASTER(cr) ? (globalValue != 0) : invalidLocally;
102     }
103 #else
104     GMX_UNUSED_VALUE(cr);
105 #endif
106     return invalidLocally;
107 }
108
109 static bool
110 get_thread_affinity_layout(FILE *fplog, const gmx::MDLogger &mdlog,
111                            const t_commrec *cr,
112                            const gmx::HardwareTopology &hwTop,
113                            int   threads,
114                            bool  automatic,
115                            int pin_offset, int * pin_stride,
116                            int **localityOrder)
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 (auto &s : hwTop.machine().sockets)
143         {
144             for (auto &c : s.cores)
145             {
146                 for (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  = NULL;
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).asParagraph().appendText(
170                 "NOTE: No information on available cores, thread pinning disabled.");
171         alreadyWarned = true;
172     }
173     bool validLayout = !invalidValue;
174
175     if (automatic)
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).asParagraph().appendText(
182                     "NOTE: The number of threads is not equal to the number of (logical) cores\n"
183                     "      and the -pin option is set to auto: will not pin thread to cores.\n"
184                     "      This can lead to significant performance degradation.\n"
185                     "      Consider using -pin on (and -pinoffset in case you run multiple jobs).");
186             alreadyWarned = true;
187         }
188         validLayout = validLayout && !invalidValue;
189     }
190
191     invalidValue = (threads > hwThreads);
192     if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
193     {
194         GMX_LOG(mdlog.warning).asParagraph().appendText(
195                 "NOTE: Oversubscribing the CPU, will not pin threads");
196         alreadyWarned = true;
197     }
198     validLayout = validLayout && !invalidValue;
199
200     invalidValue = (pin_offset + threads > hwThreads);
201     if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
202     {
203         GMX_LOG(mdlog.warning).asParagraph().appendText(
204                 "WARNING: Requested offset too large for available cores, thread pinning disabled.");
205         alreadyWarned = true;
206
207     }
208     validLayout = validLayout && !invalidValue;
209
210     invalidValue   = false;
211     /* do we need to choose the pinning stride? */
212     bPickPinStride = (*pin_stride == 0);
213
214     if (bPickPinStride)
215     {
216         if (haveTopology && pin_offset + threads*hwThreadsPerCore <= hwThreads)
217         {
218             /* Put one thread on each physical core */
219             *pin_stride = hwThreadsPerCore;
220         }
221         else
222         {
223             /* We don't know if we have SMT, and if we do, we don't know
224              * if hw threads in the same physical core are consecutive.
225              * Without SMT the pinning layout should not matter too much.
226              * so we assume a consecutive layout and maximally spread out"
227              * the threads at equal threads per core.
228              * Note that IBM is the major non-x86 case with cpuid support
229              * and probably threads are already pinned by the queuing system,
230              * so we wouldn't end up here in the first place.
231              */
232             *pin_stride = (hwThreads - pin_offset)/threads;
233         }
234     }
235     else
236     {
237         /* Check the placement of the thread with the largest index to make sure
238          * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
239         invalidValue = (pin_offset + (threads-1)*(*pin_stride) >= hwThreads);
240     }
241     if (invalidWithinSimulation(cr, invalidValue) && !alreadyWarned)
242     {
243         /* We are oversubscribing, don't pin */
244         GMX_LOG(mdlog.warning).asParagraph().appendText(
245                 "WARNING: Requested stride too large for available cores, thread pinning disabled.");
246     }
247     validLayout = validLayout && !invalidValue;
248
249     if (validLayout && fplog != NULL)
250     {
251         fprintf(fplog, "Pinning threads with a%s logical core stride of %d\n",
252                 bPickPinStride ? "n auto-selected" : " user-specified",
253                 *pin_stride);
254     }
255
256     return validLayout;
257 }
258
259 static bool set_affinity(const t_commrec *cr, int nthread_local, int thread0_id_node,
260                          int offset, int core_pinning_stride, int *localityOrder,
261                          gmx::IThreadAffinityAccess *affinityAccess)
262 {
263     // Set the per-thread affinity. In order to be able to check the success
264     // of affinity settings, we will set nth_affinity_set to 1 on threads
265     // where the affinity setting succeded and to 0 where it failed.
266     // Reducing these 0/1 values over the threads will give the total number
267     // of threads on which we succeeded.
268
269     // To avoid warnings from the static analyzer we initialize nth_affinity_set
270     // to zero outside the OpenMP block, and then add to it inside the block.
271     // The value will still always be 0 or 1 from each thread.
272     int nth_affinity_set = 0;
273 #pragma omp parallel num_threads(nthread_local) reduction(+:nth_affinity_set)
274     {
275         try
276         {
277             int      thread_id, thread_id_node;
278             int      index, core;
279
280             thread_id      = gmx_omp_get_thread_num();
281             thread_id_node = thread0_id_node + thread_id;
282             index          = offset + thread_id_node*core_pinning_stride;
283             if (localityOrder != nullptr)
284             {
285                 core = localityOrder[index];
286             }
287             else
288             {
289                 core = index;
290             }
291
292             const bool ret = affinityAccess->setCurrentThreadAffinityToCore(core);
293
294             /* store the per-thread success-values of the setaffinity */
295             nth_affinity_set += (ret ? 1 : 0);
296
297             if (debug)
298             {
299                 fprintf(debug, "On rank %2d, thread %2d, index %2d, core %2d the affinity setting returned %d\n",
300                         cr->nodeid, gmx_omp_get_thread_num(), index, core, ret ? 1 : 0);
301             }
302         }
303         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
304     }
305
306     if (nth_affinity_set > nthread_local)
307     {
308         char msg[STRLEN];
309
310         sprintf(msg, "Looks like we have set affinity for more threads than "
311                 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
312         gmx_incons(msg);
313     }
314
315     /* check & warn if some threads failed to set their affinities */
316     const bool allAffinitiesSet = (nth_affinity_set == nthread_local);
317     if (!allAffinitiesSet)
318     {
319         char sbuf1[STRLEN], sbuf2[STRLEN];
320
321         /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
322         sbuf1[0] = sbuf2[0] = '\0';
323         /* Only add rank info if we have more than one rank. */
324         if (cr->nnodes > 1)
325         {
326 #if GMX_MPI
327 #if GMX_THREAD_MPI
328             sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
329 #else           /* GMX_LIB_MPI */
330             sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
331 #endif
332 #endif          /* GMX_MPI */
333         }
334
335         if (nthread_local > 1)
336         {
337             sprintf(sbuf2, "for %d/%d thread%s ",
338                     nthread_local - nth_affinity_set, nthread_local,
339                     nthread_local > 1 ? "s" : "");
340         }
341
342         fprintf(stderr, "NOTE: %sAffinity setting %sfailed.\n", sbuf1, sbuf2);
343     }
344     return allAffinitiesSet;
345 }
346
347 /* Set CPU affinity. Can be important for performance.
348    On some systems (e.g. Cray) CPU Affinity is set by default.
349    But default assigning doesn't work (well) with only some ranks
350    having threads. This causes very low performance.
351    External tools have cumbersome syntax for setting affinity
352    in the case that only some ranks have threads.
353    Thus it is important that GROMACS sets the affinity internally
354    if only PME is using threads.
355  */
356 void
357 gmx_set_thread_affinity(FILE                        *fplog,
358                         const gmx::MDLogger         &mdlog,
359                         const t_commrec             *cr,
360                         const gmx_hw_opt_t          *hw_opt,
361                         const gmx::HardwareTopology &hwTop,
362                         int                          nthread_local,
363                         gmx::IThreadAffinityAccess  *affinityAccess)
364 {
365     int        thread0_id_node, nthread_node;
366     int *      localityOrder = nullptr;
367
368     if (hw_opt->thread_affinity == threadaffOFF)
369     {
370         /* Nothing to do */
371         return;
372     }
373
374     if (affinityAccess == nullptr)
375     {
376         affinityAccess = &g_defaultAffinityAccess;
377     }
378
379     /* If the tMPI thread affinity setting is not supported encourage the user
380      * to report it as it's either a bug or an exotic platform which we might
381      * want to support. */
382     if (!affinityAccess->isThreadAffinitySupported())
383     {
384         /* we know Mac OS & BlueGene do not support setting thread affinity, so there's
385            no point in warning the user in that case. In any other case
386            the user might be able to do something about it. */
387 #if !defined(__APPLE__) && !defined(__bg__)
388         GMX_LOG(mdlog.warning).asParagraph().appendText(
389                 "NOTE: Cannot set thread affinities on the current platform.");
390 #endif  /* __APPLE__ */
391         return;
392     }
393
394     /* map the current process to cores */
395     thread0_id_node = 0;
396     nthread_node    = nthread_local;
397 #if GMX_MPI
398     if (PAR(cr) || MULTISIM(cr))
399     {
400         /* We need to determine a scan of the thread counts in this
401          * compute node.
402          */
403         MPI_Comm comm_intra;
404
405         MPI_Comm_split(MPI_COMM_WORLD,
406                        gmx_physicalnode_id_hash(), cr->rank_intranode,
407                        &comm_intra);
408         MPI_Scan(&nthread_local, &thread0_id_node, 1, MPI_INT, MPI_SUM, comm_intra);
409         /* MPI_Scan is inclusive, but here we need exclusive */
410         thread0_id_node -= nthread_local;
411         /* Get the total number of threads on this physical node */
412         MPI_Allreduce(&nthread_local, &nthread_node, 1, MPI_INT, MPI_SUM, comm_intra);
413         MPI_Comm_free(&comm_intra);
414     }
415 #endif
416
417     int  offset              = hw_opt->core_pinning_offset;
418     int  core_pinning_stride = hw_opt->core_pinning_stride;
419     if (offset != 0)
420     {
421         GMX_LOG(mdlog.warning).appendTextFormatted("Applying core pinning offset %d", offset);
422     }
423
424     bool automatic = (hw_opt->thread_affinity == threadaffAUTO);
425     bool validLayout
426         = get_thread_affinity_layout(fplog, mdlog, cr, hwTop, nthread_node, automatic,
427                                      offset, &core_pinning_stride, &localityOrder);
428     const gmx::sfree_guard  localityOrderGuard(localityOrder);
429
430     bool                    allAffinitiesSet;
431     if (validLayout)
432     {
433         allAffinitiesSet = set_affinity(cr, nthread_local, thread0_id_node,
434                                         offset, core_pinning_stride, localityOrder,
435                                         affinityAccess);
436     }
437     else
438     {
439         // Produce the warning if any rank fails.
440         allAffinitiesSet = false;
441     }
442     if (invalidWithinSimulation(cr, !allAffinitiesSet))
443     {
444         GMX_LOG(mdlog.warning).asParagraph().appendText(
445                 "NOTE: Thread affinity setting failed. This can cause performance degradation.\n"
446                 "      If you think your settings are correct, ask on the gmx-users list.");
447     }
448 }
449
450 /* Check the process affinity mask and if it is found to be non-zero,
451  * will honor it and disable mdrun internal affinity setting.
452  * Note that this will only work on Linux as we use a GNU feature.
453  */
454 void
455 gmx_check_thread_affinity_set(const gmx::MDLogger &mdlog,
456                               const t_commrec     *cr,
457                               gmx_hw_opt_t        *hw_opt,
458                               int  gmx_unused      nthreads_hw_avail,
459                               gmx_bool             bAfterOpenmpInit)
460 {
461     GMX_RELEASE_ASSERT(hw_opt, "hw_opt must be a non-NULL pointer");
462
463     if (!bAfterOpenmpInit)
464     {
465         /* Check for externally set OpenMP affinity and turn off internal
466          * pinning if any is found. We need to do this check early to tell
467          * thread-MPI whether it should do pinning when spawning threads.
468          * TODO: the above no longer holds, we should move these checks later
469          */
470         if (hw_opt->thread_affinity != threadaffOFF)
471         {
472             char *message;
473             if (!gmx_omp_check_thread_affinity(&message))
474             {
475                 /* TODO: with -pin auto we should only warn when using all cores */
476                 GMX_LOG(mdlog.warning).asParagraph().appendText(message);
477                 sfree(message);
478                 hw_opt->thread_affinity = threadaffOFF;
479             }
480         }
481
482         /* With thread-MPI this is needed as pinning might get turned off,
483          * which needs to be known before starting thread-MPI.
484          * With thread-MPI hw_opt is processed here on the master rank
485          * and passed to the other ranks later, so we only do this on master.
486          */
487         if (!SIMMASTER(cr))
488         {
489             return;
490         }
491 #if !GMX_THREAD_MPI
492         return;
493 #endif
494     }
495
496 #if HAVE_SCHED_AFFINITY
497     int       ret;
498     cpu_set_t mask_current;
499
500     if (hw_opt->thread_affinity == threadaffOFF)
501     {
502         /* internal affinity setting is off, don't bother checking process affinity */
503         return;
504     }
505
506     CPU_ZERO(&mask_current);
507     if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
508     {
509         /* failed to query affinity mask, will just return */
510         if (debug)
511         {
512             fprintf(debug, "Failed to query affinity mask (error %d)", ret);
513         }
514         return;
515     }
516
517     /* Before proceeding with the actual check, make sure that the number of
518      * detected CPUs is >= the CPUs in the current set.
519      * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
520 #ifdef CPU_COUNT
521     if (nthreads_hw_avail < CPU_COUNT(&mask_current))
522     {
523         if (debug)
524         {
525             fprintf(debug, "%d hardware threads detected, but %d was returned by CPU_COUNT",
526                     nthreads_hw_avail, CPU_COUNT(&mask_current));
527         }
528         return;
529     }
530 #endif /* CPU_COUNT */
531
532     gmx_bool bAllSet = TRUE;
533     for (int i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
534     {
535         bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
536     }
537
538 #if GMX_MPI
539     int isInitialized;
540     MPI_Initialized(&isInitialized);
541     /* Before OpenMP initialization, thread-MPI is not yet initialized.
542      * With thread-MPI bAllSet will be the same on all MPI-ranks, but the
543      * MPI_Allreduce then functions as a barrier before setting affinities.
544      */
545     if (isInitialized)
546     {
547         gmx_bool  bAllSet_All;
548
549         MPI_Allreduce(&bAllSet, &bAllSet_All, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
550         bAllSet = bAllSet_All;
551     }
552 #endif
553
554     if (!bAllSet)
555     {
556         if (hw_opt->thread_affinity == threadaffAUTO)
557         {
558             if (!bAfterOpenmpInit)
559             {
560                 GMX_LOG(mdlog.warning).asParagraph().appendText(
561                         "Non-default thread affinity set, disabling internal thread affinity");
562             }
563             else
564             {
565                 GMX_LOG(mdlog.warning).asParagraph().appendText(
566                         "Non-default thread affinity set probably by the OpenMP library,\n"
567                         "disabling internal thread affinity");
568             }
569             hw_opt->thread_affinity = threadaffOFF;
570         }
571         else
572         {
573             /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
574             if (bAfterOpenmpInit)
575             {
576                 GMX_LOG(mdlog.warning).asParagraph().appendTextFormatted(
577                         "Overriding thread affinity set outside %s",
578                         gmx::getProgramContext().displayName());
579             }
580         }
581
582         if (debug)
583         {
584             fprintf(debug, "Non-default affinity mask found\n");
585         }
586     }
587     else
588     {
589         if (debug)
590         {
591             fprintf(debug, "Default affinity mask found\n");
592         }
593     }
594 #endif /* HAVE_SCHED_AFFINITY */
595 }