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