Improve Gromacs portability
[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)
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 "typedefs.h"
53 #include "types/commrec.h"
54 #include "types/hw_info.h"
55 #include "copyrite.h"
56 #include "gmx_cpuid.h"
57 #include "gmx_omp_nthreads.h"
58 #include "md_logging.h"
59 #include "gmx_thread_affinity.h"
60
61 #include "gmx_fatal.h"
62 #include "gromacs/utility/gmxomp.h"
63
64 static int
65 get_thread_affinity_layout(FILE *fplog,
66                            const t_commrec *cr,
67                            const gmx_hw_info_t * hwinfo,
68                            int nthreads,
69                            int pin_offset, int * pin_stride,
70                            const int **locality_order)
71 {
72     int         nhwthreads, npkg, ncores, nhwthreads_per_core, rc;
73     const int * pkg_id;
74     const int * core_id;
75     const int * hwthread_id;
76     gmx_bool    bPickPinStride;
77
78     if (pin_offset < 0)
79     {
80         gmx_fatal(FARGS, "Negative thread pinning offset requested");
81     }
82     if (*pin_stride < 0)
83     {
84         gmx_fatal(FARGS, "Negative thread pinning stride requested");
85     }
86
87     rc = gmx_cpuid_topology(hwinfo->cpuid_info, &nhwthreads, &npkg, &ncores,
88                             &nhwthreads_per_core,
89                             &pkg_id, &core_id, &hwthread_id, locality_order);
90
91     if (rc != 0)
92     {
93         /* topology information not available or invalid, ignore it */
94         nhwthreads      = hwinfo->nthreads_hw_avail;
95         *locality_order = NULL;
96
97         if (nhwthreads <= 0)
98         {
99             /* We don't know anything about the hardware, don't pin */
100             md_print_warn(cr, fplog,
101                           "NOTE: We don't know how many logical cores we have, will not pin threads");
102
103             return -1;
104         }
105     }
106
107     if (nthreads > nhwthreads)
108     {
109         /* We are oversubscribing, don't pin */
110         md_print_warn(NULL, fplog,
111                       "WARNING: Oversubscribing the CPU, will not pin threads");
112
113         return -1;
114     }
115
116     if (pin_offset + nthreads > nhwthreads)
117     {
118         /* We are oversubscribing, don't pin */
119         md_print_warn(NULL, fplog,
120                       "WARNING: The requested pin offset is too large for the available logical cores,\n"
121                       "         will not pin threads");
122
123         return -1;
124     }
125
126
127     /* do we need to choose the pinning stride? */
128     bPickPinStride = (*pin_stride == 0);
129
130     if (bPickPinStride)
131     {
132         if (rc == 0 && pin_offset + nthreads*nhwthreads_per_core <= nhwthreads)
133         {
134             /* Put one thread on each physical core */
135             *pin_stride = nhwthreads_per_core;
136         }
137         else
138         {
139             /* We don't know if we have SMT, and if we do, we don't know
140              * if hw threads in the same physical core are consecutive.
141              * Without SMT the pinning layout should not matter too much.
142              * so we assume a consecutive layout and maximally spread out"
143              * the threads at equal threads per core.
144              * Note that IBM is the major non-x86 case with cpuid support
145              * and probably threads are already pinned by the queuing system,
146              * so we wouldn't end up here in the first place.
147              */
148             *pin_stride = (nhwthreads - pin_offset)/nthreads;
149         }
150     }
151     else
152     {
153         /* Check the placement of the thread with the largest index to make sure
154          * that the offset & stride doesn't cause pinning beyond the last hardware thread. */
155         if (pin_offset + (nthreads-1)*(*pin_stride) >= nhwthreads)
156         {
157             /* We are oversubscribing, don't pin */
158             md_print_warn(NULL, fplog,
159                           "WARNING: The requested pinning stride is too large for the available logical cores,\n"
160                           "         will not pin threads");
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 doesn't 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 #ifndef __APPLE__
212         md_print_warn(NULL, fplog,
213                       "Can not set thread affinities on the current platform. On NUMA systems this\n"
214                       "can cause performance degradation. If you think your platform should support\n"
215                       "setting affinities, contact the GROMACS developers.");
216 #endif  /* __APPLE__ */
217         return;
218     }
219
220     /* threads on this MPI process or TMPI thread */
221     if (cr->duty & DUTY_PP)
222     {
223         nthread_local = gmx_omp_nthreads_get(emntNonbonded);
224     }
225     else
226     {
227         nthread_local = gmx_omp_nthreads_get(emntPME);
228     }
229
230     /* map the current process to cores */
231     thread0_id_node = 0;
232     nthread_node    = nthread_local;
233 #ifdef GMX_MPI
234     if (PAR(cr) || MULTISIM(cr))
235     {
236         /* We need to determine a scan of the thread counts in this
237          * compute node.
238          */
239         MPI_Comm comm_intra;
240
241         MPI_Comm_split(MPI_COMM_WORLD,
242                        gmx_physicalnode_id_hash(), cr->rank_intranode,
243                        &comm_intra);
244         MPI_Scan(&nthread_local, &thread0_id_node, 1, MPI_INT, MPI_SUM, comm_intra);
245         /* MPI_Scan is inclusive, but here we need exclusive */
246         thread0_id_node -= nthread_local;
247         /* Get the total number of threads on this physical node */
248         MPI_Allreduce(&nthread_local, &nthread_node, 1, MPI_INT, MPI_SUM, comm_intra);
249         MPI_Comm_free(&comm_intra);
250     }
251 #endif
252
253     if (hw_opt->thread_affinity == threadaffAUTO &&
254         nthread_node != hwinfo->nthreads_hw_avail)
255     {
256         if (nthread_node > 1 && nthread_node < hwinfo->nthreads_hw_avail)
257         {
258             md_print_warn(cr, fplog,
259                           "NOTE: The number of threads is not equal to the number of (logical) cores\n"
260                           "      and the -pin option is set to auto: will not pin thread to cores.\n"
261                           "      This can lead to significant performance degradation.\n"
262                           "      Consider using -pin on (and -pinoffset in case you run multiple jobs).\n");
263         }
264
265         return;
266     }
267
268     offset = 0;
269     if (hw_opt->core_pinning_offset != 0)
270     {
271         offset = hw_opt->core_pinning_offset;
272         md_print_info(cr, fplog, "Applying core pinning offset %d\n", offset);
273     }
274
275     rc = get_thread_affinity_layout(fplog, cr, hwinfo,
276                                     nthread_node,
277                                     offset, &hw_opt->core_pinning_stride,
278                                     &locality_order);
279
280     if (rc != 0)
281     {
282         /* Incompatible layout, don't pin, warning was already issued */
283         return;
284     }
285
286     /* Set the per-thread affinity. In order to be able to check the success
287      * of affinity settings, we will set nth_affinity_set to 1 on threads
288      * where the affinity setting succeded and to 0 where it failed.
289      * Reducing these 0/1 values over the threads will give the total number
290      * of threads on which we succeeded.
291      */
292     nth_affinity_set = 0;
293 #pragma omp parallel num_threads(nthread_local) reduction(+:nth_affinity_set)
294     {
295         int      thread_id, thread_id_node;
296         int      index, core;
297         gmx_bool setaffinity_ret;
298
299         thread_id      = gmx_omp_get_thread_num();
300         thread_id_node = thread0_id_node + thread_id;
301         index          = offset + thread_id_node*hw_opt->core_pinning_stride;
302         if (locality_order != NULL)
303         {
304             core = locality_order[index];
305         }
306         else
307         {
308             core = index;
309         }
310
311         setaffinity_ret = tMPI_Thread_setaffinity_single(tMPI_Thread_self(), core);
312
313         /* store the per-thread success-values of the setaffinity */
314         nth_affinity_set = (setaffinity_ret == 0);
315
316         if (debug)
317         {
318             fprintf(debug, "On rank %2d, thread %2d, index %2d, core %2d the affinity setting returned %d\n",
319                     cr->nodeid, gmx_omp_get_thread_num(), index, core, setaffinity_ret);
320         }
321     }
322
323     if (nth_affinity_set > nthread_local)
324     {
325         char msg[STRLEN];
326
327         sprintf(msg, "Looks like we have set affinity for more threads than "
328                 "we have (%d > %d)!\n", nth_affinity_set, nthread_local);
329         gmx_incons(msg);
330     }
331     else
332     {
333         /* check & warn if some threads failed to set their affinities */
334         if (nth_affinity_set != nthread_local)
335         {
336             char sbuf1[STRLEN], sbuf2[STRLEN];
337
338             /* sbuf1 contains rank info, while sbuf2 OpenMP thread info */
339             sbuf1[0] = sbuf2[0] = '\0';
340             /* Only add rank info if we have more than one rank. */
341             if (cr->nnodes > 1)
342             {
343 #ifdef GMX_MPI
344 #ifdef GMX_THREAD_MPI
345                 sprintf(sbuf1, "In tMPI thread #%d: ", cr->nodeid);
346 #else           /* GMX_LIB_MPI */
347                 sprintf(sbuf1, "In MPI process #%d: ", cr->nodeid);
348 #endif
349 #endif          /* GMX_MPI */
350             }
351
352             if (nthread_local > 1)
353             {
354                 sprintf(sbuf2, "for %d/%d thread%s ",
355                         nthread_local - nth_affinity_set, nthread_local,
356                         nthread_local > 1 ? "s" : "");
357             }
358
359             md_print_warn(NULL, fplog,
360                           "WARNING: %sAffinity setting %sfailed.\n"
361                           "         This can cause performance degradation! If you think your setting are\n"
362                           "         correct, contact the GROMACS developers.",
363                           sbuf1, sbuf2);
364         }
365     }
366     return;
367 }
368
369 /* Check the process affinity mask and if it is found to be non-zero,
370  * will honor it and disable mdrun internal affinity setting.
371  * Note that this will only work on Linux as we use a GNU feature.
372  */
373 void
374 gmx_check_thread_affinity_set(FILE            gmx_unused *fplog,
375                               const t_commrec gmx_unused *cr,
376                               gmx_hw_opt_t    gmx_unused *hw_opt,
377                               int             gmx_unused  nthreads_hw_avail,
378                               gmx_bool        gmx_unused  bAfterOpenmpInit)
379 {
380 #ifdef HAVE_SCHED_AFFINITY
381     cpu_set_t mask_current;
382     int       i, ret, cpu_count, cpu_set;
383     gmx_bool  bAllSet;
384
385     assert(hw_opt);
386     if (hw_opt->thread_affinity == threadaffOFF)
387     {
388         /* internal affinity setting is off, don't bother checking process affinity */
389         return;
390     }
391
392     CPU_ZERO(&mask_current);
393     if ((ret = sched_getaffinity(0, sizeof(cpu_set_t), &mask_current)) != 0)
394     {
395         /* failed to query affinity mask, will just return */
396         if (debug)
397         {
398             fprintf(debug, "Failed to query affinity mask (error %d)", ret);
399         }
400         return;
401     }
402
403     /* Before proceeding with the actual check, make sure that the number of
404      * detected CPUs is >= the CPUs in the current set.
405      * We need to check for CPU_COUNT as it was added only in glibc 2.6. */
406 #ifdef CPU_COUNT
407     if (nthreads_hw_avail < CPU_COUNT(&mask_current))
408     {
409         if (debug)
410         {
411             fprintf(debug, "%d hardware threads detected, but %d was returned by CPU_COUNT",
412                     nthreads_hw_avail, CPU_COUNT(&mask_current));
413         }
414         return;
415     }
416 #endif /* CPU_COUNT */
417
418     bAllSet = TRUE;
419     for (i = 0; (i < nthreads_hw_avail && i < CPU_SETSIZE); i++)
420     {
421         bAllSet = bAllSet && (CPU_ISSET(i, &mask_current) != 0);
422     }
423
424     if (!bAllSet)
425     {
426         if (hw_opt->thread_affinity == threadaffAUTO)
427         {
428             if (!bAfterOpenmpInit)
429             {
430                 md_print_warn(cr, fplog,
431                               "Non-default thread affinity set, disabling internal thread affinity");
432             }
433             else
434             {
435                 md_print_warn(cr, fplog,
436                               "Non-default thread affinity set probably by the OpenMP library,\n"
437                               "disabling internal thread affinity");
438             }
439             hw_opt->thread_affinity = threadaffOFF;
440         }
441         else
442         {
443             /* Only warn once, at the last check (bAfterOpenmpInit==TRUE) */
444             if (bAfterOpenmpInit)
445             {
446                 md_print_warn(cr, fplog,
447                               "Overriding thread affinity set outside %s\n",
448                               ShortProgram());
449             }
450         }
451
452         if (debug)
453         {
454             fprintf(debug, "Non-default affinity mask found\n");
455         }
456     }
457     else
458     {
459         if (debug)
460         {
461             fprintf(debug, "Default affinity mask found\n");
462         }
463     }
464 #endif /* HAVE_SCHED_AFFINITY */
465 }