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