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