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