New quotes
[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,2015, 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/gmxlib/gpu_utils/gpu_utils.h"
60 #include "gromacs/legacyheaders/copyrite.h"
61 #include "gromacs/legacyheaders/gmx_cpuid.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 #  ifdef GMX_USE_OPENCL
82 const char    *gpu_implementation        = "OpenCL";
83 /* Our current OpenCL implementation only supports using exactly one
84  * GPU per PP rank, so sharing is impossible */
85 const gmx_bool bGpuSharingSupported      = FALSE;
86 /* Our current OpenCL implementation is not known to handle
87  * concurrency correctly (at context creation, JIT compilation, or JIT
88  * cache-management stages). OpenCL runtimes need not support it
89  * either; library MPI segfaults when creating OpenCL contexts;
90  * thread-MPI seems to work but is not yet known to be safe. */
91 const gmx_bool bMultiGpuPerNodeSupported = FALSE;
92 #  else
93 const char    *gpu_implementation        = "CUDA";
94 const gmx_bool bGpuSharingSupported      = TRUE;
95 const gmx_bool bMultiGpuPerNodeSupported = TRUE;
96 #  endif
97 #else
98 const gmx_bool bGPUBinary                = FALSE;
99 const char    *gpu_implementation        = "non-GPU";
100 const gmx_bool bGpuSharingSupported      = FALSE;
101 const gmx_bool bMultiGpuPerNodeSupported = FALSE;
102 #endif
103
104 /* Names of the GPU detection/check results (see e_gpu_detect_res_t in hw_info.h). */
105 const char * const gpu_detect_res_str[egpuNR] =
106 {
107     "compatible", "inexistent", "incompatible", "insane"
108 };
109
110 static const char * invalid_gpuid_hint =
111     "A delimiter-free sequence of valid numeric IDs of available GPUs is expected.";
112
113 /* The globally shared hwinfo structure. */
114 static gmx_hw_info_t      *hwinfo_g;
115 /* A reference counter for the hwinfo structure */
116 static int                 n_hwinfo = 0;
117 /* A lock to protect the hwinfo structure */
118 static tMPI_Thread_mutex_t hw_info_lock = TMPI_THREAD_MUTEX_INITIALIZER;
119
120 #define HOSTNAMELEN 80
121
122 /* FW decl. */
123 static void set_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank);
124 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
125                                     const gmx_gpu_opt_t  *gpu_opt);
126
127 static void sprint_gpus(char *sbuf, const gmx_gpu_info_t *gpu_info)
128 {
129     int      i, ndev;
130     char     stmp[STRLEN];
131
132     ndev = gpu_info->n_dev;
133
134     sbuf[0] = '\0';
135     for (i = 0; i < ndev; i++)
136     {
137         get_gpu_device_info_string(stmp, gpu_info, i);
138         strcat(sbuf, "    ");
139         strcat(sbuf, stmp);
140         if (i < ndev - 1)
141         {
142             strcat(sbuf, "\n");
143         }
144     }
145 }
146
147 static void print_gpu_detection_stats(FILE                 *fplog,
148                                       const gmx_gpu_info_t *gpu_info,
149                                       const t_commrec      *cr)
150 {
151     char onhost[HOSTNAMELEN+10], stmp[STRLEN];
152     int  ngpu;
153
154     if (!gpu_info->bDetectGPUs)
155     {
156         /* We skipped the detection, so don't print detection stats */
157         return;
158     }
159
160     ngpu = gpu_info->n_dev;
161
162 #if defined GMX_MPI && !defined GMX_THREAD_MPI
163     /* We only print the detection on one, of possibly multiple, nodes */
164     strncpy(onhost, " on host ", 10);
165     gmx_gethostname(onhost + 9, HOSTNAMELEN);
166 #else
167     /* We detect all relevant GPUs */
168     strncpy(onhost, "", 1);
169 #endif
170
171     if (ngpu > 0)
172     {
173         sprint_gpus(stmp, gpu_info);
174         md_print_warn(cr, fplog, "%d GPU%s detected%s:\n%s\n",
175                       ngpu, (ngpu > 1) ? "s" : "", onhost, stmp);
176     }
177     else
178     {
179         md_print_warn(cr, fplog, "No GPUs detected%s\n", onhost);
180     }
181 }
182
183 /*! \brief Helper function for reporting GPU usage information
184  * in the mdrun log file
185  *
186  * \param[in] gpu_info       Pointer to per-node GPU info struct
187  * \param[in] gpu_opt        Pointer to per-node GPU options struct
188  * \param[in] numPpRanks     Number of PP ranks per node
189  * \param[in] bPrintHostName Print the hostname in the usage information
190  * \return                   String to write to the log file
191  * \throws                   std::bad_alloc if out of memory */
192 static std::string
193 makeGpuUsageReport(const gmx_gpu_info_t *gpu_info,
194                    const gmx_gpu_opt_t  *gpu_opt,
195                    size_t                numPpRanks,
196                    bool                  bPrintHostName)
197 {
198     int  ngpu_use  = gpu_opt->n_dev_use;
199     int  ngpu_comp = gpu_info->n_dev_compatible;
200     char host[HOSTNAMELEN];
201
202     if (bPrintHostName)
203     {
204         gmx_gethostname(host, HOSTNAMELEN);
205     }
206
207     /* Issue a note if GPUs are available but not used */
208     if (ngpu_comp > 0 && ngpu_use < 1)
209     {
210         return gmx::formatString("%d compatible GPU%s detected in the system, but none will be used.\n"
211                                  "Consider trying GPU acceleration with the Verlet scheme!\n",
212                                  ngpu_comp, (ngpu_comp > 1) ? "s" : "");
213     }
214
215     std::string output;
216     if (!gpu_opt->bUserSet)
217     {
218         // gpu_opt->dev_compatible is only populated during auto-selection
219         std::string gpuIdsString =
220             formatAndJoin(gmx::constArrayRefFromArray(gpu_opt->dev_compatible,
221                                                       gpu_opt->n_dev_compatible),
222                           ",", gmx::StringFormatter("%d"));
223         bool bPluralGpus = gpu_opt->n_dev_compatible > 1;
224
225         if (bPrintHostName)
226         {
227             output += gmx::formatString("On host %s ", host);
228         }
229         output += gmx::formatString("%d compatible GPU%s %s present, with ID%s %s\n",
230                                     gpu_opt->n_dev_compatible,
231                                     bPluralGpus ? "s" : "",
232                                     bPluralGpus ? "are" : "is",
233                                     bPluralGpus ? "s" : "",
234                                     gpuIdsString.c_str());
235     }
236
237     {
238         std::vector<int> gpuIdsInUse;
239         for (int i = 0; i < ngpu_use; i++)
240         {
241             gpuIdsInUse.push_back(get_gpu_device_id(gpu_info, gpu_opt, i));
242         }
243         std::string gpuIdsString =
244             formatAndJoin(gpuIdsInUse, ",", gmx::StringFormatter("%d"));
245         int         numGpusInUse = gmx_count_gpu_dev_unique(gpu_info, gpu_opt);
246         bool        bPluralGpus  = numGpusInUse > 1;
247
248         if (bPrintHostName)
249         {
250             output += gmx::formatString("On host %s ", host);
251         }
252         output += gmx::formatString("%d GPU%s %sselected for this run.\n"
253                                     "Mapping of GPU ID%s to the %d PP rank%s in this node: %s\n",
254                                     numGpusInUse, bPluralGpus ? "s" : "",
255                                     gpu_opt->bUserSet ? "user-" : "auto-",
256                                     bPluralGpus ? "s" : "",
257                                     numPpRanks,
258                                     (numPpRanks > 1) ? "s" : "",
259                                     gpuIdsString.c_str());
260     }
261
262     return output;
263 }
264
265 /* Give a suitable fatal error or warning if the build configuration
266    and runtime CPU do not match. */
267 static void
268 check_use_of_rdtscp_on_this_cpu(FILE                *fplog,
269                                 const t_commrec     *cr,
270                                 const gmx_hw_info_t *hwinfo)
271 {
272     gmx_bool bCpuHasRdtscp, bBinaryUsesRdtscp;
273 #ifdef HAVE_RDTSCP
274     bBinaryUsesRdtscp = TRUE;
275 #else
276     bBinaryUsesRdtscp = FALSE;
277 #endif
278
279     bCpuHasRdtscp = gmx_cpuid_feature(hwinfo->cpuid_info, GMX_CPUID_FEATURE_X86_RDTSCP);
280
281     if (!bCpuHasRdtscp && bBinaryUsesRdtscp)
282     {
283         gmx_fatal(FARGS, "The %s executable was compiled to use the rdtscp CPU instruction. "
284                   "However, this is not supported by the current hardware and continuing would lead to a crash. "
285                   "Please rebuild GROMACS with the GMX_USE_RDTSCP=OFF CMake option.",
286                   ShortProgram());
287     }
288
289     if (bCpuHasRdtscp && !bBinaryUsesRdtscp)
290     {
291         md_print_warn(cr, fplog, "The current CPU can measure timings more accurately than the code in\n"
292                       "%s was configured to use. This might affect your simulation\n"
293                       "speed as accurate timings are needed for load-balancing.\n"
294                       "Please consider rebuilding %s with the GMX_USE_RDTSCP=ON CMake option.\n",
295                       ShortProgram(), ShortProgram());
296     }
297 }
298
299 void gmx_check_hw_runconf_consistency(FILE                *fplog,
300                                       const gmx_hw_info_t *hwinfo,
301                                       const t_commrec     *cr,
302                                       const gmx_hw_opt_t  *hw_opt,
303                                       gmx_bool             bUseGPU)
304 {
305     int      npppn;
306     char     th_or_proc[STRLEN], th_or_proc_plural[STRLEN], pernode[STRLEN];
307     gmx_bool btMPI, bMPI, bNthreadsAuto, bEmulateGPU;
308
309     assert(hwinfo);
310     assert(cr);
311
312     /* Below we only do consistency checks for PP and GPUs,
313      * this is irrelevant for PME only nodes, so in that case we return
314      * here.
315      */
316     if (!(cr->duty & DUTY_PP))
317     {
318         return;
319     }
320
321 #if defined(GMX_THREAD_MPI)
322     bMPI          = FALSE;
323     btMPI         = TRUE;
324     bNthreadsAuto = (hw_opt->nthreads_tmpi < 1);
325 #elif defined(GMX_LIB_MPI)
326     bMPI          = TRUE;
327     btMPI         = FALSE;
328     bNthreadsAuto = FALSE;
329 #else
330     bMPI          = FALSE;
331     btMPI         = FALSE;
332     bNthreadsAuto = FALSE;
333 #endif
334
335     /* GPU emulation detection is done later, but we need here as well
336      * -- uncool, but there's no elegant workaround */
337     bEmulateGPU       = (getenv("GMX_EMULATE_GPU") != NULL);
338
339     if (hwinfo->gpu_info.n_dev_compatible > 0)
340     {
341         std::string gpuUseageReport;
342         try
343         {
344             gpuUseageReport = makeGpuUsageReport(&hwinfo->gpu_info,
345                                                  &hw_opt->gpu_opt,
346                                                  cr->nrank_pp_intranode,
347                                                  bMPI && cr->nnodes > 1);
348         }
349         GMX_CATCH_ALL_AND_EXIT_WITH_FATAL_ERROR;
350
351         /* NOTE: this print is only for and on one physical node */
352         md_print_info(cr, fplog, "%s\n", gpuUseageReport.c_str());
353     }
354
355     /* Need to ensure that we have enough GPUs:
356      * - need one GPU per PP node
357      * - no GPU oversubscription with tMPI
358      * */
359     /* number of PP processes per node */
360     npppn = cr->nrank_pp_intranode;
361
362     pernode[0]           = '\0';
363     th_or_proc_plural[0] = '\0';
364     if (btMPI)
365     {
366         sprintf(th_or_proc, "thread-MPI thread");
367         if (npppn > 1)
368         {
369             sprintf(th_or_proc_plural, "s");
370         }
371     }
372     else if (bMPI)
373     {
374         sprintf(th_or_proc, "MPI process");
375         if (npppn > 1)
376         {
377             sprintf(th_or_proc_plural, "es");
378         }
379         sprintf(pernode, " per node");
380     }
381     else
382     {
383         /* neither MPI nor tMPI */
384         sprintf(th_or_proc, "process");
385     }
386
387     if (bUseGPU && hwinfo->gpu_info.n_dev_compatible > 0 &&
388         !bEmulateGPU)
389     {
390         int  ngpu_comp, ngpu_use;
391         char gpu_comp_plural[2], gpu_use_plural[2];
392
393         ngpu_comp = hwinfo->gpu_info.n_dev_compatible;
394         ngpu_use  = hw_opt->gpu_opt.n_dev_use;
395
396         sprintf(gpu_comp_plural, "%s", (ngpu_comp > 1) ? "s" : "");
397         sprintf(gpu_use_plural,  "%s", (ngpu_use > 1) ? "s" : "");
398
399         /* number of tMPI threads auto-adjusted */
400         if (btMPI && bNthreadsAuto)
401         {
402             if (hw_opt->gpu_opt.bUserSet && npppn < ngpu_use)
403             {
404                 /* The user manually provided more GPUs than threads we
405                    could automatically start. */
406                 gmx_fatal(FARGS,
407                           "%d GPU%s provided, but only %d PP thread-MPI thread%s coud be started.\n"
408                           "%s requires one PP tread-MPI thread per GPU; use fewer GPUs.",
409                           ngpu_use, gpu_use_plural,
410                           npppn, th_or_proc_plural,
411                           ShortProgram());
412             }
413
414             if (!hw_opt->gpu_opt.bUserSet && npppn < ngpu_comp)
415             {
416                 /* There are more GPUs than tMPI threads; we have
417                    limited the number GPUs used. */
418                 md_print_warn(cr, fplog,
419                               "NOTE: %d GPU%s were detected, but only %d PP thread-MPI thread%s can be started.\n"
420                               "      %s can use one GPU per PP tread-MPI thread, so only %d GPU%s will be used.\n",
421                               ngpu_comp, gpu_comp_plural,
422                               npppn, th_or_proc_plural,
423                               ShortProgram(), npppn,
424                               npppn > 1 ? "s" : "");
425             }
426         }
427
428         if (hw_opt->gpu_opt.bUserSet)
429         {
430             if (ngpu_use != npppn)
431             {
432                 gmx_fatal(FARGS,
433                           "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
434                           "%s was started with %d PP %s%s%s, but you provided %d GPU%s.",
435                           th_or_proc, btMPI ? "s" : "es", pernode,
436                           ShortProgram(), npppn, th_or_proc,
437                           th_or_proc_plural, pernode,
438                           ngpu_use, gpu_use_plural);
439             }
440         }
441         else
442         {
443             if (ngpu_comp > npppn)
444             {
445                 md_print_warn(cr, fplog,
446                               "NOTE: potentially sub-optimal launch configuration, %s started with less\n"
447                               "      PP %s%s%s than GPU%s available.\n"
448                               "      Each PP %s can use only one GPU, %d GPU%s%s will be used.\n",
449                               ShortProgram(), th_or_proc,
450                               th_or_proc_plural, pernode, gpu_comp_plural,
451                               th_or_proc, npppn, gpu_use_plural, pernode);
452             }
453
454             if (ngpu_use != npppn)
455             {
456                 /* Avoid duplicate error messages.
457                  * Unfortunately we can only do this at the physical node
458                  * level, since the hardware setup and MPI process count
459                  * might differ between physical nodes.
460                  */
461                 if (cr->rank_pp_intranode == 0)
462                 {
463                     gmx_fatal(FARGS,
464                               "Incorrect launch configuration: mismatching number of PP %s%s and GPUs%s.\n"
465                               "%s was started with %d PP %s%s%s, but only %d GPU%s were detected.",
466                               th_or_proc, btMPI ? "s" : "es", pernode,
467                               ShortProgram(), npppn, th_or_proc,
468                               th_or_proc_plural, pernode,
469                               ngpu_use, gpu_use_plural);
470                 }
471             }
472         }
473
474         {
475             int      same_count;
476
477             same_count = gmx_count_gpu_dev_shared(&hw_opt->gpu_opt);
478
479             if (same_count > 0)
480             {
481                 md_print_info(cr, fplog,
482                               "NOTE: You assigned %s to multiple %s%s.\n",
483                               same_count > 1 ? "GPUs" : "a GPU", th_or_proc, btMPI ? "s" : "es");
484             }
485         }
486     }
487
488 #ifdef GMX_MPI
489     if (PAR(cr))
490     {
491         /* Avoid other ranks to continue after
492            inconsistency */
493         MPI_Barrier(cr->mpi_comm_mygroup);
494     }
495 #endif
496
497 }
498
499 /* Return 0 if none of the GPU (per node) are shared among PP ranks.
500  *
501  * Sharing GPUs among multiple PP ranks is possible when the user passes
502  * GPU IDs. Here we check for sharing and return a non-zero value when
503  * this is detected. Note that the return value represents the number of
504  * PP rank pairs that share a device.
505  */
506 int gmx_count_gpu_dev_shared(const gmx_gpu_opt_t *gpu_opt)
507 {
508     int      same_count    = 0;
509     int      ngpu          = gpu_opt->n_dev_use;
510
511     if (gpu_opt->bUserSet)
512     {
513         int      i, j;
514
515         for (i = 0; i < ngpu - 1; i++)
516         {
517             for (j = i + 1; j < ngpu; j++)
518             {
519                 same_count      += (gpu_opt->dev_use[i] ==
520                                     gpu_opt->dev_use[j]);
521             }
522         }
523     }
524
525     return same_count;
526 }
527
528 /* Count and return the number of unique GPUs (per node) selected.
529  *
530  * As sharing GPUs among multiple PP ranks is possible when the user passes
531  * GPU IDs, the number of GPUs user (per node) can be different from the
532  * number of GPU IDs selected.
533  */
534 static int gmx_count_gpu_dev_unique(const gmx_gpu_info_t *gpu_info,
535                                     const gmx_gpu_opt_t  *gpu_opt)
536 {
537     int  i, uniq_count, ngpu;
538     int *uniq_ids;
539
540     assert(gpu_info);
541     assert(gpu_opt);
542
543     ngpu = gpu_info->n_dev;
544
545     uniq_count  = 0;
546
547     snew(uniq_ids, ngpu);
548
549     /* Each element in uniq_ids will be set to 0 or 1. The n-th element set
550      * to 1 indicates that the respective GPU was selected to be used. */
551     for (i = 0; i < gpu_opt->n_dev_use; i++)
552     {
553         int device_id;
554
555         device_id           = bGpuSharingSupported ? get_gpu_device_id(gpu_info, gpu_opt, i) : i;
556         uniq_ids[device_id] = 1;
557     }
558     /* Count the devices used. */
559     for (i = 0; i < ngpu; i++)
560     {
561         uniq_count += uniq_ids[i];
562     }
563
564     sfree(uniq_ids);
565
566     return uniq_count;
567 }
568
569 static int get_ncores(gmx_cpuid_t cpuid)
570 {
571     int        nprocessors, npackages, ncores_per_package, nhwthreads_per_core;
572     const int *package_id, *core_id, *hwthread_id, *locality_order;
573     int        rc;
574
575     rc = gmx_cpuid_topology(cpuid,
576                             &nprocessors, &npackages,
577                             &ncores_per_package, &nhwthreads_per_core,
578                             &package_id, &core_id,
579                             &hwthread_id, &locality_order);
580
581     if (rc == 0)
582     {
583         return npackages*ncores_per_package;
584     }
585     else
586     {
587         /* We don't have cpuid topology info, return 0 core count */
588         return 0;
589     }
590 }
591
592 /* Return the number of hardware threads supported by the current CPU.
593  * We assume that this is equal with the number of "processors"
594  * reported to be online by the OS at the time of the call. The
595  * definition of "processor" is according to an old POSIX standard.
596  *
597  * Note that the number of hardware threads is generally greater than
598  * the number of cores (e.g. x86 hyper-threading, Power). Managing the
599  * mapping of software threads to hardware threads is managed
600  * elsewhere. */
601 static int get_nthreads_hw_avail(FILE gmx_unused *fplog, const t_commrec gmx_unused *cr)
602 {
603     int ret = 0;
604
605 #if ((defined(WIN32) || defined( _WIN32 ) || defined(WIN64) || defined( _WIN64 )) && !(defined (__CYGWIN__) || defined (__CYGWIN32__)))
606     /* Windows */
607     SYSTEM_INFO sysinfo;
608     GetSystemInfo( &sysinfo );
609     ret = sysinfo.dwNumberOfProcessors;
610 #elif defined HAVE_SYSCONF
611     /* We are probably on Unix.
612      * Now check if we have the argument to use before executing the call
613      */
614 #if defined(_SC_NPROCESSORS_ONLN)
615     ret = sysconf(_SC_NPROCESSORS_ONLN);
616 #elif defined(_SC_NPROC_ONLN)
617     ret = sysconf(_SC_NPROC_ONLN);
618 #elif defined(_SC_NPROCESSORS_CONF)
619     ret = sysconf(_SC_NPROCESSORS_CONF);
620 #elif defined(_SC_NPROC_CONF)
621     ret = sysconf(_SC_NPROC_CONF);
622 #else
623 #warning "No valid sysconf argument value found. Executables will not be able to determine the number of logical cores: mdrun will use 1 thread by default!"
624 #endif /* End of check for sysconf argument values */
625
626 #else
627     /* Neither windows nor Unix. No fscking idea how many hardware threads we have! */
628     ret = -1;
629 #endif
630
631     if (debug)
632     {
633         fprintf(debug, "Detected %d hardware threads to use.\n", ret);
634     }
635
636 #ifdef GMX_OPENMP
637     if (ret != gmx_omp_get_num_procs())
638     {
639         md_print_warn(cr, fplog,
640                       "Number of logical cores detected (%d) does not match the number reported by OpenMP (%d).\n"
641                       "Consider setting the launch configuration manually!",
642                       ret, gmx_omp_get_num_procs());
643     }
644 #endif
645
646     return ret;
647 }
648
649 static void gmx_detect_gpus(FILE *fplog, const t_commrec *cr)
650 {
651 #ifdef GMX_LIB_MPI
652     int              rank_world;
653     MPI_Comm         physicalnode_comm;
654 #endif
655     int              rank_local;
656
657     /* Under certain circumstances MPI ranks on the same physical node
658      * can not simultaneously access the same GPU(s). Therefore we run
659      * the detection only on one MPI rank per node and broadcast the info.
660      * Note that with thread-MPI only a single thread runs this code.
661      *
662      * TODO: We should also do CPU hardware detection only once on each
663      * physical node and broadcast it, instead of do it on every MPI rank.
664      */
665 #ifdef GMX_LIB_MPI
666     /* A split of MPI_COMM_WORLD over physical nodes is only required here,
667      * so we create and destroy it locally.
668      */
669     MPI_Comm_rank(MPI_COMM_WORLD, &rank_world);
670     MPI_Comm_split(MPI_COMM_WORLD, gmx_physicalnode_id_hash(),
671                    rank_world, &physicalnode_comm);
672     MPI_Comm_rank(physicalnode_comm, &rank_local);
673 #else
674     /* Here there should be only one process, check this */
675     assert(cr->nnodes == 1 && cr->sim_nodeid == 0);
676
677     rank_local = 0;
678 #endif
679
680     if (rank_local == 0)
681     {
682         char detection_error[STRLEN] = "", sbuf[STRLEN];
683
684         if (detect_gpus(&hwinfo_g->gpu_info, detection_error) != 0)
685         {
686             if (detection_error[0] != '\0')
687             {
688                 sprintf(sbuf, ":\n      %s\n", detection_error);
689             }
690             else
691             {
692                 sprintf(sbuf, ".");
693             }
694             md_print_warn(cr, fplog,
695                           "NOTE: Error occurred during GPU detection%s"
696                           "      Can not use GPU acceleration, will fall back to CPU kernels.\n",
697                           sbuf);
698         }
699     }
700
701 #ifdef GMX_LIB_MPI
702     /* Broadcast the GPU info to the other ranks within this node */
703     MPI_Bcast(&hwinfo_g->gpu_info.n_dev, 1, MPI_INT, 0, physicalnode_comm);
704
705     if (hwinfo_g->gpu_info.n_dev > 0)
706     {
707         int dev_size;
708
709         dev_size = hwinfo_g->gpu_info.n_dev*sizeof_gpu_dev_info();
710
711         if (rank_local > 0)
712         {
713             hwinfo_g->gpu_info.gpu_dev =
714                 (struct gmx_device_info_t *)malloc(dev_size);
715         }
716         MPI_Bcast(hwinfo_g->gpu_info.gpu_dev, dev_size, MPI_BYTE,
717                   0, physicalnode_comm);
718         MPI_Bcast(&hwinfo_g->gpu_info.n_dev_compatible, 1, MPI_INT,
719                   0, physicalnode_comm);
720     }
721
722     MPI_Comm_free(&physicalnode_comm);
723 #endif
724 }
725
726 static void gmx_collect_hardware_mpi()
727 {
728 #ifdef GMX_LIB_MPI
729     int  rank_id;
730     int  nrank, rank, ncore, nhwthread, ngpu, i;
731     int  gpu_hash;
732     int *buf, *all;
733
734     rank_id   = gmx_physicalnode_id_hash();
735     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
736     MPI_Comm_size(MPI_COMM_WORLD, &nrank);
737     ncore     = hwinfo_g->ncore;
738     nhwthread = hwinfo_g->nthreads_hw_avail;
739     ngpu      = hwinfo_g->gpu_info.n_dev_compatible;
740     /* Create a unique hash of the GPU type(s) in this node */
741     gpu_hash  = 0;
742     /* Here it might be better to only loop over the compatible GPU, but we
743      * don't have that information available and it would also require
744      * removing the device ID from the device info string.
745      */
746     for (i = 0; i < hwinfo_g->gpu_info.n_dev; i++)
747     {
748         char stmp[STRLEN];
749
750         /* Since the device ID is incorporated in the hash, the order of
751          * the GPUs affects the hash. Also two identical GPUs won't give
752          * a gpu_hash of zero after XORing.
753          */
754         get_gpu_device_info_string(stmp, &hwinfo_g->gpu_info, i);
755         gpu_hash ^= gmx_string_fullhash_func(stmp, gmx_string_hash_init);
756     }
757
758     snew(buf, nrank);
759     snew(all, nrank);
760     buf[rank] = rank_id;
761
762     MPI_Allreduce(buf, all, nrank, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
763
764     gmx_bool bFound;
765     int      nnode0, ncore0, nhwthread0, ngpu0, r;
766
767     bFound     = FALSE;
768     ncore0     = 0;
769     nnode0     = 0;
770     nhwthread0 = 0;
771     ngpu0      = 0;
772     for (r = 0; r < nrank; r++)
773     {
774         if (all[r] == rank_id)
775         {
776             if (!bFound && r == rank)
777             {
778                 /* We are the first rank in this physical node */
779                 nnode0     = 1;
780                 ncore0     = ncore;
781                 nhwthread0 = nhwthread;
782                 ngpu0      = ngpu;
783             }
784             bFound = TRUE;
785         }
786     }
787
788     sfree(buf);
789     sfree(all);
790
791     int sum[4], maxmin[10];
792
793     {
794         int buf[4];
795
796         /* Sum values from only intra-rank 0 so we get the sum over all nodes */
797         buf[0] = nnode0;
798         buf[1] = ncore0;
799         buf[2] = nhwthread0;
800         buf[3] = ngpu0;
801
802         MPI_Allreduce(buf, sum, 4, MPI_INT, MPI_SUM, MPI_COMM_WORLD);
803     }
804
805     {
806         int buf[10];
807
808         /* Store + and - values for all ranks,
809          * so we can get max+min with one MPI call.
810          */
811         buf[0] = ncore;
812         buf[1] = nhwthread;
813         buf[2] = ngpu;
814         buf[3] = gmx_cpuid_simd_suggest(hwinfo_g->cpuid_info);
815         buf[4] = gpu_hash;
816         buf[5] = -buf[0];
817         buf[6] = -buf[1];
818         buf[7] = -buf[2];
819         buf[8] = -buf[3];
820         buf[9] = -buf[4];
821
822         MPI_Allreduce(buf, maxmin, 10, MPI_INT, MPI_MAX, MPI_COMM_WORLD);
823     }
824
825     hwinfo_g->nphysicalnode       = sum[0];
826     hwinfo_g->ncore_tot           = sum[1];
827     hwinfo_g->ncore_min           = -maxmin[5];
828     hwinfo_g->ncore_max           = maxmin[0];
829     hwinfo_g->nhwthread_tot       = sum[2];
830     hwinfo_g->nhwthread_min       = -maxmin[6];
831     hwinfo_g->nhwthread_max       = maxmin[1];
832     hwinfo_g->ngpu_compatible_tot = sum[3];
833     hwinfo_g->ngpu_compatible_min = -maxmin[7];
834     hwinfo_g->ngpu_compatible_max = maxmin[2];
835     hwinfo_g->simd_suggest_min    = static_cast<enum gmx_cpuid_simd>(-maxmin[8]);
836     hwinfo_g->simd_suggest_max    = static_cast<enum gmx_cpuid_simd>(maxmin[3]);
837     hwinfo_g->bIdenticalGPUs      = (maxmin[4] == -maxmin[9]);
838 #else
839     /* All ranks use the same pointer, protect it with a mutex */
840     tMPI_Thread_mutex_lock(&hw_info_lock);
841     hwinfo_g->nphysicalnode       = 1;
842     hwinfo_g->ncore_tot           = hwinfo_g->ncore;
843     hwinfo_g->ncore_min           = hwinfo_g->ncore;
844     hwinfo_g->ncore_max           = hwinfo_g->ncore;
845     hwinfo_g->nhwthread_tot       = hwinfo_g->nthreads_hw_avail;
846     hwinfo_g->nhwthread_min       = hwinfo_g->nthreads_hw_avail;
847     hwinfo_g->nhwthread_max       = hwinfo_g->nthreads_hw_avail;
848     hwinfo_g->ngpu_compatible_tot = hwinfo_g->gpu_info.n_dev_compatible;
849     hwinfo_g->ngpu_compatible_min = hwinfo_g->gpu_info.n_dev_compatible;
850     hwinfo_g->ngpu_compatible_max = hwinfo_g->gpu_info.n_dev_compatible;
851     hwinfo_g->simd_suggest_min    = gmx_cpuid_simd_suggest(hwinfo_g->cpuid_info);
852     hwinfo_g->simd_suggest_max    = gmx_cpuid_simd_suggest(hwinfo_g->cpuid_info);
853     hwinfo_g->bIdenticalGPUs      = TRUE;
854     tMPI_Thread_mutex_unlock(&hw_info_lock);
855 #endif
856 }
857
858 gmx_hw_info_t *gmx_detect_hardware(FILE *fplog, const t_commrec *cr,
859                                    gmx_bool bDetectGPUs)
860 {
861     int ret;
862
863     /* make sure no one else is doing the same thing */
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     /* only initialize the hwinfo structure if it is not already initalized */
871     if (n_hwinfo == 0)
872     {
873         snew(hwinfo_g, 1);
874
875         /* detect CPUID info; no fuss, we don't detect system-wide
876          * -- sloppy, but that's it for now */
877         if (gmx_cpuid_init(&hwinfo_g->cpuid_info) != 0)
878         {
879             gmx_fatal_collective(FARGS, cr, NULL, "CPUID detection failed!");
880         }
881
882         /* get the number of cores, will be 0 when not detected */
883         hwinfo_g->ncore             = get_ncores(hwinfo_g->cpuid_info);
884
885         /* detect number of hardware threads */
886         hwinfo_g->nthreads_hw_avail = get_nthreads_hw_avail(fplog, cr);
887
888         /* detect GPUs */
889         hwinfo_g->gpu_info.n_dev            = 0;
890         hwinfo_g->gpu_info.n_dev_compatible = 0;
891         hwinfo_g->gpu_info.gpu_dev          = NULL;
892
893         /* Run the detection if the binary was compiled with GPU support
894          * and we requested detection.
895          */
896         hwinfo_g->gpu_info.bDetectGPUs =
897             (bGPUBinary && bDetectGPUs &&
898              getenv("GMX_DISABLE_GPU_DETECTION") == NULL);
899         if (hwinfo_g->gpu_info.bDetectGPUs)
900         {
901             gmx_detect_gpus(fplog, cr);
902         }
903     }
904     /* increase the reference counter */
905     n_hwinfo++;
906
907     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
908     if (ret != 0)
909     {
910         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
911     }
912
913     gmx_collect_hardware_mpi();
914
915     return hwinfo_g;
916 }
917
918 static std::string detected_hardware_string(const gmx_hw_info_t *hwinfo,
919                                             bool                 bFullCpuInfo)
920 {
921     std::string s;
922
923     s  = gmx::formatString("\n");
924     s += gmx::formatString("Running on %d node%s with total",
925                            hwinfo->nphysicalnode,
926                            hwinfo->nphysicalnode == 1 ? "" : "s");
927     if (hwinfo->ncore_tot > 0)
928     {
929         s += gmx::formatString(" %d cores,", hwinfo->ncore_tot);
930     }
931     s += gmx::formatString(" %d logical cores", hwinfo->nhwthread_tot);
932     if (hwinfo->gpu_info.bDetectGPUs)
933     {
934         s += gmx::formatString(", %d compatible GPU%s",
935                                hwinfo->ngpu_compatible_tot,
936                                hwinfo->ngpu_compatible_tot == 1 ? "" : "s");
937     }
938     else if (bGPUBinary)
939     {
940         s += gmx::formatString(" (GPU detection deactivated)");
941     }
942     s += gmx::formatString("\n");
943
944     if (hwinfo->nphysicalnode > 1)
945     {
946         /* Print per node hardware feature counts */
947         if (hwinfo->ncore_max > 0)
948         {
949             s += gmx::formatString("  Cores per node:           %2d", hwinfo->ncore_min);
950             if (hwinfo->ncore_max > hwinfo->ncore_min)
951             {
952                 s += gmx::formatString(" - %2d", hwinfo->ncore_max);
953             }
954             s += gmx::formatString("\n");
955         }
956         s += gmx::formatString("  Logical cores per node:   %2d", hwinfo->nhwthread_min);
957         if (hwinfo->nhwthread_max > hwinfo->nhwthread_min)
958         {
959             s += gmx::formatString(" - %2d", hwinfo->nhwthread_max);
960         }
961         s += gmx::formatString("\n");
962         if (bGPUBinary)
963         {
964             s += gmx::formatString("  Compatible GPUs per node: %2d",
965                                    hwinfo->ngpu_compatible_min);
966             if (hwinfo->ngpu_compatible_max > hwinfo->ngpu_compatible_min)
967             {
968                 s += gmx::formatString(" - %2d", hwinfo->ngpu_compatible_max);
969             }
970             s += gmx::formatString("\n");
971             if (hwinfo->ngpu_compatible_tot > 0)
972             {
973                 if (hwinfo->bIdenticalGPUs)
974                 {
975                     s += gmx::formatString("  All nodes have identical type(s) of GPUs\n");
976                 }
977                 else
978                 {
979                     /* This message will also appear with identical GPU types
980                      * when at least one node has no GPU.
981                      */
982                     s += gmx::formatString("  Different nodes have different type(s) and/or order of GPUs\n");
983                 }
984             }
985         }
986     }
987
988 #ifdef GMX_LIB_MPI
989     char host[HOSTNAMELEN];
990     int  rank;
991
992     gmx_gethostname(host, HOSTNAMELEN);
993     MPI_Comm_rank(MPI_COMM_WORLD, &rank);
994
995     s += gmx::formatString("Hardware detected on host %s (the node of MPI rank %d):\n",
996                            host, rank);
997 #else
998     s += gmx::formatString("Hardware detected:\n");
999 #endif
1000     s += gmx::formatString("  CPU info:\n");
1001     if (bFullCpuInfo)
1002     {
1003         char buf[1024];
1004
1005         gmx_cpuid_formatstring(hwinfo->cpuid_info, buf, 1023);
1006         buf[1023] = '\0';
1007
1008         s += gmx::formatString("%s", buf);
1009     }
1010     else
1011     {
1012         s += gmx::formatString("    Vendor: %s\n",
1013                                gmx_cpuid_vendor_string[gmx_cpuid_vendor(hwinfo->cpuid_info)]);
1014         s += gmx::formatString("    Brand:  %s\n",
1015                                gmx_cpuid_brand(hwinfo->cpuid_info));
1016     }
1017     s += gmx::formatString("    SIMD instructions most likely to fit this hardware: %s",
1018                            gmx_cpuid_simd_string[hwinfo->simd_suggest_min]);
1019     if (hwinfo->simd_suggest_max > hwinfo->simd_suggest_min)
1020     {
1021         s += gmx::formatString(" - %s",
1022                                gmx_cpuid_simd_string[hwinfo->simd_suggest_max]);
1023     }
1024     s += gmx::formatString("\n");
1025     s += gmx::formatString("    SIMD instructions selected at GROMACS compile time: %s\n",
1026                            gmx_cpuid_simd_string[gmx_compiled_simd()]);
1027     if (bGPUBinary && (hwinfo->ngpu_compatible_tot > 0 ||
1028                        hwinfo->gpu_info.n_dev > 0))
1029     {
1030         s += gmx::formatString("  GPU info:\n");
1031         s += gmx::formatString("    Number of GPUs detected: %d\n",
1032                                hwinfo->gpu_info.n_dev);
1033         if (hwinfo->gpu_info.n_dev > 0)
1034         {
1035             char buf[STRLEN];
1036
1037             sprint_gpus(buf, &hwinfo->gpu_info);
1038             s += gmx::formatString("%s\n", buf);
1039         }
1040     }
1041
1042     return s;
1043 }
1044
1045 void gmx_print_detected_hardware(FILE *fplog, const t_commrec *cr,
1046                                  const gmx_hw_info_t *hwinfo)
1047 {
1048     if (fplog != NULL)
1049     {
1050         std::string detected;
1051
1052         detected = detected_hardware_string(hwinfo, TRUE);
1053
1054         fprintf(fplog, "%s\n", detected.c_str());
1055     }
1056
1057     if (MULTIMASTER(cr))
1058     {
1059         std::string detected;
1060
1061         detected = detected_hardware_string(hwinfo, FALSE);
1062
1063         fprintf(stderr, "%s\n", detected.c_str());
1064     }
1065
1066     /* Check the compiled SIMD instruction set against that of the node
1067      * with the lowest SIMD level support.
1068      */
1069     gmx_cpuid_simd_check(hwinfo->simd_suggest_min, fplog, MULTIMASTER(cr));
1070
1071     /* For RDTSCP we only check on our local node and skip the MPI reduction */
1072     check_use_of_rdtscp_on_this_cpu(fplog, cr, hwinfo);
1073 }
1074
1075 //! \brief Return if any GPU ID (e.g in a user-supplied string) is repeated
1076 static gmx_bool anyGpuIdIsRepeated(const gmx_gpu_opt_t *gpu_opt)
1077 {
1078     /* Loop over IDs in the string */
1079     for (int i = 0; i < gpu_opt->n_dev_use - 1; ++i)
1080     {
1081         /* Look for the ID in location i in the following part of the
1082            string */
1083         for (int j = i + 1; j < gpu_opt->n_dev_use; ++j)
1084         {
1085             if (gpu_opt->dev_use[i] == gpu_opt->dev_use[j])
1086             {
1087                 /* Same ID found in locations i and j */
1088                 return TRUE;
1089             }
1090         }
1091     }
1092
1093     return FALSE;
1094 }
1095
1096 void gmx_parse_gpu_ids(gmx_gpu_opt_t *gpu_opt)
1097 {
1098     char *env;
1099
1100     if (gpu_opt->gpu_id != NULL && !bGPUBinary)
1101     {
1102         gmx_fatal(FARGS, "GPU ID string set, but %s was compiled without GPU support!", ShortProgram());
1103     }
1104
1105     env = getenv("GMX_GPU_ID");
1106     if (env != NULL && gpu_opt->gpu_id != NULL)
1107     {
1108         gmx_fatal(FARGS, "GMX_GPU_ID and -gpu_id can not be used at the same time");
1109     }
1110     if (env == NULL)
1111     {
1112         env = gpu_opt->gpu_id;
1113     }
1114
1115     /* parse GPU IDs if the user passed any */
1116     if (env != NULL)
1117     {
1118         /* Parse a "plain" GPU ID string which contains a sequence of
1119          * digits corresponding to GPU IDs; the order will indicate
1120          * the process/tMPI thread - GPU assignment. */
1121         parse_digits_from_plain_string(env,
1122                                        &gpu_opt->n_dev_use,
1123                                        &gpu_opt->dev_use);
1124         if (!bMultiGpuPerNodeSupported && 1 < gpu_opt->n_dev_use)
1125         {
1126             gmx_fatal(FARGS, "The %s implementation only supports using exactly one PP rank per node", gpu_implementation);
1127         }
1128         if (!bGpuSharingSupported && anyGpuIdIsRepeated(gpu_opt))
1129         {
1130             gmx_fatal(FARGS, "The %s implementation only supports using exactly one PP rank per GPU", gpu_implementation);
1131         }
1132         if (gpu_opt->n_dev_use == 0)
1133         {
1134             gmx_fatal(FARGS, "Empty GPU ID string encountered.\n%s\n",
1135                       invalid_gpuid_hint);
1136         }
1137
1138         gpu_opt->bUserSet = TRUE;
1139     }
1140 }
1141
1142 void gmx_select_gpu_ids(FILE *fplog, const t_commrec *cr,
1143                         const gmx_gpu_info_t *gpu_info,
1144                         gmx_bool bForceUseGPU,
1145                         gmx_gpu_opt_t *gpu_opt)
1146 {
1147     int              i;
1148     char             sbuf[STRLEN], stmp[STRLEN];
1149
1150     /* Bail if binary is not compiled with GPU acceleration, but this is either
1151      * explicitly (-nb gpu) or implicitly (gpu ID passed) requested. */
1152     if (bForceUseGPU && !bGPUBinary)
1153     {
1154         gmx_fatal(FARGS, "GPU acceleration requested, but %s was compiled without GPU support!", ShortProgram());
1155     }
1156
1157     if (!(cr->duty & DUTY_PP))
1158     {
1159         /* Our rank is not doing PP, we don't use a GPU */
1160         return;
1161     }
1162
1163     if (gpu_opt->bUserSet)
1164     {
1165         /* Check the GPU IDs passed by the user.
1166          * (GPU IDs have been parsed by gmx_parse_gpu_ids before)
1167          */
1168         int *checkres;
1169         int  res;
1170
1171         snew(checkres, gpu_opt->n_dev_use);
1172
1173         res = check_selected_gpus(checkres, gpu_info, gpu_opt);
1174
1175         if (!res)
1176         {
1177             print_gpu_detection_stats(fplog, gpu_info, cr);
1178
1179             sprintf(sbuf, "Some of the requested GPUs do not exist, behave strangely, or are not compatible:\n");
1180             for (i = 0; i < gpu_opt->n_dev_use; i++)
1181             {
1182                 if (checkres[i] != egpuCompatible)
1183                 {
1184                     sprintf(stmp, "    GPU #%d: %s\n",
1185                             gpu_opt->dev_use[i],
1186                             gpu_detect_res_str[checkres[i]]);
1187                     strcat(sbuf, stmp);
1188                 }
1189             }
1190             gmx_fatal(FARGS, "%s", sbuf);
1191         }
1192
1193         sfree(checkres);
1194     }
1195     else if (getenv("GMX_EMULATE_GPU") == NULL)
1196     {
1197         pick_compatible_gpus(&hwinfo_g->gpu_info, gpu_opt);
1198         set_gpu_ids(gpu_opt, cr->nrank_pp_intranode, cr->rank_pp_intranode);
1199     }
1200
1201     /* If the user asked for a GPU, check whether we have a GPU */
1202     if (bForceUseGPU && gpu_info->n_dev_compatible == 0)
1203     {
1204         gmx_fatal(FARGS, "GPU acceleration requested, but no compatible GPUs were detected.");
1205     }
1206 }
1207
1208 /* Select the GPUs we will use. This is an operation local to each physical
1209  * node. If we have less MPI ranks than GPUs, we will waste some GPUs.
1210  * nrank and rank are the rank count and id for PP processes in our node.
1211  */
1212 static void set_gpu_ids(gmx_gpu_opt_t *gpu_opt, int nrank, int rank)
1213 {
1214     GMX_RELEASE_ASSERT(gpu_opt, "Invalid gpu_opt pointer passed");
1215     GMX_RELEASE_ASSERT(nrank >= 1,
1216                        gmx::formatString("Invalid limit (%d) for the number of GPUs (detected %d compatible GPUs)",
1217                                          rank, gpu_opt->n_dev_compatible).c_str());
1218
1219     if (gpu_opt->n_dev_compatible == 0)
1220     {
1221         char host[HOSTNAMELEN];
1222
1223         gmx_gethostname(host, HOSTNAMELEN);
1224         gmx_fatal(FARGS, "A GPU was requested on host %s, but no compatible GPUs were detected. All nodes with PP ranks need to have GPUs. If you intended to use GPU acceleration in a parallel run, you can either avoid using the nodes that don't have GPUs or place PME ranks on these nodes.", host);
1225     }
1226
1227     int nshare;
1228
1229     nshare = 1;
1230     if (nrank > gpu_opt->n_dev_compatible)
1231     {
1232         if (nrank % gpu_opt->n_dev_compatible == 0)
1233         {
1234             nshare = bGpuSharingSupported ? nrank/gpu_opt->n_dev_compatible : 1;
1235         }
1236         else
1237         {
1238             if (rank == 0)
1239             {
1240                 gmx_fatal(FARGS, "The number of MPI ranks (%d) in a physical node is not a multiple of the number of GPUs (%d). Select a different number of MPI ranks or use the -gpu_id option to manually specify the GPU to be used.",
1241                           nrank, gpu_opt->n_dev_compatible);
1242             }
1243
1244 #ifdef GMX_MPI
1245             /* We use a global barrier to prevent ranks from continuing with
1246              * an invalid setup.
1247              */
1248             MPI_Barrier(MPI_COMM_WORLD);
1249 #endif
1250         }
1251     }
1252
1253     /* Here we will waste GPUs when nrank < gpu_opt->n_dev_compatible */
1254     gpu_opt->n_dev_use = std::min(gpu_opt->n_dev_compatible*nshare, nrank);
1255     if (!bMultiGpuPerNodeSupported)
1256     {
1257         gpu_opt->n_dev_use = std::min(gpu_opt->n_dev_use, 1);
1258     }
1259     snew(gpu_opt->dev_use, gpu_opt->n_dev_use);
1260     for (int i = 0; i != gpu_opt->n_dev_use; ++i)
1261     {
1262         /* TODO: improve this implementation: either sort GPUs or remove the weakest here */
1263         gpu_opt->dev_use[i] = gpu_opt->dev_compatible[i/nshare];
1264     }
1265 }
1266
1267 void gmx_hardware_info_free(gmx_hw_info_t *hwinfo)
1268 {
1269     int ret;
1270
1271     ret = tMPI_Thread_mutex_lock(&hw_info_lock);
1272     if (ret != 0)
1273     {
1274         gmx_fatal(FARGS, "Error locking hwinfo mutex: %s", strerror(errno));
1275     }
1276
1277     /* decrease the reference counter */
1278     n_hwinfo--;
1279
1280
1281     if (hwinfo != hwinfo_g)
1282     {
1283         gmx_incons("hwinfo < hwinfo_g");
1284     }
1285
1286     if (n_hwinfo < 0)
1287     {
1288         gmx_incons("n_hwinfo < 0");
1289     }
1290
1291     if (n_hwinfo == 0)
1292     {
1293         gmx_cpuid_done(hwinfo_g->cpuid_info);
1294         free_gpu_info(&hwinfo_g->gpu_info);
1295         sfree(hwinfo_g);
1296     }
1297
1298     ret = tMPI_Thread_mutex_unlock(&hw_info_lock);
1299     if (ret != 0)
1300     {
1301         gmx_fatal(FARGS, "Error unlocking hwinfo mutex: %s", strerror(errno));
1302     }
1303 }