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