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