b43444468007852d8029d6c8dbcfa883ee807977
[alexxy/gromacs.git] / src / gromacs / gmxlib / gmx_thread_affinity.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
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 "gromacs/legacyheaders/gmx_thread_affinity.h"
38
39 #include "config.h"
40
41 #include <assert.h>
42 #include <errno.h>
43 #include <stdio.h>
44 #include <string.h>
45
46 #ifdef HAVE_SCHED_AFFINITY
47 #  include <sched.h>
48 #  include <sys/syscall.h>
49 #endif
50
51 #include "thread_mpi/threads.h"
52
53 #include "gromacs/legacyheaders/copyrite.h"
54 #include "gromacs/legacyheaders/gmx_cpuid.h"
55 #include "gromacs/legacyheaders/gmx_omp_nthreads.h"
56 #include "gromacs/legacyheaders/md_logging.h"
57 #include "gromacs/legacyheaders/typedefs.h"
58 #include "gromacs/legacyheaders/types/commrec.h"
59 #include "gromacs/legacyheaders/types/hw_info.h"
60 #include "gromacs/utility/basenetwork.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/fatalerror.h"
63 #include "gromacs/utility/gmxomp.h"
64 #include "gromacs/utility/smalloc.h"
65
66 static int
67 get_thread_affinity_layout(FILE *fplog,
68                            const t_commrec *cr,
69                            const gmx_hw_info_t * hwinfo,
70                            int nthreads,
71                            int pin_offset, int * pin_stride,
72                            const int **locality_order)
73 {
74     int         nhwthreads, npkg, ncores, nhwthreads_per_core, rc;
75     const int * pkg_id;
76     const int * core_id;
77     const int * hwthread_id;
78     gmx_bool    bPickPinStride;
79
80     if (pin_offset < 0)
81     {
82         gmx_fatal(FARGS, "Negative thread pinning offset requested");
83     }
84     if (*pin_stride < 0)
85     {
86         gmx_fatal(FARGS, "Negative thread pinning stride requested");
87     }
88
89     rc = gmx_cpuid_topology(hwinfo->cpuid_info, &nhwthreads, &npkg, &ncores,
90                             &nhwthreads_per_core,
91                             &pkg_id, &core_id, &hwthread_id, locality_order);
92
93     if (rc != 0)
94     {
95         /* topology information not available or invalid, ignore it */
96         nhwthreads      = hwinfo->nthreads_hw_avail;
97         *locality_order = NULL;
98
99         if (nhwthreads <= 0)
100         {
101             /* We don't know anything about the hardware, don't pin */
102             md_print_warn(cr, fplog,
103                           "NOTE: No information on available cores, thread pinning disabled.");
104
105             return -1;
106         }
107     }
108
109     if (nthreads > nhwthreads)
110     {
111         /* We are oversubscribing, don't pin */
112         md_print_warn(NULL, fplog,
113                       "NOTE: Oversubscribing the CPU, will not pin threads");
114
115         return -1;
116     }
117
118     if (pin_offset + nthreads > nhwthreads)
119     {
120         /* We are oversubscribing, don't pin */
121         md_print_warn(NULL, fplog,
122                       "WARNING: Requested offset too large for available cores, thread pinning disabled.");
123
124         return -1;
125     }
126
127
128     /* do we need to choose the pinning stride? */
129     bPickPinStride = (*pin_stride == 0);
130
131     if (bPickPinStride)
132     {
133         if (rc == 0 && pin_offset + nthreads*nhwthreads_per_core <= nhwthreads)
134         {
135             /* Put one thread on each physical core */
136             *pin_stride = nhwthreads_per_core;
137         }
138         else
139         {
140             /* We don't know if we have SMT, and if we do, we don't know
141              * if hw threads in the same physical core are consecutive.
142              * Without SMT the pinning layout should not matter too much.
143              * so we assume a consecutive layout and maximally spread out"
144              * the threads at equal threads per core.
145              * Note that IBM is the major non-x86 case with cpuid support
146              * and probably threads are already pinned by the queuing system,
147              * so we wouldn't end up here in the first place.
148              */
149             *pin_stride = (nhwthreads - pin_offset)/nthreads;
150         }
151     }
152     else
153     {
154         /* Check the placement of the thread with the largest index to make sure
155          * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
156         if (pin_offset + (nthreads-1)*(*pin_stride) >= nhwthreads)
157         {
158             /* We are oversubscribing, don't pin */
159             md_print_warn(NULL, fplog,
160                           "WARNING: Requested stride too large for available cores, thread pinning disabled.");
161
162             return -1;
163         }
164     }
165
166     if (fplog != NULL)
167     {
168         fprintf(fplog, "Pinning threads with a%s logical core stride of %d\n",
169                 bPickPinStride ? "n auto-selected" : " user-specified",
170                 *pin_stride);
171     }
172
173     return 0;
174 }
175
176 /* Set CPU affinity. Can be important for performance.
177    On some systems (e.g. Cray) CPU Affinity is set by default.
178    But default assigning doesn't work (well) with only some ranks
179    having threads. This causes very low performance.
180    External tools have cumbersome syntax for setting affinity
181    in the case that only some ranks have threads.
182    Thus it is important that GROMACS sets the affinity internally
183    if only PME is using threads.
184  */
185 void
186 gmx_set_thread_affinity(FILE                *fplog,
187                         const t_commrec     *cr,
188                         gmx_hw_opt_t        *hw_opt,
189                         const gmx_hw_info_t *hwinfo)
190 {
191     int        nth_affinity_set, thread0_id_node,
192                nthread_local, nthread_node, nthread_hw_max, nphyscore;
193     int        offset;
194     const int *locality_order;
195     int        rc;
196
197     if (hw_opt->thread_affinity == threadaffOFF)
198     {
199         /* Nothing to do */
200         return;
201     }
202
203     /* If the tMPI thread affinity setting is not supported encourage the user
204      * to report it as it's either a bug or an exotic platform which we might
205      * want to support. */
206     if (tMPI_Thread_setaffinity_support() != TMPI_SETAFFINITY_SUPPORT_YES)
207     {
208         /* we know Mac OS & BlueGene do not support setting thread affinity, so there's
209            no point in warning the user in that case. In any other case
210            the user might be able to do something about it. */
211 #if !defined(__APPLE__) && !defined(__bg__)
212         md_print_warn(NULL, fplog,
213                       "NOTE: Cannot set thread affinities on the current platform.");
214 #endif  /* __APPLE__ */
215         return;
216     }
217
218     /* threads on this MPI process or TMPI thread */
219     if (cr->duty & DUTY_PP)
220     {
221         nthread_local = gmx_omp_nthreads_get(emntNonbonded);
222     }
223     else
224     {
225         nthread_local = gmx_omp_nthreads_get(emntPME);
226     }
227
228     /* map the current process to cores */
229     thread0_id_node = 0;
230     nthread_node    = nthread_local;
231 #ifdef GMX_MPI
232     if (PAR(cr) || MULTISIM(cr))
233     {
234         /* We need to determine a scan of the thread counts in this
235          * compute node.
236          */
237         MPI_Comm comm_intra;
238
239         MPI_Comm_split(MPI_COMM_WORLD,
240                        gmx_physicalnode_id_hash(), cr->rank_intranode,
241                        &comm_intra);
242         MPI_Scan(&nthread_local, &thread0_id_node, 1, MPI_INT, MPI_SUM, comm_intra);
243         /* MPI_Scan is inclusive, but here we need exclusive */
244         thread0_id_node -= nthread_local;
245         /* Get the total number of threads on this physical node */
246         MPI_Allreduce(&nthread_local, &nthread_node, 1, MPI_INT, MPI_SUM, comm_intra);
247         MPI_Comm_free(&comm_intra);
248     }
249 #endif
250
251     if (hw_opt->thread_affinity == threadaffAUTO &&
252         nthread_node != hwinfo->nthreads_hw_avail)
253     {
254         if (nthread_node > 1 && nthread_node < hwinfo->nthreads_hw_avail)
255         {
256             md_print_warn(cr, fplog,
257                           "NOTE: The number of threads is not equal to the number of (logical) cores\n"
258                           "      and the -pin option is set to auto: will not pin thread to cores.\n"
259                           "      This can lead to significant performance degradation.\n"
260                           "      Consider using -pin on (and -pinoffset in case you run multiple jobs).\n");
261         }
262
263         return;
264     }
265
266     offset = 0;
267     if (hw_opt->core_pinning_offset != 0)
268     {
269         offset = hw_opt->core_pinning_offset;
270         md_print_info(cr, fplog, "Applying core pinning offset %d\n", offset);
271     }
272
273     rc = get_thread_affinity_layout(fplog, cr, hwinfo,
274                                     nthread_node,
275                                     offset, &hw_opt->core_pinning_stride,
276                                     &locality_order);
277
278     if (rc != 0)
279     {
280         /* Incompatible layout, don't pin, warning was already issued */
281         return;
282     }
283
284     /* Set the per-thread affinity. In order to be able to check the success
285      * of affinity settings, we will set nth_affinity_set to 1 on threads
286      * where the affinity setting succeded and to 0 where it failed.
287      * Reducing these 0/1 values over the threads will give the total number
288      * of threads on which we succeeded.
289      */
290     nth_affinity_set = 0;
291 #pragma omp parallel num_threads(nthread_local) reduction(+:nth_affinity_set)
292     {
293         int      thread_id, thread_id_node;
294         int      index, core;
295         gmx_bool setaffinity_ret;
296
297         thread_id      = gmx_omp_get_thread_num();
298         thread_id_node = thread0_id_node + thread_id;
299         index          = offset + thread_id_node*hw_opt->core_pinning_stride;
300         if (locality_order != NULL)
301         {
302             core = locality_order[index];
303         }
304         else
305         {
306             core = index;
307         }
308
309         setaffinity_ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
310
311         /* store the per-thread success-values of the setaffinity */
312         nth_affinity_set = (setaffinity_ret == 0);
313
314         if (debug)
315         {
316             fprintf(debug, "On rank %2d, thread %2d, index %2d, core %2d the affinity setting returned %d\n",
317                     cr->nodeid, gmx_omp_get_thread_num(), index, core, setaffinity_ret);
318         }
319     }
320
321     if (nth_affinity_set > nthread_local)
322     {
323         char msg[STRLEN];
324
325         sprintf(msg, "Looks like we have set affinity for more threads than "
326                 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
327         gmx_incons(msg);
328     }
329     else
330     {
331         /* check & warn if some threads failed to set their affinities */
332         if (nth_affinity_set != nthread_local)
333         {
334             char sbuf1[STRLEN], sbuf2[STRLEN];
335
336             /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
337             sbuf1[0] = sbuf2[0] = '\0';
338             /* Only add rank info if we have more than one rank. */
339             if (cr->nnodes > 1)
340             {
341 #ifdef GMX_MPI
342 #ifdef GMX_THREAD_MPI
343                 sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
344 #else           /* GMX_LIB_MPI */
345                 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
346 #endif
347 #endif          /* GMX_MPI */
348             }
349
350             if (nthread_local > 1)
351             {
352                 sprintf(sbuf2, "for %d/%d thread%s ",
353                         nthread_local - nth_affinity_set, nthread_local,
354                         nthread_local > 1 ? "s" : "");
355             }
356
357             md_print_warn(NULL, fplog,
358                           "NOTE: %sAffinity setting %sfailed. This can cause performance degradation.\n"
359                           "      If you think your settings are correct, ask on the gmx-users list.",
360                           sbuf1, sbuf2);
361         }
362     }
363     return;
364 }
365
366 /* Check the process affinity mask and if it is found to be non-zero,
367  * will honor it and disable mdrun internal affinity setting.
368  * Note that this will only work on Linux as we use a GNU feature.
369  */
370 void
371 gmx_check_thread_affinity_set(FILE            *fplog,
372                               const t_commrec *cr,
373                               gmx_hw_opt_t    *hw_opt,
374                               int  gmx_unused  nthreads_hw_avail,
375                               gmx_bool         bAfterOpenmpInit)
376 {
377 #ifdef HAVE_SCHED_AFFINITY
378     cpu_set_t mask_current;
379     int       i, ret, cpu_count, cpu_set;
380     gmx_bool  bAllSet;
381 #endif
382 #ifdef GMX_LIB_MPI
383     gmx_bool  bAllSet_All;
384 #endif
385
386     assert(hw_opt);
387     if (!bAfterOpenmpInit)
388     {
389         /* Check for externally set OpenMP affinity and turn off internal
390          * pinning if any is found. We need to do this check early to tell
391          * thread-MPI whether it should do pinning when spawning threads.
392          * TODO: the above no longer holds, we should move these checks later
393          */
394         if (hw_opt->thread_affinity != threadaffOFF)
395         {
396             char *message;
397             if (!gmx_omp_check_thread_affinity(&message))
398             {
399                 /* TODO: with -pin auto we should only warn when using all cores */
400                 md_print_warn(cr, fplog, "%s", message);
401                 sfree(message);
402                 hw_opt->thread_affinity = threadaffOFF;
403             }
404         }
405
406         /* With thread-MPI this is needed as pinning might get turned off,
407          * which needs to be known before starting thread-MPI.
408          * With thread-MPI hw_opt is processed here on the master rank
409          * and passed to the other ranks later, so we only do this on master.
410          */
411         if (!SIMMASTER(cr))
412         {
413             return;
414         }
415 #ifndef GMX_THREAD_MPI
416         return;
417 #endif
418     }
419
420 #ifdef HAVE_SCHED_GETAFFINITY
421     if (hw_opt->thread_affinity == threadaffOFF)
422     {
423         /* internal affinity setting is off, don't bother checking process affinity */
424         return;
425     }
426
427     CPU_ZERO(&mask_current);
428     if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
429     {
430         /* failed to query affinity mask, will just return */
431         if (debug)
432         {
433             fprintf(debug, "Failed to query affinity mask (error %d)", ret);
434         }
435         return;
436     }
437
438     /* Before proceeding with the actual check, make sure that the number of
439      * detected CPUs is >= the CPUs in the current set.
440      * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
441 #ifdef CPU_COUNT
442     if (nthreads_hw_avail < CPU_COUNT(&mask_current))
443     {
444         if (debug)
445         {
446             fprintf(debug, "%d hardware threads detected, but %d was returned by CPU_COUNT",
447                     nthreads_hw_avail, CPU_COUNT(&mask_current));
448         }
449         return;
450     }
451 #endif /* CPU_COUNT */
452
453     bAllSet = TRUE;
454     for (i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
455     {
456         bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
457     }
458
459 #ifdef GMX_LIB_MPI
460     MPI_Allreduce(&bAllSet, &bAllSet_All, 1, MPI_INT, MPI_LAND, MPI_COMM_WORLD);
461     bAllSet = bAllSet_All;
462 #endif
463
464     if (!bAllSet)
465     {
466         if (hw_opt->thread_affinity == threadaffAUTO)
467         {
468             if (!bAfterOpenmpInit)
469             {
470                 md_print_warn(cr, fplog,
471                               "Non-default thread affinity set, disabling internal thread affinity");
472             }
473             else
474             {
475                 md_print_warn(cr, fplog,
476                               "Non-default thread affinity set probably by the OpenMP library,\n"
477                               "disabling internal thread affinity");
478             }
479             hw_opt->thread_affinity = threadaffOFF;
480         }
481         else
482         {
483             /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
484             if (bAfterOpenmpInit)
485             {
486                 md_print_warn(cr, fplog,
487                               "Overriding thread affinity set outside %s\n",
488                               ShortProgram());
489             }
490         }
491
492         if (debug)
493         {
494             fprintf(debug, "Non-default affinity mask found\n");
495         }
496     }
497     else
498     {
499         if (debug)
500         {
501             fprintf(debug, "Default affinity mask found\n");
502         }
503     }
504 #endif /* HAVE_SCHED_AFFINITY */
505 }