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