StandardPreconditioners_serial.hpp
Go to the documentation of this file.
1/*
2 Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
3 Copyright 2019 SINTEF Digital, Mathematics and Cybernetics.
4
5 This file is part of the Open Porous Media project (OPM).
6
7 OPM is free software: you can redistribute it and/or modify
8 it under the terms of the GNU General Public License as published by
9 the Free Software Foundation, either version 3 of the License, or
10 (at your option) any later version.
11
12 OPM is distributed in the hope that it will be useful,
13 but WITHOUT ANY WARRANTY; without even the implied warranty of
14 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
15 GNU General Public License for more details.
16
17 You should have received a copy of the GNU General Public License
18 along with OPM. If not, see <http://www.gnu.org/licenses/>.
19*/
20
21#ifndef OPM_STANDARDPRECONDITIONERS_SERIAL_HPP
22#define OPM_STANDARDPRECONDITIONERS_SERIAL_HPP
23
24#if HAVE_CUDA
26#endif
27
28namespace Opm {
29
30
31template <class Operator>
32struct StandardPreconditioners<Operator, Dune::Amg::SequentialInformation, typename std::enable_if_t<!Opm::is_gpu_operator_v<Operator>>>
33{
34 static void add()
35 {
36 using namespace Dune;
37 using O = Operator;
38 using C = Dune::Amg::SequentialInformation;
40 using M = typename F::Matrix;
41 using V = typename F::Vector;
42 using P = PropertyTree;
43 F::addCreator("ilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
44 const double w = prm.get<double>("relaxation", 1.0);
45 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
46 op.getmat(), 0, w, MILU_VARIANT::ILU);
47 });
48 F::addCreator("duneilu", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
49 const double w = prm.get<double>("relaxation", 1.0);
50 const int n = prm.get<int>("ilulevel", 0);
51 const bool resort = prm.get<bool>("resort", false);
52 return getRebuildOnUpdateWrapper<Dune::SeqILU<M, V, V>>(op.getmat(), n, w, resort);
53 });
54 F::addCreator("paroverilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
55 const double w = prm.get<double>("relaxation", 1.0);
56 const int n = prm.get<int>("ilulevel", 0);
57 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
58 op.getmat(), n, w, MILU_VARIANT::ILU);
59 });
60 F::addCreator("ilun", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
61 const int n = prm.get<int>("ilulevel", 0);
62 const double w = prm.get<double>("relaxation", 1.0);
63 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
64 op.getmat(), n, w, MILU_VARIANT::ILU);
65 });
66 F::addCreator("dilu", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
67 DUNE_UNUSED_PARAMETER(prm);
68 return std::make_shared<MultithreadDILU<M, V, V>>(op.getmat());
69 });
70 F::addCreator("jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
71 const int n = prm.get<int>("repeats", 1);
72 const double w = prm.get<double>("relaxation", 1.0);
73 return getDummyUpdateWrapper<SeqJac<M, V, V>>(op.getmat(), n, w);
74 });
75 F::addCreator("gs", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
76 const int n = prm.get<int>("repeats", 1);
77 const double w = prm.get<double>("relaxation", 1.0);
78 return getDummyUpdateWrapper<SeqGS<M, V, V>>(op.getmat(), n, w);
79 });
80 F::addCreator("sor", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
81 const int n = prm.get<int>("repeats", 1);
82 const double w = prm.get<double>("relaxation", 1.0);
83 return getDummyUpdateWrapper<SeqSOR<M, V, V>>(op.getmat(), n, w);
84 });
85 F::addCreator("ssor", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
86 const int n = prm.get<int>("repeats", 1);
87 const double w = prm.get<double>("relaxation", 1.0);
88 return getDummyUpdateWrapper<SeqSSOR<M, V, V>>(op.getmat(), n, w);
89 });
90
91 // Only add AMG preconditioners to the factory if the operator
92 // is an actual matrix operator.
93 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
94 F::addCreator("amg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
95 std::string smoother = prm.get<std::string>("smoother", "paroverilu0");
96 // Make the smoother type lowercase for internal canonical representation
97 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
98 if (smoother == "ilu0" || smoother == "paroverilu0") {
99 using Smoother = SeqILU<M, V, V>;
100 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
101 } else if (smoother == "jac") {
102 using Smoother = SeqJac<M, V, V>;
103 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
104 } else if (smoother == "gs") {
105 using Smoother = SeqGS<M, V, V>;
106 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
107 } else if (smoother == "dilu") {
108 using Smoother = MultithreadDILU<M, V, V>;
109 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
110 } else if (smoother == "sor") {
111 using Smoother = SeqSOR<M, V, V>;
112 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
113 } else if (smoother == "ssor") {
114 using Smoother = SeqSSOR<M, V, V>;
115 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
116 } else if (smoother == "ilun") {
117 using Smoother = SeqILU<M, V, V>;
118 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
119 } else {
120 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
121 }
122 });
123 F::addCreator("kamg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
124 std::string smoother = prm.get<std::string>("smoother", "paroverilu0");
125 // Make the smoother type lowercase for internal canonical representation
126 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
127 if (smoother == "ilu0" || smoother == "paroverilu0") {
128 using Smoother = SeqILU<M, V, V>;
129 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
130 } else if (smoother == "jac") {
131 using Smoother = SeqJac<M, V, V>;
132 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
133 } else if (smoother == "sor") {
134 using Smoother = SeqSOR<M, V, V>;
135 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
136 } else if (smoother == "gs") {
137 using Smoother = SeqGS<M, V, V>;
138 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
139 } else if (smoother == "ssor") {
140 using Smoother = SeqSSOR<M, V, V>;
141 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
142 } else if (smoother == "ilun") {
143 using Smoother = SeqILU<M, V, V>;
144 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
145 } else {
146 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
147 }
148 });
149 F::addCreator("famg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
150 if constexpr (std::is_same_v<typename V::field_type, float>) {
151 OPM_THROW(std::logic_error, "famg requires UMFPack which is not available for floats");
152 return nullptr;
153 } else {
154 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
155 Dune::Amg::Parameters parms;
156 parms.setNoPreSmoothSteps(1);
157 parms.setNoPostSmoothSteps(1);
158 return getRebuildOnUpdateWrapper<Dune::Amg::FastAMG<O, V>>(op, crit, parms);
159 }
160 });
161
162#if HAVE_AMGX
163 // Only add AMGX for scalar matrices
164 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1) {
165 F::addCreator("amgx", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
166 auto prm_copy = prm;
167 prm_copy.put("setup_frequency", Opm::Parameters::Get<Opm::Parameters::CprReuseInterval>());
168 return std::make_shared<Amgx::AmgxPreconditioner<M, V, V>>(op.getmat(), prm_copy);
169 });
170 }
171#endif
172
173#if HAVE_HYPRE
174 // Only add Hypre for scalar matrices
175 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1 &&
176 std::is_same_v<HYPRE_Real, typename V::field_type>) {
177 F::addCreator("hypre", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
178 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, Dune::Amg::SequentialInformation>>(op.getmat(), prm, Dune::Amg::SequentialInformation());
179 });
180 }
181#endif
182 }
183
184 // Add CPRW only for the WellModelMatrixAdapter, as the method requires that the operator
185 // has the addWellPressureEquations() method (and a few more) it can not be combined with
186 // a well-less operator such as Dune::MatrixAdapter. For OPM Flow this corresponds to
187 // requiring --matrix-add-well-contributions=false (which is the default).
188 if constexpr (std::is_same_v<O, WellModelMatrixAdapter<M, V, V>>) {
189 F::addCreator(
190 "cprw",
191 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
192 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
193 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
194 }
195 using Scalar = typename V::field_type;
196 using LevelTransferPolicy
198 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
199 op, prm, weightsCalculator, pressureIndex);
200 });
201 }
202
203 F::addCreator(
204 "cpr",
205 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
206 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
207 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
208 }
209 using Scalar = typename V::field_type;
211 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
212 op, prm, weightsCalculator, pressureIndex);
213 });
214 F::addCreator(
215 "cprt",
216 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
217 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
218 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
219 }
220 using Scalar = typename V::field_type;
222 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
223 op, prm, weightsCalculator, pressureIndex);
224 });
225
226#if HAVE_CUDA
227 // Here we create the *wrapped* GPU preconditioners
228 // meaning they will act as CPU preconditioners on the outside,
229 // but copy data back and forth to the GPU as needed.
230
231 // TODO: Make this use the GPU preconditioner factory once that is up and running.
232 F::addCreator("gpuilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
233 const double w = prm.get<double>("relaxation", 1.0);
234 using field_type = typename V::field_type;
235 using GpuILU0 = typename gpuistl::
236 GpuSeqILU0<M, gpuistl::GpuVector<field_type>, gpuistl::GpuVector<field_type>>;
237 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(
238 std::make_shared<GpuILU0>(op.getmat(), w));
239 });
240
241 F::addCreator("gpuilu0float", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
242 const double w = prm.get<double>("relaxation", 1.0);
243 using block_type = typename V::block_type;
244 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
245 using matrix_type_to =
246 typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
247 using GpuILU0 = typename gpuistl::
248 GpuSeqILU0<matrix_type_to, gpuistl::GpuVector<float>, gpuistl::GpuVector<float>>;
251 auto converted = std::make_shared<Converter>(op.getmat());
252 auto adapted = std::make_shared<Adapter>(std::make_shared<GpuILU0>(converted->getConvertedMatrix(), w));
253 converted->setUnderlyingPreconditioner(adapted);
254 return converted;
255 });
256
257 F::addCreator("gpujac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
258 const double w = prm.get<double>("relaxation", 1.0);
259 using field_type = typename V::field_type;
260 using GPUJac =
263
266 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
267 std::make_shared<MatrixOwner>(op.getmat(), w));
268 });
269
270 F::addCreator("opmgpuilu0", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
271 const bool split_matrix = prm.get<bool>("split_matrix", true);
272 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
273 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
274
275 using field_type = typename V::field_type;
279 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
280 // Note: op.getmat() is passed twice, because the ILU0 needs both the CPU and GPU matrix.
281 // The first argument will be converted to a GPU matrix, and the second one is used as a CPU matrix.
282 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme));
283 });
284
285
286 F::addCreator("gpudilu", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
287 const bool split_matrix = prm.get<bool>("split_matrix", true);
288 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
289 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
290 const bool reorder = prm.get<bool>("reorder", true);
291 using field_type = typename V::field_type;
295 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
296 // Note: op.getmat() is passed twice, because the DILU needs both the CPU and GPU matrix.
297 // The first argument will be converted to a GPU matrix, and the second one is used as a CPU matrix.
298 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
299 });
300
301 F::addCreator("gpudilufloat", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
302 const bool split_matrix = prm.get<bool>("split_matrix", true);
303 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
304 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
305 const bool reorder = prm.get<bool>("reorder", true);
306
307 using block_type = typename V::block_type;
308 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
309 using matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
312 gpuistl::GpuVector<float>, GpuDILU, matrix_type_to>;
315
316
317 auto converted = std::make_shared<Converter>(op.getmat());
318 // Note: converted->getConvertedMatrix() is passed twice, because the DILU needs both the CPU and GPU matrix.
319 // The first argument will be converted to a GPU matrix, and the second one is used as a CPU matrix.
320 auto adapted = std::make_shared<Adapter>(std::make_shared<MatrixOwner>(
321 converted->getConvertedMatrix(), converted->getConvertedMatrix(),
322 split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
323 converted->setUnderlyingPreconditioner(adapted);
324 return converted;
325 });
326#endif // HAVE_CUDA
327 }
328};
329
330
331} // namespace Opm
332
333
334#endif // OPM_STANDARDPRECONDITIONERS_SERIAL_HPP
The OpenMP thread parallelized DILU preconditioner.
Definition: DILU.hpp:53
Definition: PreconditionerFactory.hpp:64
Definition: PressureBhpTransferPolicy.hpp:99
Definition: PressureTransferPolicy.hpp:55
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
DILU preconditioner on the GPU.
Definition: GpuDILU.hpp:53
Jacobi preconditioner on the GPU.
Definition: GpuJac.hpp:47
ILU0 preconditioner on the GPU.
Definition: OpmGpuILU0.hpp:51
Makes a CUDA preconditioner available to a CPU simulator.
Definition: PreconditionerAdapter.hpp:43
Convert a CPU matrix to a GPU matrix and use a CUDA preconditioner on the GPU.
Definition: PreconditionerCPUMatrixToGPUMatrix.hpp:42
Converts the field type (eg. double to float) to benchmark single precision preconditioners.
Definition: PreconditionerConvertFieldTypeAdapter.hpp:86
Definition: fvbaseprimaryvariables.hh:141
Definition: blackoilboundaryratevector.hh:39
@ ILU
Do not perform modified ILU.
Definition: PreconditionerFactory.hpp:43
static Criterion criterion(const PropertyTree &prm)
Definition: StandardPreconditioners_mpi.hpp:83
Definition: StandardPreconditioners_mpi.hpp:128