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