3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
initialize.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_COMMON_INITIALIZE_HH
25#define DUMUX_COMMON_INITIALIZE_HH
26
27#include <string>
28#include <algorithm>
29#include <cstdlib>
30
31#include <dune/common/parallel/mpihelper.hh>
32
33#if HAVE_TBB
34#include <oneapi/tbb/info.h>
35#include <oneapi/tbb/global_control.h>
36
37#ifndef DOXYGEN
38namespace Dumux::Detail {
39
40class TBBGlobalControl
41{
42public:
43 static oneapi::tbb::global_control& instance(int& argc, char* argv[])
44 {
45 int maxNumThreads = oneapi::tbb::info::default_concurrency();
46 if (const char* dumuxNumThreads = std::getenv("DUMUX_NUM_THREADS"))
47 maxNumThreads = std::max(1, std::stoi(std::string{ dumuxNumThreads }));
48
49 static oneapi::tbb::global_control global_limit(
50 oneapi::tbb::global_control::max_allowed_parallelism, maxNumThreads
51 );
52
53 return global_limit;
54 }
55};
56
57} // namespace Dumux::Detail
58#endif // DOXYGEN
59
60#endif // HAVE_TBB
61
62
63#if HAVE_OPENMP
64#include <omp.h>
65#endif // HAVE_OPENMP
66
67
68#if HAVE_KOKKOS
69#include <Kokkos_Core.hpp>
70
71#ifndef DOXYGEN
72namespace Dumux::Detail {
73
74class KokkosScopeGuard
75{
76public:
77 static Kokkos::ScopeGuard& instance(int& argc, char* argv[])
78 {
79 Kokkos::InitArguments arguments;
80 if (const char* dumuxNumThreads = std::getenv("DUMUX_NUM_THREADS"))
81 arguments.num_threads = std::max(1, std::stoi(std::string{ dumuxNumThreads }));
82
83 static Kokkos::ScopeGuard guard(arguments);
84 return guard;
85 }
86};
87
88} // namespace Dumux::Detail
89#endif // DOXYGEN
90
91#endif // HAVE_KOKKOS
92
93namespace Dumux {
94
95void initialize(int& argc, char* argv[])
96{
97 // initialize MPI if available
98 // otherwise this will create a sequential (fake) helper
99 Dune::MPIHelper::instance(argc, argv);
100
101#if HAVE_TBB
102 // initialize TBB and keep global control alive
103 Detail::TBBGlobalControl::instance(argc, argv);
104#endif
105
106#if HAVE_OPENMP
107 if (const char* dumuxNumThreads = std::getenv("DUMUX_NUM_THREADS"))
108 omp_set_num_threads(
109 std::max(1, std::stoi(std::string{ dumuxNumThreads }))
110 );
111#endif
112
113#if HAVE_KOKKOS
114 // initialize Kokkos (command line / environmental variable interface)
115 Detail::KokkosScopeGuard::instance(argc, argv);
116#endif
117}
118
119} // end namespace Dumux
120
121#endif
Definition: adapt.hh:29
void initialize(int &argc, char *argv[])
Definition: initialize.hh:95
Distance implementation details.
Definition: fclocalassembler.hh:42