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