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