StandardPreconditioners_gpu_mpi.hpp
Go to the documentation of this file.
1/*
2 Copyright 2025 Equinor ASA
3
4 This file is part of the Open Porous Media project (OPM).
5 OPM is free software: you can redistribute it and/or modify
6 it under the terms of the GNU General Public License as published by
7 the Free Software Foundation, either version 3 of the License, or
8 (at your option) any later version.
9 OPM is distributed in the hope that it will be useful,
10 but WITHOUT ANY WARRANTY; without even the implied warranty of
11 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
12 GNU General Public License for more details.
13 You should have received a copy of the GNU General Public License
14 along with OPM. If not, see <http://www.gnu.org/licenses/>.
15*/
16
17#ifndef OPM_STANDARDPRECONDITIONERS_GPU_MPI_HEADER
18#define OPM_STANDARDPRECONDITIONERS_GPU_MPI_HEADER
19
20#include <dune/istl/bcrsmatrix.hh>
21#include <string>
22
26
27#if HAVE_CUDA
28#if USE_HIP
29#include <opm/simulators/linalg/gpuistl_hip/GpuDILU.hpp>
30#include <opm/simulators/linalg/gpuistl_hip/GpuSeqILU0.hpp>
31#include <opm/simulators/linalg/gpuistl_hip/GpuJac.hpp>
32#include <opm/simulators/linalg/gpuistl_hip/OpmGpuILU0.hpp>
33#include <opm/simulators/linalg/gpuistl_hip/GpuBlockPreconditioner.hpp>
34#include <opm/simulators/linalg/gpuistl_hip/detail/gpu_preconditioner_utils.hpp>
35#else
42#endif
43#endif
44
45namespace Opm {
46
47template <class Operator, class Comm>
48struct StandardPreconditioners<Operator, Comm, typename std::enable_if_t<Opm::is_gpu_operator_v<Operator>>>
49{
50 using O = Operator;
51 using C = Comm;
53 using M = typename F::Matrix;
54 using V = typename F::Vector;
55 using P = PropertyTree;
56 using PrecPtr = typename F::PrecPtr;
57
58 using field_type = typename V::field_type;
59
60 static constexpr int maxblocksize = 6;
61
62 static void add()
63 {
64#if HAVE_CUDA
65 F::addCreator("ilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
66 const double w = prm.get<double>("relaxation", 1.0);
67 using GpuILU0 = typename gpuistl::GpuSeqILU0<M, V, V>;
68 auto gpuPrec = std::make_shared<GpuILU0>(op.getmat(), w);
69 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
70 });
71
72 F::addCreator("jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
73 const double w = prm.get<double>("relaxation", 1.0);
74 using GpuJac = typename gpuistl::GpuJac<M, V, V>;
75 auto gpuPrec = std::make_shared<GpuJac>(op.getmat(), w);
76 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
77 });
78
79 F::addCreator("dilu", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t, const C& comm) -> PrecPtr {
80 return op.getmat().dispatchOnBlocksize([&](auto blockSizeVal) -> PrecPtr {
81 constexpr int blockSize = decltype(blockSizeVal)::value;
82 const auto cpuMatrix = gpuistl::detail::makeCPUMatrix<O, Dune::FieldMatrix<field_type, blockSize, blockSize>>(op);
83 const bool split_matrix = prm.get<bool>("split_matrix", true);
84 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
85 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
86 const bool reorder = prm.get<bool>("reorder", true);
87 using CPUMatrixType = std::remove_const_t<std::remove_reference_t<decltype(cpuMatrix)>>;
88 using GPUDILU = typename gpuistl::GpuDILU<CPUMatrixType, V, V>;
89 auto gpuPrec = std::make_shared<GPUDILU>(op.getmat(), cpuMatrix, split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder);
90 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
91 });
92 });
93
94 F::addCreator("opmilu0", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t, const C& comm) -> PrecPtr {
95 return op.getmat().dispatchOnBlocksize([&](auto blockSizeVal) -> PrecPtr {
96 constexpr int blockSize = decltype(blockSizeVal)::value;
97 const auto cpuMatrix = gpuistl::detail::makeCPUMatrix<O, Dune::FieldMatrix<field_type, blockSize, blockSize>>(op);
98 const bool split_matrix = prm.get<bool>("split_matrix", true);
99 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
100 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
101 using CPUMatrixType = std::remove_const_t<std::remove_reference_t<decltype(cpuMatrix)>>;
102 using GPUILU0 = typename gpuistl::OpmGpuILU0<CPUMatrixType, V, V>;
103 auto gpuPrec = std::make_shared<GPUILU0>(op.getmat(), cpuMatrix, split_matrix, tune_gpu_kernels, mixed_precision_scheme);
104 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
105 });
106 });
107#endif // HAVE_CUDA
108 }
109
110};
111
112
113}// namespace Opm
114
115
116#endif // OPM_STANDARDPRECONDITIONERS_GPU_MPI_HEADER
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:325
Definition: PreconditionerFactory.hpp:64
typename Operator::domain_type Vector
Definition: PreconditionerFactory.hpp:68
std::shared_ptr< Dune::PreconditionerWithUpdate< Vector, Vector > > PrecPtr
The type of pointer returned by create().
Definition: PreconditionerFactory.hpp:71
typename Operator::matrix_type Matrix
Linear algebra types.
Definition: PreconditionerFactory.hpp:67
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
T get(const std::string &key) const
DILU preconditioner on the GPU.
Definition: GpuDILU.hpp:53
Jacobi preconditioner on the GPU.
Definition: GpuJac.hpp:47
Sequential ILU0 preconditioner on the GPU through the CuSparse library.
Definition: GpuSeqILU0.hpp:52
ILU0 preconditioner on the GPU.
Definition: OpmGpuILU0.hpp:51
Definition: blackoilbioeffectsmodules.hh:45
typename V::field_type field_type
Definition: StandardPreconditioners_gpu_mpi.hpp:58
Definition: StandardPreconditioners_mpi.hpp:135