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