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