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