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