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