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
34#endif
35
36namespace Opm {
37
38template <class Operator, class Comm>
39struct StandardPreconditioners<Operator, Comm, typename std::enable_if_t<Opm::is_gpu_operator_v<Operator>>>
40{
41 using O = Operator;
42 using C = Comm;
44 using M = typename F::Matrix;
45 using V = typename F::Vector;
46 using P = PropertyTree;
47 using PrecPtr = typename F::PrecPtr;
48
49 using field_type = typename V::field_type;
50
51 static constexpr int maxblocksize = 6;
52
53 static void add()
54 {
55#if HAVE_CUDA
56 F::addCreator("ilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
57 const double w = prm.get<double>("relaxation", 1.0);
58 using GpuILU0 = typename gpuistl::GpuSeqILU0<M, V, V>;
59 auto gpuPrec = std::make_shared<GpuILU0>(op.getmat(), w);
60 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
61 });
62
63 F::addCreator("jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
64 const double w = prm.get<double>("relaxation", 1.0);
65 using GpuJac = typename gpuistl::GpuJac<M, V, V>;
66 auto gpuPrec = std::make_shared<GpuJac>(op.getmat(), w);
67 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
68 });
69
70 F::addCreator("dilu", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t, const C& comm) -> PrecPtr {
71 return op.getmat().dispatchOnBlocksize([&](auto blockSizeVal) -> PrecPtr {
72 constexpr int blockSize = decltype(blockSizeVal)::value;
73 const auto cpuMatrix = gpuistl::detail::makeCPUMatrix<O, Dune::FieldMatrix<field_type, blockSize, blockSize>>(op);
74 const bool split_matrix = prm.get<bool>("split_matrix", true);
75 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
76 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
77 const bool reorder = prm.get<bool>("reorder", true);
78 using CPUMatrixType = std::remove_const_t<std::remove_reference_t<decltype(cpuMatrix)>>;
79 using GPUDILU = typename gpuistl::GpuDILU<CPUMatrixType, V, V>;
80 auto gpuPrec = std::make_shared<GPUDILU>(op.getmat(), cpuMatrix, split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder);
81 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
82 });
83 });
84
85 F::addCreator("opmilu0", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t, const C& comm) -> PrecPtr {
86 return op.getmat().dispatchOnBlocksize([&](auto blockSizeVal) -> PrecPtr {
87 constexpr int blockSize = decltype(blockSizeVal)::value;
88 const auto cpuMatrix = gpuistl::detail::makeCPUMatrix<O, Dune::FieldMatrix<field_type, blockSize, blockSize>>(op);
89 const bool split_matrix = prm.get<bool>("split_matrix", true);
90 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
91 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
92 using CPUMatrixType = std::remove_const_t<std::remove_reference_t<decltype(cpuMatrix)>>;
93 using GPUILU0 = typename gpuistl::OpmGpuILU0<CPUMatrixType, V, V>;
94 auto gpuPrec = std::make_shared<GPUILU0>(op.getmat(), cpuMatrix, split_matrix, tune_gpu_kernels, mixed_precision_scheme);
95 return std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, C>>(gpuPrec, comm);
96 });
97 });
98#endif // HAVE_CUDA
99 }
100
101};
102
103
104}// namespace Opm
105
106
107#endif // OPM_STANDARDPRECONDITIONERS_GPU_MPI_HEADER
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:304
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:43
typename V::field_type field_type
Definition: StandardPreconditioners_gpu_mpi.hpp:49
Definition: StandardPreconditioners_mpi.hpp:132