Corrected RDTSCP warning text
[alexxy/gromacs.git] / src / gromacs / gmxlib / gmx_detect_hardware.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015, 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 h=ttp://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
39 #include <stdlib.h>
40 #include <assert.h>
41 #include <string.h>
42
43 #ifdef HAVE_UNISTD_H
44 /* For sysconf */
45 #include <unistd.h>
46 #endif
47
48 #include "types/enums.h"
49 #include "types/hw_info.h"
50 #include "types/commrec.h"
51 #include "gmx_fatal.h"
52 #include "gmx_fatal_collective.h"
53 #include "md_logging.h"
54 #include "gmx_cpuid.h"
55 #include "gromacs/utility/smalloc.h"
56 #include "gpu_utils.h"
57 #include "copyrite.h"
58 #include "gmx_detect_hardware.h"
59 #include "main.h"
60 #include "md_logging.h"
61 #include "gromacs/utility/gmxomp.h"
62 #include "gromacs/utility/cstringutil.h"
63
64 #include "thread_mpi/threads.h"
65
66 #ifdef GMX_NATIVE_WINDOWS
67 #include <windows.h>
68 #endif
69
70 #ifdef GMX_GPU
71 const gmx_bool bGPUBinary = TRUE;
72 #else
73 const gmx_bool bGPUBinary = FALSE;
74 #endif
75
76 static const char * invalid_gpuid_hint =
77     "A delimiter-free sequence of valid numeric IDs of available GPUs is expected.";
78
79 /* The globally shared hwinfo structure. */
80 static gmx_hw_info_t      *hwinfo_g;
81 /* A reference counter for the hwinfo structure */
82 static int                 n_hwinfo = 0;
83 /* A lock to protect the hwinfo structure */
84 static tMPI_Thread_mutex_t hw_info_lock = TMPI_THREAD_MUTEX_INITIALIZER;
85
86
87 /* FW decl. */
88 static void limit_num_gpus_used(gmx_gpu_opt_t *gpu_opt, int count);
89 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
90                                     const gmx_gpu_opt_t  *gpu_opt);
91
92 static void sprint_gpus(char *sbuf, const gmx_gpu_info_t *gpu_info)
93 {
94     int      i, ndev;
95     char     stmp[STRLEN];
96
97     ndev = gpu_info->ncuda_dev;
98
99     sbuf[0] = '\0';
100     for (i = 0; i < ndev; i++)
101     {
102         get_gpu_device_info_string(stmp, gpu_info, i);
103         strcat(sbuf, "  ");
104         strcat(sbuf, stmp);
105         if (i < ndev - 1)
106         {
107             strcat(sbuf, "\n");
108         }
109     }
110 }
111
112 static void print_gpu_detection_stats(FILE                 *fplog,
113                                       const gmx_gpu_info_t *gpu_info,
114                                       const t_commrec      *cr)
115 {
116     char onhost[266], stmp[STRLEN];
117     int  ngpu;
118
119     if (!gpu_info->bDetectGPUs)
120     {
121         /* We skipped the detection, so don't print detection stats */
122         return;
123     }
124
125     ngpu = gpu_info->ncuda_dev;
126
127 #if defined GMX_MPI && !defined GMX_THREAD_MPI
128     /* We only print the detection on one, of possibly multiple, nodes */
129     strncpy(onhost, " on host ", 10);
130     gmx_gethostname(onhost+9, 256);
131 #else
132     /* We detect all relevant GPUs */
133     strncpy(onhost, "", 1);
134 #endif
135
136     if (ngpu > 0)
137     {
138         sprint_gpus(stmp, gpu_info);
139         md_print_warn(cr, fplog, "%d GPU%s detected%s:\n%s\n",
140                       ngpu, (ngpu > 1) ? "s" : "", onhost, stmp);
141     }
142     else
143     {
144         md_print_warn(cr, fplog, "No GPUs detected%s\n", onhost);
145     }
146 }
147
148 static void print_gpu_use_stats(FILE                 *fplog,
149                                 const gmx_gpu_info_t *gpu_info,
150                                 const gmx_gpu_opt_t  *gpu_opt,
151                                 const t_commrec      *cr)
152 {
153     char sbuf[STRLEN], stmp[STRLEN];
154     int  i, ngpu_comp, ngpu_use;
155
156     ngpu_comp = gpu_info->ncuda_dev_compatible;
157     ngpu_use  = gpu_opt->ncuda_dev_use;
158
159     /* Issue a note if GPUs are available but not used */
160     if (ngpu_comp > 0 && ngpu_use < 1)
161     {
162         sprintf(sbuf,
163                 "%d compatible GPU%s detected in the system, but none will be used.\n"
164                 "Consider trying GPU acceleration with the Verlet scheme!",
165                 ngpu_comp, (ngpu_comp > 1) ? "s" : "");
166     }
167     else
168     {
169         int ngpu_use_uniq;
170
171         ngpu_use_uniq = gmx_count_gpu_dev_unique(gpu_info, gpu_opt);
172
173         sprintf(sbuf, "%d GPU%s %sselected for this run.\n"
174                 "Mapping of GPU%s to the %d PP rank%s in this node: ",
175                 ngpu_use_uniq, (ngpu_use_uniq > 1) ? "s" : "",
176                 gpu_opt->bUserSet ? "user-" : "auto-",
177                 (ngpu_use > 1) ? "s" : "",
178                 cr->nrank_pp_intranode,
179                 (cr->nrank_pp_intranode > 1) ? "s" : "");
180
181         for (i = 0; i < ngpu_use; i++)
182         {
183             sprintf(stmp, "#%d", get_gpu_device_id(gpu_info, gpu_opt, i));
184             if (i < ngpu_use - 1)
185             {
186                 strcat(stmp, ", ");
187             }
188             strcat(sbuf, stmp);
189         }
190     }
191     md_print_info(cr, fplog, "%s\n\n", sbuf);
192 }
193
194 /* Give a suitable fatal error or warning if the build configuration
195    and runtime CPU do not match. */
196 static void
197 check_use_of_rdtscp_on_this_cpu(FILE                *fplog,
198                                 const t_commrec     *cr,
199                                 const gmx_hw_info_t *hwinfo)
200 {
201     gmx_bool bCpuHasRdtscp, bBinaryUsesRdtscp;
202 #ifdef HAVE_RDTSCP
203     bBinaryUsesRdtscp = TRUE;
204 #else
205     bBinaryUsesRdtscp = FALSE;
206 #endif
207
208     bCpuHasRdtscp = gmx_cpuid_feature(hwinfo->cpuid_info, GMX_CPUID_FEATURE_X86_RDTSCP);
209
210     if (!bCpuHasRdtscp && bBinaryUsesRdtscp)
211     {
212         gmx_fatal(FARGS, "The %s executable was compiled to use the rdtscp CPU instruction. "
213                   "However, this is not supported by the current hardware and continuing would lead to a crash. "
214                   "Please rebuild GROMACS with the GMX_USE_RDTSCP=OFF CMake option.",
215                   ShortProgram());
216     }
217
218     if (bCpuHasRdtscp && !bBinaryUsesRdtscp)
219     {
220         md_print_warn(cr, fplog, "The current CPU can measure timings more accurately than the code in\n"
221                       "%s was configured to use. This might affect your simulation\n"
222                       "speed as accurate timings are needed for load-balancing.\n"
223                       "Please consider rebuilding %s with the GMX_USE_RDTSCP=ON CMake option.\n",
224                       ShortProgram(), ShortProgram());
225     }
226 }
227
228 void gmx_check_hw_runconf_consistency(FILE                *fplog,
229                                       const gmx_hw_info_t *hwinfo,
230                                       const t_commrec     *cr,
231                                       const gmx_hw_opt_t  *hw_opt,
232                                       gmx_bool             bUseGPU)
233 {
234     int      npppn, ntmpi_pp;
235     char     sbuf[STRLEN], th_or_proc[STRLEN], th_or_proc_plural[STRLEN], pernode[STRLEN];
236     gmx_bool btMPI, bMPI, bMaxMpiThreadsSet, bNthreadsAuto, bEmulateGPU;
237
238     assert(hwinfo);
239     assert(cr);
240
241     /* Below we only do consistency checks for PP and GPUs,
242      * this is irrelevant for PME only nodes, so in that case we return
243      * here.
244      */
245     if (!(cr->duty & DUTY_PP))
246     {
247         return;
248     }
249
250     btMPI         = bMPI = FALSE;
251     bNthreadsAuto = FALSE;
252 #if defined(GMX_THREAD_MPI)
253     btMPI         = TRUE;
254     bNthreadsAuto = (hw_opt->nthreads_tmpi < 1);
255 #elif defined(GMX_LIB_MPI)
256     bMPI  = TRUE;
257 #endif
258
259     /* GPU emulation detection is done later, but we need here as well
260      * -- uncool, but there's no elegant workaround */
261     bEmulateGPU       = (getenv("GMX_EMULATE_GPU") != NULL);
262     bMaxMpiThreadsSet = (getenv("GMX_MAX_MPI_THREADS") != NULL);
263
264     /* check the SIMD level mdrun is compiled with against hardware
265        capabilities */
266     /* TODO: Here we assume homogeneous hardware which is not necessarily
267              the case! Might not hurt to add an extra check over MPI. */
268     gmx_cpuid_simd_check(hwinfo->cpuid_info, fplog, SIMMASTER(cr));
269
270     check_use_of_rdtscp_on_this_cpu(fplog, cr, hwinfo);
271
272     /* NOTE: this print is only for and on one physical node */
273     print_gpu_detection_stats(fplog, &hwinfo->gpu_info, cr);
274
275     if (hwinfo->gpu_info.ncuda_dev_compatible > 0)
276     {
277         /* NOTE: this print is only for and on one physical node */
278         print_gpu_use_stats(fplog, &hwinfo->gpu_info, &hw_opt->gpu_opt, cr);
279     }
280
281     /* Need to ensure that we have enough GPUs:
282      * - need one GPU per PP node
283      * - no GPU oversubscription with tMPI
284      * */
285     /* number of PP processes per node */
286     npppn = cr->nrank_pp_intranode;
287
288     pernode[0]           = '\0';
289     th_or_proc_plural[0] = '\0';
290     if (btMPI)
291     {
292         sprintf(th_or_proc, "thread-MPI thread");
293         if (npppn > 1)
294         {
295             sprintf(th_or_proc_plural, "s");
296         }
297     }
298     else if (bMPI)
299     {
300         sprintf(th_or_proc, "MPI process");
301         if (npppn > 1)
302         {
303             sprintf(th_or_proc_plural, "es");
304         }
305         sprintf(pernode, " per node");
306     }
307     else
308     {
309         /* neither MPI nor tMPI */
310         sprintf(th_or_proc, "process");
311     }
312
313     if (bUseGPU && hwinfo->gpu_info.ncuda_dev_compatible > 0 &&
314         !bEmulateGPU)
315     {
316         int  ngpu_comp, ngpu_use;
317         char gpu_comp_plural[2], gpu_use_plural[2];
318
319         ngpu_comp = hwinfo->gpu_info.ncuda_dev_compatible;
320         ngpu_use  = hw_opt->gpu_opt.ncuda_dev_use;
321
322         sprintf(gpu_comp_plural, "%s", (ngpu_comp > 1) ? "s" : "");
323         sprintf(gpu_use_plural,  "%s", (ngpu_use > 1) ? "s" : "");
324
325         /* number of tMPI threads auto-adjusted */
326         if (btMPI && bNthreadsAuto)
327         {
328             if (hw_opt->gpu_opt.bUserSet && npppn < ngpu_use)
329             {
330                 /* The user manually provided more GPUs than threads we
331                    could automatically start. */
332                 gmx_fatal(FARGS,
333                           "%d GPU%s provided, but only %d PP thread-MPI thread%s coud be started.\n"
334                           "%s requires one PP tread-MPI thread per GPU; use fewer GPUs%s.",
335                           ngpu_use, gpu_use_plural,
336                           npppn, th_or_proc_plural,
337                           ShortProgram(), bMaxMpiThreadsSet ? "\nor allow more threads to be used" : "");
338             }
339
340             if (!hw_opt->gpu_opt.bUserSet && npppn < ngpu_comp)
341             {
342                 /* There are more GPUs than tMPI threads; we have
343                    limited the number GPUs used. */
344                 md_print_warn(cr, fplog,
345                               "NOTE: %d GPU%s were detected, but only %d PP thread-MPI thread%s can be started.\n"
346                               "      %s can use one GPU per PP tread-MPI thread, so only %d GPU%s will be used.%s\n",
347                               ngpu_comp, gpu_comp_plural,
348                               npppn, th_or_proc_plural,
349                               ShortProgram(), npppn,
350                               npppn > 1 ? "s" : "",
351                               bMaxMpiThreadsSet ? "\n      Also, you can allow more threads to be used by increasing GMX_MAX_MPI_THREADS" : "");
352             }
353         }
354
355         if (hw_opt->gpu_opt.bUserSet)
356         {
357             if (ngpu_use != npppn)
358             {
359                 gmx_fatal(FARGS,
360                           "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
361                           "%s was started with %d PP %s%s%s, but you provided %d GPU%s.",
362                           th_or_proc, btMPI ? "s" : "es", pernode,
363                           ShortProgram(), npppn, th_or_proc,
364                           th_or_proc_plural, pernode,
365                           ngpu_use, gpu_use_plural);
366             }
367         }
368         else
369         {
370             if (ngpu_comp > npppn)
371             {
372                 md_print_warn(cr, fplog,
373                               "NOTE: potentially sub-optimal launch configuration, %s started with less\n"
374                               "      PP %s%s%s than GPU%s available.\n"
375                               "      Each PP %s can use only one GPU, %d GPU%s%s will be used.\n",
376                               ShortProgram(), th_or_proc,
377                               th_or_proc_plural, pernode, gpu_comp_plural,
378                               th_or_proc, npppn, gpu_use_plural, pernode);
379             }
380
381             if (ngpu_use != npppn)
382             {
383                 /* Avoid duplicate error messages.
384                  * Unfortunately we can only do this at the physical node
385                  * level, since the hardware setup and MPI process count
386                  * might differ between physical nodes.
387                  */
388                 if (cr->rank_pp_intranode == 0)
389                 {
390                     gmx_fatal(FARGS,
391                               "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
392                               "%s was started with %d PP %s%s%s, but only %d GPU%s were detected.",
393                               th_or_proc, btMPI ? "s" : "es", pernode,
394                               ShortProgram(), npppn, th_or_proc,
395                               th_or_proc_plural, pernode,
396                               ngpu_use, gpu_use_plural);
397                 }
398             }
399         }
400
401         {
402             int      same_count;
403
404             same_count = gmx_count_gpu_dev_shared(&hw_opt->gpu_opt);
405
406             if (same_count > 0)
407             {
408                 md_print_info(cr, fplog,
409                               "NOTE: You assigned %s to multiple %s%s.\n",
410                               same_count > 1 ? "GPUs" : "a GPU", th_or_proc, btMPI ? "s" : "es");
411             }
412         }
413     }
414
415 #ifdef GMX_MPI
416     if (PAR(cr))
417     {
418         /* Avoid other ranks to continue after
419            inconsistency */
420         MPI_Barrier(cr->mpi_comm_mygroup);
421     }
422 #endif
423
424 }
425
426 /* Return 0 if none of the GPU (per node) are shared among PP ranks.
427  *
428  * Sharing GPUs among multiple PP ranks is possible when the user passes
429  * GPU IDs. Here we check for sharing and return a non-zero value when
430  * this is detected. Note that the return value represents the number of
431  * PP rank pairs that share a device.
432  */
433 int gmx_count_gpu_dev_shared(const gmx_gpu_opt_t *gpu_opt)
434 {
435     int      same_count    = 0;
436     int      ngpu          = gpu_opt->ncuda_dev_use;
437
438     if (gpu_opt->bUserSet)
439     {
440         int      i, j;
441
442         for (i = 0; i < ngpu - 1; i++)
443         {
444             for (j = i + 1; j < ngpu; j++)
445             {
446                 same_count      += (gpu_opt->cuda_dev_use[i] ==
447                                     gpu_opt->cuda_dev_use[j]);
448             }
449         }
450     }
451
452     return same_count;
453 }
454
455 /* Count and return the number of unique GPUs (per node) selected.
456  *
457  * As sharing GPUs among multiple PP ranks is possible when the user passes
458  * GPU IDs, the number of GPUs user (per node) can be different from the
459  * number of GPU IDs selected.
460  */
461 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
462                                     const gmx_gpu_opt_t  *gpu_opt)
463 {
464     int  i, uniq_count, ngpu;
465     int *uniq_ids;
466
467     assert(gpu_info);
468     assert(gpu_opt);
469
470     ngpu        = gpu_info->ncuda_dev;
471     uniq_count  = 0;
472
473     snew(uniq_ids, ngpu);
474
475     /* Each element in uniq_ids will be set to 0 or 1. The n-th element set
476      * to 1 indicates that the respective GPU was selected to be used. */
477     for (i = 0; i < gpu_opt->ncuda_dev_use; i++)
478     {
479         uniq_ids[get_gpu_device_id(gpu_info, gpu_opt, i)] = 1;
480     }
481     /* Count the devices used. */
482     for (i = 0; i < ngpu; i++)
483     {
484         uniq_count += uniq_ids[i];
485     }
486
487     sfree(uniq_ids);
488
489     return uniq_count;
490 }
491
492
493 /* Return the number of hardware threads supported by the current CPU.
494  * We assume that this is equal with the number of "processors"
495  * reported to be online by the OS at the time of the call. The
496  * definition of "processor" is according to an old POSIX standard.
497  *
498  * Note that the number of hardware threads is generally greater than
499  * the number of cores (e.g. x86 hyper-threading, Power). Managing the
500  * mapping of software threads to hardware threads is managed
501  * elsewhere. */
502 static int get_nthreads_hw_avail(FILE gmx_unused *fplog, const t_commrec gmx_unused *cr)
503 {
504     int ret = 0;
505
506 #if ((defined(WIN32) || defined( _WIN32 ) || defined(WIN64) || defined( _WIN64 )) && !(defined (__CYGWIN__) || defined (__CYGWIN32__)))
507     /* Windows */
508     SYSTEM_INFO sysinfo;
509     GetSystemInfo( &sysinfo );
510     ret = sysinfo.dwNumberOfProcessors;
511 #elif defined HAVE_SYSCONF
512     /* We are probably on Unix.
513      * Now check if we have the argument to use before executing the call
514      */
515 #if defined(_SC_NPROCESSORS_ONLN)
516     ret = sysconf(_SC_NPROCESSORS_ONLN);
517 #elif defined(_SC_NPROC_ONLN)
518     ret = sysconf(_SC_NPROC_ONLN);
519 #elif defined(_SC_NPROCESSORS_CONF)
520     ret = sysconf(_SC_NPROCESSORS_CONF);
521 #elif defined(_SC_NPROC_CONF)
522     ret = sysconf(_SC_NPROC_CONF);
523 #else
524 #warning "No valid sysconf argument value found. Executables will not be able to determine the number of hardware threads: mdrun will use 1 thread by default!"
525 #endif /* End of check for sysconf argument values */
526
527 #else
528     /* Neither windows nor Unix. No fscking idea how many hardware threads we have! */
529     ret = -1;
530 #endif
531
532     if (debug)
533     {
534         fprintf(debug, "Detected %d hardware threads to use.\n", ret);
535     }
536
537 #ifdef GMX_OPENMP
538     if (ret != gmx_omp_get_num_procs())
539     {
540         md_print_warn(cr, fplog,
541                       "Number of hardware threads detected (%d) does not match the number reported by OpenMP (%d).\n"
542                       "Consider setting the launch configuration manually!",
543                       ret, gmx_omp_get_num_procs());
544     }
545 #endif
546
547     return ret;
548 }
549
550 static void gmx_detect_gpus(FILE *fplog, const t_commrec *cr)
551 {
552 #ifdef GMX_LIB_MPI
553     int              rank_world;
554     MPI_Comm         physicalnode_comm;
555 #endif
556     int              rank_local;
557
558     /* Under certain circumstances MPI ranks on the same physical node
559      * can not simultaneously access the same GPU(s). Therefore we run
560      * the detection only on one MPI rank per node and broadcast the info.
561      * Note that with thread-MPI only a single thread runs this code.
562      *
563      * TODO: We should also do CPU hardware detection only once on each
564      * physical node and broadcast it, instead of do it on every MPI rank.
565      */
566 #ifdef GMX_LIB_MPI
567     /* A split of MPI_COMM_WORLD over physical nodes is only required here,
568      * so we create and destroy it locally.
569      */
570     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
571     MPI_Comm_split(MPI_COMM_WORLD, gmx_physicalnode_id_hash(),
572                    rank_world, &physicalnode_comm);
573     MPI_Comm_rank(physicalnode_comm, &rank_local);
574 #else
575     /* Here there should be only one process, check this */
576     assert(cr->nnodes == 1 && cr->sim_nodeid == 0);
577
578     rank_local = 0;
579 #endif
580
581     if (rank_local == 0)
582     {
583         char detection_error[STRLEN] = "", sbuf[STRLEN];
584
585         if (detect_cuda_gpus(&hwinfo_g->gpu_info, detection_error) != 0)
586         {
587             if (detection_error[0] != '\0')
588             {
589                 sprintf(sbuf, ":\n      %s\n", detection_error);
590             }
591             else
592             {
593                 sprintf(sbuf, ".");
594             }
595             md_print_warn(cr, fplog,
596                           "NOTE: Error occurred during GPU detection%s"
597                           "      Can not use GPU acceleration, will fall back to CPU kernels.\n",
598                           sbuf);
599         }
600     }
601
602 #ifdef GMX_LIB_MPI
603     /* Broadcast the GPU info to the other ranks within this node */
604     MPI_Bcast(&hwinfo_g->gpu_info.ncuda_dev, 1, MPI_INT, 0, physicalnode_comm);
605
606     if (hwinfo_g->gpu_info.ncuda_dev > 0)
607     {
608         int cuda_dev_size;
609
610         cuda_dev_size = hwinfo_g->gpu_info.ncuda_dev*sizeof_cuda_dev_info();
611
612         if (rank_local > 0)
613         {
614             hwinfo_g->gpu_info.cuda_dev =
615                 (cuda_dev_info_ptr_t)malloc(cuda_dev_size);
616         }
617         MPI_Bcast(hwinfo_g->gpu_info.cuda_dev, cuda_dev_size, MPI_BYTE,
618                   0, physicalnode_comm);
619         MPI_Bcast(&hwinfo_g->gpu_info.ncuda_dev_compatible, 1, MPI_INT,
620                   0, physicalnode_comm);
621     }
622
623     MPI_Comm_free(&physicalnode_comm);
624 #endif
625 }
626
627 gmx_hw_info_t *gmx_detect_hardware(FILE *fplog, const t_commrec *cr,
628                                    gmx_bool bDetectGPUs)
629 {
630     gmx_hw_info_t   *hw;
631     int              ret;
632
633     /* make sure no one else is doing the same thing */
634     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
635     if (ret != 0)
636     {
637         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
638     }
639
640     /* only initialize the hwinfo structure if it is not already initalized */
641     if (n_hwinfo == 0)
642     {
643         snew(hwinfo_g, 1);
644
645         /* detect CPUID info; no fuss, we don't detect system-wide
646          * -- sloppy, but that's it for now */
647         if (gmx_cpuid_init(&hwinfo_g->cpuid_info) != 0)
648         {
649             gmx_fatal_collective(FARGS, cr, NULL, "CPUID detection failed!");
650         }
651
652         /* detect number of hardware threads */
653         hwinfo_g->nthreads_hw_avail = get_nthreads_hw_avail(fplog, cr);
654
655         /* detect GPUs */
656         hwinfo_g->gpu_info.ncuda_dev            = 0;
657         hwinfo_g->gpu_info.cuda_dev             = NULL;
658         hwinfo_g->gpu_info.ncuda_dev_compatible = 0;
659
660         /* Run the detection if the binary was compiled with GPU support
661          * and we requested detection.
662          */
663         hwinfo_g->gpu_info.bDetectGPUs =
664             (bGPUBinary && bDetectGPUs &&
665              getenv("GMX_DISABLE_GPU_DETECTION") == NULL);
666         if (hwinfo_g->gpu_info.bDetectGPUs)
667         {
668             gmx_detect_gpus(fplog, cr);
669         }
670     }
671     /* increase the reference counter */
672     n_hwinfo++;
673
674     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
675     if (ret != 0)
676     {
677         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
678     }
679
680     return hwinfo_g;
681 }
682
683 void gmx_parse_gpu_ids(gmx_gpu_opt_t *gpu_opt)
684 {
685     char *env;
686
687     if (gpu_opt->gpu_id != NULL && !bGPUBinary)
688     {
689         gmx_fatal(FARGS, "GPU ID string set, but %s was compiled without GPU support!", ShortProgram());
690     }
691
692     env = getenv("GMX_GPU_ID");
693     if (env != NULL && gpu_opt->gpu_id != NULL)
694     {
695         gmx_fatal(FARGS, "GMX_GPU_ID and -gpu_id can not be used at the same time");
696     }
697     if (env == NULL)
698     {
699         env = gpu_opt->gpu_id;
700     }
701
702     /* parse GPU IDs if the user passed any */
703     if (env != NULL)
704     {
705         /* Parse a "plain" GPU ID string which contains a sequence of
706          * digits corresponding to GPU IDs; the order will indicate
707          * the process/tMPI thread - GPU assignment. */
708         parse_digits_from_plain_string(env,
709                                        &gpu_opt->ncuda_dev_use,
710                                        &gpu_opt->cuda_dev_use);
711
712         if (gpu_opt->ncuda_dev_use == 0)
713         {
714             gmx_fatal(FARGS, "Empty GPU ID string encountered.\n%s\n",
715                       invalid_gpuid_hint);
716         }
717
718         gpu_opt->bUserSet = TRUE;
719     }
720 }
721
722 void gmx_select_gpu_ids(FILE *fplog, const t_commrec *cr,
723                         const gmx_gpu_info_t *gpu_info,
724                         gmx_bool bForceUseGPU,
725                         gmx_gpu_opt_t *gpu_opt)
726 {
727     int              i;
728     const char      *env;
729     char             sbuf[STRLEN], stmp[STRLEN];
730
731     /* Bail if binary is not compiled with GPU acceleration, but this is either
732      * explicitly (-nb gpu) or implicitly (gpu ID passed) requested. */
733     if (bForceUseGPU && !bGPUBinary)
734     {
735         gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!", ShortProgram());
736     }
737
738     if (gpu_opt->bUserSet)
739     {
740         /* Check the GPU IDs passed by the user.
741          * (GPU IDs have been parsed by gmx_parse_gpu_ids before)
742          */
743         int *checkres;
744         int  res;
745
746         snew(checkres, gpu_opt->ncuda_dev_use);
747
748         res = check_selected_cuda_gpus(checkres, gpu_info, gpu_opt);
749
750         if (!res)
751         {
752             print_gpu_detection_stats(fplog, gpu_info, cr);
753
754             sprintf(sbuf, "Some of the requested GPUs do not exist, behave strangely, or are not compatible:\n");
755             for (i = 0; i < gpu_opt->ncuda_dev_use; i++)
756             {
757                 if (checkres[i] != egpuCompatible)
758                 {
759                     sprintf(stmp, "    GPU #%d: %s\n",
760                             gpu_opt->cuda_dev_use[i],
761                             gpu_detect_res_str[checkres[i]]);
762                     strcat(sbuf, stmp);
763                 }
764             }
765             gmx_fatal(FARGS, "%s", sbuf);
766         }
767
768         sfree(checkres);
769     }
770     else
771     {
772         pick_compatible_gpus(&hwinfo_g->gpu_info, gpu_opt);
773
774         if (gpu_opt->ncuda_dev_use > cr->nrank_pp_intranode)
775         {
776             /* We picked more GPUs than we can use: limit the number.
777              * We print detailed messages about this later in
778              * gmx_check_hw_runconf_consistency.
779              */
780             limit_num_gpus_used(gpu_opt, cr->nrank_pp_intranode);
781         }
782
783         gpu_opt->bUserSet = FALSE;
784     }
785
786     /* If the user asked for a GPU, check whether we have a GPU */
787     if (bForceUseGPU && gpu_info->ncuda_dev_compatible == 0)
788     {
789         gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
790     }
791 }
792
793 static void limit_num_gpus_used(gmx_gpu_opt_t *gpu_opt, int count)
794 {
795     int ndev_use;
796
797     assert(gpu_opt);
798
799     ndev_use = gpu_opt->ncuda_dev_use;
800
801     if (count > ndev_use)
802     {
803         /* won't increase the # of GPUs */
804         return;
805     }
806
807     if (count < 1)
808     {
809         char sbuf[STRLEN];
810         sprintf(sbuf, "Limiting the number of GPUs to <1 doesn't make sense (detected %d, %d requested)!",
811                 ndev_use, count);
812         gmx_incons(sbuf);
813     }
814
815     /* TODO: improve this implementation: either sort GPUs or remove the weakest here */
816     gpu_opt->ncuda_dev_use = count;
817 }
818
819 void gmx_hardware_info_free(gmx_hw_info_t *hwinfo)
820 {
821     int ret;
822
823     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
824     if (ret != 0)
825     {
826         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
827     }
828
829     /* decrease the reference counter */
830     n_hwinfo--;
831
832
833     if (hwinfo != hwinfo_g)
834     {
835         gmx_incons("hwinfo < hwinfo_g");
836     }
837
838     if (n_hwinfo < 0)
839     {
840         gmx_incons("n_hwinfo < 0");
841     }
842
843     if (n_hwinfo == 0)
844     {
845         gmx_cpuid_done(hwinfo_g->cpuid_info);
846         free_gpu_info(&hwinfo_g->gpu_info);
847         sfree(hwinfo_g);
848     }
849
850     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
851     if (ret != 0)
852     {
853         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
854     }
855 }