Update mdrun message to reflect that Cuda CC>=3.5 is supported.
[alexxy/gromacs.git] / src / gromacs / hardware / hardwaretopology.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2016 by the GROMACS development team.
5  * Copyright (c) 2017,2018,2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36
37 /*! \internal \file
38  * \brief
39  * Implements gmx::HardwareTopology.
40  *
41  * \author Erik Lindahl <erik.lindahl@gmail.com>
42  * \ingroup module_hardware
43  */
44
45 #include "gmxpre.h"
46
47 #include "hardwaretopology.h"
48
49 #include "config.h"
50
51 #include <cstdio>
52
53 #include <algorithm>
54 #include <functional>
55 #include <limits>
56 #include <utility>
57 #include <vector>
58
59 #if GMX_USE_HWLOC
60 #    include <hwloc.h>
61 #endif
62
63 #include "gromacs/hardware/cpuinfo.h"
64 #include "gromacs/utility/gmxassert.h"
65
66 #ifdef HAVE_UNISTD_H
67 #    include <unistd.h> // sysconf()
68 #endif
69 #if GMX_NATIVE_WINDOWS
70 #    include <windows.h> // GetSystemInfo()
71 #endif
72
73 //! Convenience macro to help us avoid ifdefs each time we use sysconf
74 #if !defined(_SC_NPROCESSORS_ONLN) && defined(_SC_NPROC_ONLN)
75 #    define _SC_NPROCESSORS_ONLN _SC_NPROC_ONLN
76 #endif
77
78 namespace gmx
79 {
80
81 namespace
82 {
83
84 /*****************************************************************************
85  *                                                                           *
86  *   Utility functions for extracting hardware topology from CpuInfo object  *
87  *                                                                           *
88  *****************************************************************************/
89
90 /*! \brief Initialize machine data from basic information in cpuinfo
91  *
92  *  \param  machine      Machine tree structure where information will be assigned
93  *                       if the cpuinfo object contains topology information.
94  *  \param  supportLevel If topology information is available in CpuInfo,
95  *                       this will be updated to reflect the amount of
96  *                       information written to the machine structure.
97  */
98 void parseCpuInfo(HardwareTopology::Machine* machine, HardwareTopology::SupportLevel* supportLevel)
99 {
100     CpuInfo cpuInfo(CpuInfo::detect());
101
102     if (!cpuInfo.logicalProcessors().empty())
103     {
104         int nSockets   = 0;
105         int nCores     = 0;
106         int nHwThreads = 0;
107
108         // Copy the logical processor information from cpuinfo
109         for (const auto& l : cpuInfo.logicalProcessors())
110         {
111             machine->logicalProcessors.push_back(
112                     { l.socketRankInMachine, l.coreRankInSocket, l.hwThreadRankInCore, -1 });
113             nSockets   = std::max(nSockets, l.socketRankInMachine);
114             nCores     = std::max(nCores, l.coreRankInSocket);
115             nHwThreads = std::max(nHwThreads, l.hwThreadRankInCore);
116         }
117
118         // Fill info form sockets/cores/hwthreads
119         int socketId   = 0;
120         int coreId     = 0;
121         int hwThreadId = 0;
122
123         machine->sockets.resize(nSockets + 1);
124         for (auto& s : machine->sockets)
125         {
126             s.id = socketId++;
127             s.cores.resize(nCores + 1);
128             for (auto& c : s.cores)
129             {
130                 c.id         = coreId++;
131                 c.numaNodeId = -1; // No numa information
132                 c.hwThreads.resize(nHwThreads + 1);
133                 for (auto& t : c.hwThreads)
134                 {
135                     t.id                 = hwThreadId++;
136                     t.logicalProcessorId = -1; // set as unassigned for now
137                 }
138             }
139         }
140
141         // Fill the logical processor id in the right place
142         for (std::size_t i = 0; i < machine->logicalProcessors.size(); i++)
143         {
144             const HardwareTopology::LogicalProcessor& l = machine->logicalProcessors[i];
145             machine->sockets[l.socketRankInMachine]
146                     .cores[l.coreRankInSocket]
147                     .hwThreads[l.hwThreadRankInCore]
148                     .logicalProcessorId = static_cast<int>(i);
149         }
150         machine->logicalProcessorCount = machine->logicalProcessors.size();
151         *supportLevel                  = HardwareTopology::SupportLevel::Basic;
152     }
153     else
154     {
155         *supportLevel = HardwareTopology::SupportLevel::None;
156     }
157 }
158
159 #if GMX_USE_HWLOC
160
161 #    if HWLOC_API_VERSION < 0x00010b00
162 #        define HWLOC_OBJ_PACKAGE HWLOC_OBJ_SOCKET
163 #        define HWLOC_OBJ_NUMANODE HWLOC_OBJ_NODE
164 #    endif
165
166 // Preprocessor variable for if hwloc api is version 1.x.x or 2.x.x
167 #    if HWLOC_API_VERSION >= 0x00020000
168 #        define GMX_HWLOC_API_VERSION_IS_2XX 1
169 #        if GMX_HWLOC_API_VERSION < 0x00020000
170 #            error "HWLOC library major version set during configuration is 1, but currently using version 2 headers"
171 #        endif
172 #    else
173 #        define GMX_HWLOC_API_VERSION_IS_2XX 0
174 #        if GMX_HWLOC_API_VERSION >= 0x00020000
175 #            error "HWLOC library major version set during configuration is 2, but currently using version 1 headers"
176 #        endif
177 #    endif
178
179 /*****************************************************************************
180  *                                                                           *
181  *   Utility functions for extracting hardware topology from hwloc library   *
182  *                                                                           *
183  *****************************************************************************/
184
185 // Compatibility function for accessing hwloc_obj_t object memory with different API versions of hwloc
186 std::size_t getHwLocObjectMemory(const hwloc_obj* obj)
187 {
188 #    if GMX_HWLOC_API_VERSION_IS_2XX
189     return obj->total_memory;
190 #    else
191     return obj->memory.total_memory;
192 #    endif
193 }
194
195 /*! \brief Return vector of all descendants of a given type in hwloc topology
196  *
197  *  \param topo  hwloc topology handle that has been initialized and loaded
198  *  \param obj   Non-null hwloc object.
199  *  \param type  hwloc object type to find. The routine will only search
200  *               on levels below obj.
201  *
202  *  \return vector containing all the objects of given type that are
203  *          descendants of the provided object. If no objects of this type
204  *          were found, the vector will be empty.
205  */
206 std::vector<const hwloc_obj*> getHwLocDescendantsByType(const hwloc_topology*  topo,
207                                                         const hwloc_obj*       obj,
208                                                         const hwloc_obj_type_t type)
209 {
210     GMX_RELEASE_ASSERT(obj, "NULL hwloc object provided to getHwLocDescendantsByType()");
211
212     std::vector<const hwloc_obj*> v;
213
214     if (obj->type == type)
215     {
216         v.push_back(obj);
217     }
218     // Go through children; if this object has no children obj->arity is 0,
219     // and we'll return an empty vector.
220     hwloc_obj_t tempNode = nullptr;
221     while ((tempNode = hwloc_get_next_child(
222                     const_cast<hwloc_topology_t>(topo), const_cast<hwloc_obj_t>(obj), tempNode))
223            != nullptr)
224     {
225         std::vector<const hwloc_obj*> v2 = getHwLocDescendantsByType(topo, tempNode, type);
226         v.insert(v.end(), v2.begin(), v2.end());
227     }
228     return v;
229 }
230
231 /*! \brief Read information about sockets, cores and threads from hwloc topology
232  *
233  *  \param topo    hwloc topology handle that has been initialized and loaded
234  *  \param machine Pointer to the machine structure in the HardwareTopology
235  *                 class, where the tree of sockets/cores/threads will be written.
236  *
237  *  \return If all the data is found
238  */
239 bool parseHwLocSocketsCoresThreads(hwloc_topology_t topo, HardwareTopology::Machine* machine)
240 {
241     const hwloc_obj*              root         = hwloc_get_root_obj(topo);
242     std::vector<const hwloc_obj*> hwlocSockets = getHwLocDescendantsByType(topo, root, HWLOC_OBJ_PACKAGE);
243
244     machine->logicalProcessorCount = hwloc_get_nbobjs_by_type(topo, HWLOC_OBJ_PU);
245     machine->logicalProcessors.resize(machine->logicalProcessorCount);
246     machine->sockets.resize(hwlocSockets.size());
247
248     bool topologyOk = !hwlocSockets.empty(); // Fail if we have no sockets in machine
249
250     for (std::size_t i = 0; i < hwlocSockets.size() && topologyOk; i++)
251     {
252         // Assign information about this socket
253         machine->sockets[i].id = hwlocSockets[i]->logical_index;
254
255         // Get children (cores)
256         std::vector<const hwloc_obj*> hwlocCores =
257                 getHwLocDescendantsByType(topo, hwlocSockets[i], HWLOC_OBJ_CORE);
258         machine->sockets[i].cores.resize(hwlocCores.size());
259
260         topologyOk = topologyOk && !hwlocCores.empty(); // Fail if we have no cores in socket
261
262         // Loop over child cores
263         for (std::size_t j = 0; j < hwlocCores.size() && topologyOk; j++)
264         {
265             // Assign information about this core
266             machine->sockets[i].cores[j].id         = hwlocCores[j]->logical_index;
267             machine->sockets[i].cores[j].numaNodeId = -1;
268
269             // Get children (hwthreads)
270             std::vector<const hwloc_obj*> hwlocPUs =
271                     getHwLocDescendantsByType(topo, hwlocCores[j], HWLOC_OBJ_PU);
272             machine->sockets[i].cores[j].hwThreads.resize(hwlocPUs.size());
273
274             topologyOk = topologyOk && !hwlocPUs.empty(); // Fail if we have no hwthreads in core
275
276             // Loop over child hwthreads
277             for (std::size_t k = 0; k < hwlocPUs.size() && topologyOk; k++)
278             {
279                 // Assign information about this hwthread
280                 std::size_t logicalProcessorId               = hwlocPUs[k]->os_index;
281                 machine->sockets[i].cores[j].hwThreads[k].id = hwlocPUs[k]->logical_index;
282                 machine->sockets[i].cores[j].hwThreads[k].logicalProcessorId = logicalProcessorId;
283
284                 if (logicalProcessorId < machine->logicalProcessors.size())
285                 {
286                     // Cross-assign data for this hwthread to the logicalprocess vector
287                     machine->logicalProcessors[logicalProcessorId].socketRankInMachine =
288                             static_cast<int>(i);
289                     machine->logicalProcessors[logicalProcessorId].coreRankInSocket =
290                             static_cast<int>(j);
291                     machine->logicalProcessors[logicalProcessorId].hwThreadRankInCore =
292                             static_cast<int>(k);
293                     machine->logicalProcessors[logicalProcessorId].numaNodeId = -1;
294                 }
295                 else
296                 {
297                     topologyOk = false;
298                 }
299             }
300         }
301     }
302
303     if (!topologyOk)
304     {
305         machine->logicalProcessors.clear();
306         machine->sockets.clear();
307     }
308     return topologyOk;
309 }
310
311 /*! \brief Read cache information from hwloc topology
312  *
313  *  \param topo    hwloc topology handle that has been initialized and loaded
314  *  \param machine Pointer to the machine structure in the HardwareTopology
315  *                 class, where cache data will be filled.
316  *
317  *  \return If any cache data is found
318  */
319 bool parseHwLocCache(hwloc_topology_t topo, HardwareTopology::Machine* machine)
320 {
321     // Parse caches up to L5
322     for (int cachelevel : { 1, 2, 3, 4, 5 })
323     {
324         int depth = hwloc_get_cache_type_depth(topo, cachelevel, HWLOC_OBJ_CACHE_DATA);
325
326         if (depth >= 0)
327         {
328             hwloc_obj_t cache = hwloc_get_next_obj_by_depth(topo, depth, nullptr);
329             if (cache != nullptr)
330             {
331                 std::vector<const hwloc_obj*> hwThreads =
332                         getHwLocDescendantsByType(topo, cache, HWLOC_OBJ_PU);
333
334                 machine->caches.push_back({ static_cast<int>(cache->attr->cache.depth),
335                                             static_cast<std::size_t>(cache->attr->cache.size),
336                                             static_cast<int>(cache->attr->cache.linesize),
337                                             static_cast<int>(cache->attr->cache.associativity),
338                                             std::max<int>(hwThreads.size(), 1) });
339             }
340         }
341     }
342     return !machine->caches.empty();
343 }
344
345 /*! \brief Read numa information from hwloc topology
346  *
347  *  \param topo    hwloc topology handle that has been initialized and loaded
348  *  \param machine Pointer to the machine structure in the HardwareTopology
349  *                 class, where numa information will be filled.
350  *
351  *  Hwloc should virtually always be able to detect numa information, but if
352  *  there is only a single numa node in the system it is not reported at all.
353  *  In this case we create a single numa node covering all cores.
354  *
355  *  This function uses the basic socket/core/thread information detected by
356  *  parseHwLocSocketsCoresThreads(), which means that routine must have
357  *  completed successfully before calling this one. If this is not the case,
358  *  you will get an error return code.
359  *
360  *  \return If the data found makes sense (either in the numa node or the
361  *          entire machine)
362  */
363 bool parseHwLocNuma(hwloc_topology_t topo, HardwareTopology::Machine* machine)
364 {
365     const hwloc_obj*              root = hwloc_get_root_obj(topo);
366     std::vector<const hwloc_obj*> hwlocNumaNodes =
367             getHwLocDescendantsByType(topo, root, HWLOC_OBJ_NUMANODE);
368     bool topologyOk = true;
369
370     if (!hwlocNumaNodes.empty())
371     {
372         machine->numa.nodes.resize(hwlocNumaNodes.size());
373
374         for (std::size_t i = 0; i < hwlocNumaNodes.size(); i++)
375         {
376             machine->numa.nodes[i].id     = hwlocNumaNodes[i]->logical_index;
377             machine->numa.nodes[i].memory = getHwLocObjectMemory(hwlocNumaNodes[i]);
378
379             machine->numa.nodes[i].logicalProcessorId.clear();
380
381             // Get list of PUs in this numa node. Get from numa node if v1.x.x, get from numa node's parent if 2.x.x
382 #    if GMX_HWLOC_API_VERSION_IS_2XX
383             std::vector<const hwloc_obj*> hwlocPUs =
384                     getHwLocDescendantsByType(topo, hwlocNumaNodes[i]->parent, HWLOC_OBJ_PU);
385 #    else
386             std::vector<const hwloc_obj*> hwlocPUs =
387                     getHwLocDescendantsByType(topo, hwlocNumaNodes[i], HWLOC_OBJ_PU);
388 #    endif
389             for (auto& p : hwlocPUs)
390             {
391                 machine->numa.nodes[i].logicalProcessorId.push_back(p->os_index);
392
393                 GMX_RELEASE_ASSERT(p->os_index < machine->logicalProcessors.size(),
394                                    "OS index of PU in hwloc larger than processor count");
395
396                 machine->logicalProcessors[p->os_index].numaNodeId = static_cast<int>(i);
397                 std::size_t s = machine->logicalProcessors[p->os_index].socketRankInMachine;
398                 std::size_t c = machine->logicalProcessors[p->os_index].coreRankInSocket;
399
400                 GMX_RELEASE_ASSERT(s < machine->sockets.size(),
401                                    "Socket index in logicalProcessors larger than socket count");
402                 GMX_RELEASE_ASSERT(c < machine->sockets[s].cores.size(),
403                                    "Core index in logicalProcessors larger than core count");
404                 // Set numaNodeId in core too
405                 machine->sockets[s].cores[c].numaNodeId = i;
406             }
407         }
408         // Getting the distance matrix
409 #    if GMX_HWLOC_API_VERSION_IS_2XX
410         // with hwloc api v. 2.x.x, distances are no longer directly accessible. Need to retrieve and release hwloc_distances_s object
411         // In addition, there can now be multiple types of distances, ie latency, bandwidth. We look only for latency, but have to check
412         // if multiple distance matrices are returned.
413
414         // If only 1 numa node exists, the v2.x.x hwloc api won't have a distances matrix, set manually
415         if (hwlocNumaNodes.size() == 1)
416         {
417             machine->numa.relativeLatency = { { 1.0 } };
418         }
419         else
420         {
421             hwloc_distances_s* dist;
422             // Set the number of distance matrices to return (1 in our case, but hwloc 2.x.x allows
423             // for multiple distances types and therefore multiple distance matrices)
424             unsigned nr = 1;
425             hwloc_distances_get(topo, &nr, &dist, HWLOC_DISTANCES_KIND_MEANS_LATENCY, 0);
426             // If no distances were found, nr will be 0, otherwise distances will be populated with
427             // 1 hwloc_distances_s object
428             if (nr > 0 && dist->nbobjs == hwlocNumaNodes.size())
429             {
430
431                 machine->numa.relativeLatency.resize(dist->nbobjs);
432                 for (std::size_t i = 0; i < dist->nbobjs; i++)
433                 {
434                     machine->numa.relativeLatency[i].resize(dist->nbobjs);
435                     for (std::size_t j = 0; j < dist->nbobjs; j++)
436                     {
437                         machine->numa.relativeLatency[i][j] = dist->values[i * dist->nbobjs + j];
438                     }
439                 }
440             }
441             else
442             {
443                 topologyOk = false;
444             }
445             hwloc_distances_release(topo, dist);
446         }
447
448         // hwloc-2.x provides latencies as integers, but to make things more similar to the case of
449         // a single numa node as well as hwloc-1.x, we rescale to relative floating-point values and
450         // also set the largest relative latency value.
451
452         // find smallest value in matrix
453         float minLatency = std::numeric_limits<float>::max(); // large number
454         float maxLatency = std::numeric_limits<float>::min(); // 0.0
455         for (const auto& v : machine->numa.relativeLatency)
456         {
457             auto result = std::minmax_element(v.begin(), v.end());
458             minLatency  = std::min(minLatency, *result.first);
459             maxLatency  = std::max(maxLatency, *result.second);
460         }
461
462         // assign stuff
463         for (auto& v : machine->numa.relativeLatency)
464         {
465             // Rescale the latencies to a relative float-point value
466             for (float& value : v)
467             {
468                 value /= minLatency;
469             }
470         }
471         machine->numa.baseLatency = 1.0; // latencies still do not have any units in hwloc-2.x
472         machine->numa.maxRelativeLatency = maxLatency / minLatency;
473
474 #    else  // GMX_HWLOC_API_VERSION_IS_2XX == false, hwloc api is 1.x.x
475         int                             depth = hwloc_get_type_depth(topo, HWLOC_OBJ_NUMANODE);
476         const struct hwloc_distances_s* dist = hwloc_get_whole_distance_matrix_by_depth(topo, depth);
477         if (dist != nullptr && dist->nbobjs == hwlocNumaNodes.size())
478         {
479             machine->numa.baseLatency        = dist->latency_base;
480             machine->numa.maxRelativeLatency = dist->latency_max;
481             machine->numa.relativeLatency.resize(dist->nbobjs);
482             for (std::size_t i = 0; i < dist->nbobjs; i++)
483             {
484                 machine->numa.relativeLatency[i].resize(dist->nbobjs);
485                 for (std::size_t j = 0; j < dist->nbobjs; j++)
486                 {
487                     machine->numa.relativeLatency[i][j] = dist->latency[i * dist->nbobjs + j];
488                 }
489             }
490         }
491         else
492         {
493             topologyOk = false;
494         }
495 #    endif // end GMX_HWLOC_API_VERSION_IS_2XX == false
496     }
497     else
498     // Deals with the case of no numa nodes found.
499 #    if GMX_HWLOC_API_VERSION_IS_2XX
500     // If the hwloc version is 2.x.x, and there is no numa node, something went wrong
501     {
502         topologyOk = false;
503     }
504 #    else
505     {
506         // No numa nodes found. Use the entire machine as a numa node.
507         // Note that this should only be the case with hwloc api v 1.x.x,
508         // a numa node is assigned to the machine by default in v 2.x.x
509         const hwloc_obj* const hwlocMachine = hwloc_get_next_obj_by_type(topo, HWLOC_OBJ_MACHINE, nullptr);
510
511         if (hwlocMachine != nullptr)
512         {
513             machine->numa.nodes.resize(1);
514             machine->numa.nodes[0].id        = 0;
515             machine->numa.nodes[0].memory    = hwlocMachine->memory.total_memory;
516             machine->numa.baseLatency        = 10;
517             machine->numa.maxRelativeLatency = 1;
518             machine->numa.relativeLatency    = { { 1.0 } };
519
520             for (int i = 0; i < machine->logicalProcessorCount; i++)
521             {
522                 machine->numa.nodes[0].logicalProcessorId.push_back(i);
523             }
524             for (auto& l : machine->logicalProcessors)
525             {
526                 l.numaNodeId = 0;
527             }
528             for (auto& s : machine->sockets)
529             {
530                 for (auto& c : s.cores)
531                 {
532                     c.numaNodeId = 0;
533                 }
534             }
535         }
536         else
537         {
538             topologyOk = false;
539         }
540     }
541 #    endif // end if not GMX_HWLOC_API_VERSION_IS_2XX
542     if (!topologyOk)
543     {
544         machine->numa.nodes.clear();
545     }
546     return topologyOk;
547 }
548
549 /*! \brief Read PCI device information from hwloc topology
550  *
551  *  \param topo    hwloc topology handle that has been initialized and loaded
552  *  \param machine Pointer to the machine structure in the HardwareTopology
553  *                 class, where PCI device information will be filled.
554  * *
555  *  \return If any devices were found
556  */
557 bool parseHwLocDevices(hwloc_topology_t topo, HardwareTopology::Machine* machine)
558 {
559     const hwloc_obj*              root = hwloc_get_root_obj(topo);
560     std::vector<const hwloc_obj*> pcidevs = getHwLocDescendantsByType(topo, root, HWLOC_OBJ_PCI_DEVICE);
561
562     for (auto& p : pcidevs)
563     {
564 #    if GMX_HWLOC_API_VERSION_IS_2XX
565         const hwloc_obj* ancestor = nullptr;
566         // Numa nodes not directly part of tree. Walk up the tree until we find an ancestor with a numa node
567         hwloc_obj_t parent = p->parent;
568         while (parent && !parent->memory_arity)
569         {
570             parent = parent->parent;
571         }
572         if (parent)
573         {
574             ancestor = parent->memory_first_child;
575         }
576 #    else  // GMX_HWLOC_API_VERSION_IS_2XX = false, api v 1.x.x
577         // numa nodes are normal part of tree, can use hwloc ancestor function
578         const hwloc_obj* const ancestor =
579                 hwloc_get_ancestor_obj_by_type(topo, HWLOC_OBJ_NUMANODE, const_cast<hwloc_obj_t>(p));
580 #    endif // end if GMX_HWLOC_API_VERSION_IS_2XX
581         int numaId;
582         if (ancestor != nullptr)
583         {
584             numaId = ancestor->logical_index;
585         }
586         else
587         {
588             // If we only have a single numa node we belong to it, otherwise set it to -1 (unknown)
589             numaId = (machine->numa.nodes.size() == 1) ? 0 : -1;
590         }
591
592         GMX_RELEASE_ASSERT(p->attr, "Attributes should not be NULL for hwloc PCI object");
593
594         machine->devices.push_back({ p->attr->pcidev.vendor_id,
595                                      p->attr->pcidev.device_id,
596                                      p->attr->pcidev.class_id,
597                                      p->attr->pcidev.domain,
598                                      p->attr->pcidev.bus,
599                                      p->attr->pcidev.dev,
600                                      p->attr->pcidev.func,
601                                      numaId });
602     }
603     return !pcidevs.empty();
604 }
605
606 void parseHwLoc(HardwareTopology::Machine* machine, HardwareTopology::SupportLevel* supportLevel, bool* isThisSystem)
607 {
608     hwloc_topology_t topo;
609
610     // Initialize a hwloc object, set flags to request IO device information too,
611     // try to load the topology, and get the root object. If either step fails,
612     // return that we do not have any support at all from hwloc.
613     if (hwloc_topology_init(&topo) != 0)
614     {
615         hwloc_topology_destroy(topo);
616         return; // SupportLevel::None.
617     }
618
619     // Flags to look for io devices
620 #    if GMX_HWLOC_API_VERSION_IS_2XX
621     GMX_RELEASE_ASSERT(
622             (hwloc_get_api_version() >= 0x20000),
623             "Mismatch between hwloc headers and library, using v2 headers with v1 library");
624     hwloc_topology_set_io_types_filter(topo, HWLOC_TYPE_FILTER_KEEP_IMPORTANT);
625 #    else
626     GMX_RELEASE_ASSERT(
627             (hwloc_get_api_version() < 0x20000),
628             "Mismatch between hwloc headers and library, using v1 headers with v2 library");
629     hwloc_topology_set_flags(topo, HWLOC_TOPOLOGY_FLAG_IO_DEVICES);
630 #    endif
631
632     if (hwloc_topology_load(topo) != 0 || hwloc_get_root_obj(topo) == nullptr)
633     {
634         hwloc_topology_destroy(topo);
635         return; // SupportLevel::None.
636     }
637
638     // If we get here, we can get a valid root object for the topology
639     *isThisSystem = hwloc_topology_is_thissystem(topo) != 0;
640
641     // Parse basic information about sockets, cores, and hardware threads
642     if (parseHwLocSocketsCoresThreads(topo, machine))
643     {
644         *supportLevel = HardwareTopology::SupportLevel::Basic;
645     }
646     else
647     {
648         hwloc_topology_destroy(topo);
649         return; // SupportLevel::None.
650     }
651
652     // Get information about cache and numa nodes
653     if (parseHwLocCache(topo, machine) && parseHwLocNuma(topo, machine))
654     {
655         *supportLevel = HardwareTopology::SupportLevel::Full;
656     }
657     else
658     {
659         hwloc_topology_destroy(topo);
660         return; // SupportLevel::Basic.
661     }
662
663     // PCI devices
664     if (parseHwLocDevices(topo, machine))
665     {
666         *supportLevel = HardwareTopology::SupportLevel::FullWithDevices;
667     }
668
669     hwloc_topology_destroy(topo);
670     // SupportLevel::Full or SupportLevel::FullWithDevices.
671 }
672
673 #endif
674
675 /*! \brief Try to detect the number of logical processors.
676  *
677  *  \return The number of hardware processing units, or 0 if it fails.
678  */
679 int detectLogicalProcessorCount()
680 {
681     int count = 0;
682
683     {
684 #if GMX_NATIVE_WINDOWS
685         // Windows
686         SYSTEM_INFO sysinfo;
687         GetSystemInfo(&sysinfo);
688         count = sysinfo.dwNumberOfProcessors;
689 #elif defined(HAVE_SYSCONF) && defined(_SC_NPROCESSORS_ONLN)
690         // We are probably on Unix. Check if we have the argument to use before executing any calls
691         count = sysconf(_SC_NPROCESSORS_ONLN);
692 #else
693         count = 0; // Neither windows nor Unix.
694 #endif
695     }
696
697     return count;
698 }
699
700 } // namespace
701
702 // static
703 HardwareTopology HardwareTopology::detect()
704 {
705     HardwareTopology result;
706
707 #if GMX_USE_HWLOC
708     parseHwLoc(&result.machine_, &result.supportLevel_, &result.isThisSystem_);
709 #endif
710
711     // If something went wrong in hwloc (or if it was not present) we might
712     // have more information in cpuInfo
713     if (result.supportLevel_ < SupportLevel::Basic)
714     {
715         // There might be topology information in cpuInfo
716         parseCpuInfo(&result.machine_, &result.supportLevel_);
717     }
718     // If we did not manage to get anything from either hwloc or cpuInfo, find the cpu count at least
719     if (result.supportLevel_ == SupportLevel::None)
720     {
721         // No topology information; try to detect the number of logical processors at least
722         result.machine_.logicalProcessorCount = detectLogicalProcessorCount();
723         if (result.machine_.logicalProcessorCount > 0)
724         {
725             result.supportLevel_ = SupportLevel::LogicalProcessorCount;
726         }
727     }
728     return result;
729 }
730
731 HardwareTopology::Machine::Machine()
732 {
733     logicalProcessorCount   = 0;
734     numa.baseLatency        = 0.0;
735     numa.maxRelativeLatency = 0.0;
736 }
737
738
739 HardwareTopology::HardwareTopology() :
740     supportLevel_(SupportLevel::None), machine_(), isThisSystem_(true)
741 {
742 }
743
744 HardwareTopology::HardwareTopology(int logicalProcessorCount) :
745     supportLevel_(SupportLevel::None), machine_(), isThisSystem_(true)
746 {
747     if (logicalProcessorCount > 0)
748     {
749         machine_.logicalProcessorCount = logicalProcessorCount;
750         supportLevel_                  = SupportLevel::LogicalProcessorCount;
751     }
752 }
753
754 int HardwareTopology::numberOfCores() const
755 {
756     if (supportLevel() >= SupportLevel::Basic)
757     {
758         // We assume all sockets have the same number of cores as socket 0.
759         // Since topology information is present, we can assume there is at least one socket.
760         return machine().sockets.size() * machine().sockets[0].cores.size();
761     }
762     else if (supportLevel() >= SupportLevel::LogicalProcessorCount)
763     {
764         return machine().logicalProcessorCount;
765     }
766     else
767     {
768         return 0;
769     }
770 }
771
772 } // namespace gmx