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