StandardPreconditioners_gpu_serial.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_SERIAL_HEADER
18#define OPM_STANDARDPRECONDITIONERS_GPU_SERIAL_HEADER
19
20#if USE_HIP
21#include <opm/simulators/linalg/gpuistl_hip/detail/gpu_preconditioner_utils.hpp>
22#else
24#endif
25
26#include <dune/istl/bcrsmatrix.hh>
27
28#include <type_traits>
29
30
31
32namespace Opm
33{
34
35
36template <class Operator>
38 Dune::Amg::SequentialInformation,
39 typename std::enable_if_t<Opm::is_gpu_operator_v<Operator>>>
40{
41
42 using O = Operator;
43 using C = Dune::Amg::SequentialInformation;
45 using M = typename F::Matrix;
46 using V = typename F::Vector;
47 using P = PropertyTree;
48 using PrecPtr = typename F::PrecPtr;
49
50 using field_type = typename V::field_type;
51
52 static constexpr int maxblocksize = 6;
53 static void add()
54 {
55 F::addCreator("ilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
56 const double w = prm.get<double>("relaxation", 1.0);
57 using GpuILU0 =
59 return std::make_shared<GpuILU0>(op.getmat(), w);
60 });
61
62 F::addCreator("jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
63 const double w = prm.get<double>("relaxation", 1.0);
65 return std::make_shared<GPUJac>(op.getmat(), w);
66 });
67
68 // The next two (OPMILU0 and DILU) are the GPU preconditioners that need CPU matrices.
69 // Since we are not storing the block size compile time for the GPU matrices
70 // **and** we need the CPU matrix for the creation, we need to
71 // dispatch the creation of the preconditioner based on the block size.
72 //
73 // Note that this dispatching is not needed in the future, since we will have a constructor taking GPU matrices directly.
74 F::addCreator("opmilu0", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) -> PrecPtr {
75 return op.getmat().dispatchOnBlocksize([&](auto blockSizeVal) -> PrecPtr {
76 constexpr int blockSize = decltype(blockSizeVal)::value;
77 const auto cpuMatrix = gpuistl::detail::makeCPUMatrix<O, Dune::FieldMatrix<field_type, blockSize, blockSize>>(op);
78 const bool split_matrix = prm.get<bool>("split_matrix", true);
79 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
80 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
81 using CPUMatrixType = std::remove_const_t<std::remove_reference_t<decltype(cpuMatrix)>>;
82 using GPUILU0 =
84
85 return std::make_shared<GPUILU0>(op.getmat(), cpuMatrix, split_matrix, tune_gpu_kernels, mixed_precision_scheme);
86 });
87 });
88
89 F::addCreator("dilu", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) -> PrecPtr {
90 return op.getmat().dispatchOnBlocksize([&](auto blockSizeVal) -> PrecPtr {
91 constexpr int blockSize = decltype(blockSizeVal)::value;
92 const auto cpuMatrix = gpuistl::detail::makeCPUMatrix<O, Dune::FieldMatrix<field_type, blockSize, blockSize>>(op);
93 const bool split_matrix = prm.get<bool>("split_matrix", true);
94 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
95 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
96 const bool reorder = prm.get<bool>("reorder", true);
97 using CPUMatrixType = std::remove_const_t<std::remove_reference_t<decltype(cpuMatrix)>>;
98 using GPUDILU =
100 return std::make_shared<GPUDILU>(op.getmat(), cpuMatrix, split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder);
101 });
102 });
103
104 // Only add AMG preconditioners to the factory if the operator
105 // is an actual matrix operator.
106 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
107#if HAVE_AMGX
108 // Only add AMGX for scalar matrices
109 F::addCreator("amgx", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
110 // Only create AMGX preconditioner for scalar matrices
111 if (op.getmat().blockSize() == 1) {
112 auto prm_copy = prm;
113 prm_copy.put("setup_frequency", Opm::Parameters::Get<Opm::Parameters::CprReuseInterval>());
114 return std::make_shared<Amgx::AmgxPreconditioner<M, V, V>>(op.getmat(), prm_copy);
115 } else {
116 OPM_THROW(std::logic_error, "AMGX preconditioner only works with scalar matrices (block size 1)");
117 }
118 });
119#endif
120
121#if HAVE_HYPRE && (HYPRE_USING_CUDA || HYPRE_USING_HIP)
122 // Only register Hypre preconditioner for double precision
123 if constexpr (std::is_same_v<HYPRE_Real, typename V::field_type>) {
124 F::addCreator("hypre", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
125 // Only create Hypre preconditioner for scalar matrices
126 if (op.getmat().blockSize() == 1) {
127 return std::make_shared<linalg::HyprePreconditioner<M, V, V, Dune::Amg::SequentialInformation>>(op.getmat(), prm, Dune::Amg::SequentialInformation());
128 } else {
129 OPM_THROW(std::logic_error, "Hypre preconditioner only works with scalar matrices (block size 1).");
130 }
131 });
132 }
133#endif
134 }
135
136 F::addCreator("cpr", [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
137 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
138 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
139 }
140 using Scalar = typename V::field_type;
141 using GpuVector = gpuistl::GpuVector<Scalar>;
143 return std::make_shared<Dune::OwningTwoLevelPreconditioner<O, GpuVector, LevelTransferPolicy>>(op, prm, weightsCalculator, pressureIndex);
144 });
145
146 F::addCreator("cprt", [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
147 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
148 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
149 }
150 using Scalar = typename V::field_type;
151 using GpuVector = gpuistl::GpuVector<Scalar>;
153 return std::make_shared<Dune::OwningTwoLevelPreconditioner<O, GpuVector, LevelTransferPolicy>>(op, prm, weightsCalculator, pressureIndex);
154 });
155 }
156
157
158};
159
160
161} // namespace Opm
162
163#endif // OPM_STANDARDPRECONDITIONERS_GPU_SERIAL_HEADER
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
Definition: GpuPressureTransferPolicy.hpp:54
Sequential ILU0 preconditioner on the GPU through the CuSparse library.
Definition: GpuSeqILU0.hpp:52
ILU0 preconditioner on the GPU.
Definition: OpmGpuILU0.hpp:51
Definition: fvbaseprimaryvariables.hh:161
Definition: blackoilbioeffectsmodules.hh:45
Definition: StandardPreconditioners_mpi.hpp:135