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