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