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
25#if USE_HIP
26#include <opm/simulators/linalg/gpuistl_hip/PreconditionerCPUMatrixToGPUMatrix.hpp>
27#else
29#endif
30#endif
31
32namespace Opm {
33
34template <class X, class Y>
36{
37 public:
39 virtual void update() override {};
40 virtual bool hasPerfectUpdate() const override {return true;}
41 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& y) override {};
42 virtual void post ([[maybe_unused]] X& x) override {};
43 virtual void apply ([[maybe_unused]] X& x, [[maybe_unused]] const Y& y) override {};
44 virtual Dune::SolverCategory::Category category() const override { return Dune::SolverCategory::sequential; };
45};
46
47
48template <class Operator>
49struct StandardPreconditioners<Operator, Dune::Amg::SequentialInformation, typename std::enable_if_t<!Opm::is_gpu_operator_v<Operator>>>
50{
51 static void add()
52 {
53 using namespace Dune;
54 using O = Operator;
55 using C = Dune::Amg::SequentialInformation;
57 using M = typename F::Matrix;
58 using V = typename F::Vector;
59 using P = PropertyTree;
60 F::addCreator("ilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
61 const double w = prm.get<double>("relaxation", 1.0);
62 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
63 op.getmat(), 0, w, MILU_VARIANT::ILU);
64 });
65 F::addCreator("duneilu", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
66 const double w = prm.get<double>("relaxation", 1.0);
67 const int n = prm.get<int>("ilulevel", 0);
68 const bool resort = prm.get<bool>("resort", false);
69 return getRebuildOnUpdateWrapper<Dune::SeqILU<M, V, V>>(op.getmat(), n, w, resort);
70 });
71 F::addCreator("paroverilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
72 const double w = prm.get<double>("relaxation", 1.0);
73 const int n = prm.get<int>("ilulevel", 0);
74 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
75 op.getmat(), n, w, MILU_VARIANT::ILU);
76 });
77 F::addCreator("ilun", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
78 const int n = prm.get<int>("ilulevel", 0);
79 const double w = prm.get<double>("relaxation", 1.0);
80 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
81 op.getmat(), n, w, MILU_VARIANT::ILU);
82 });
83 F::addCreator("dilu", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
84 DUNE_UNUSED_PARAMETER(prm);
85 return std::make_shared<MultithreadDILU<M, V, V>>(op.getmat());
86 });
87 F::addCreator("mixed-ilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
88 DUNE_UNUSED_PARAMETER(prm);
89 DUNE_UNUSED_PARAMETER(op);
90 return std::make_shared<TrivialPreconditioner<V,V>>();
91 });
92 F::addCreator("mixed-dilu", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
93 DUNE_UNUSED_PARAMETER(prm);
94 DUNE_UNUSED_PARAMETER(op);
95 return std::make_shared<TrivialPreconditioner<V,V>>();
96 });
97 F::addCreator("jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
98 const int n = prm.get<int>("repeats", 1);
99 const double w = prm.get<double>("relaxation", 1.0);
100 return getDummyUpdateWrapper<SeqJac<M, V, V>>(op.getmat(), n, w);
101 });
102 F::addCreator("gs", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
103 const int n = prm.get<int>("repeats", 1);
104 const double w = prm.get<double>("relaxation", 1.0);
105 return getDummyUpdateWrapper<SeqGS<M, V, V>>(op.getmat(), n, w);
106 });
107 F::addCreator("sor", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
108 const int n = prm.get<int>("repeats", 1);
109 const double w = prm.get<double>("relaxation", 1.0);
110 return getDummyUpdateWrapper<SeqSOR<M, V, V>>(op.getmat(), n, w);
111 });
112 F::addCreator("ssor", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
113 const int n = prm.get<int>("repeats", 1);
114 const double w = prm.get<double>("relaxation", 1.0);
115 return getDummyUpdateWrapper<SeqSSOR<M, V, V>>(op.getmat(), n, w);
116 });
117
118 // Only add AMG preconditioners to the factory if the operator
119 // is an actual matrix operator.
120 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
121 F::addCreator("amg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
122 std::string smoother = prm.get<std::string>("smoother", "paroverilu0");
123 // Make the smoother type lowercase for internal canonical representation
124 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
125 if (smoother == "ilu0" || smoother == "paroverilu0") {
126 using Smoother = SeqILU<M, V, V>;
127 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
128 } else if (smoother == "jac") {
129 using Smoother = SeqJac<M, V, V>;
130 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
131 } else if (smoother == "gs") {
132 using Smoother = SeqGS<M, V, V>;
133 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
134 } else if (smoother == "dilu") {
135 using Smoother = MultithreadDILU<M, V, V>;
136 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
137 } else if (smoother == "sor") {
138 using Smoother = SeqSOR<M, V, V>;
139 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
140 } else if (smoother == "ssor") {
141 using Smoother = SeqSSOR<M, V, V>;
142 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
143 } else if (smoother == "ilun") {
144 using Smoother = SeqILU<M, V, V>;
145 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
146 } else {
147 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
148 }
149 });
150 F::addCreator("kamg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
151 std::string smoother = prm.get<std::string>("smoother", "paroverilu0");
152 // Make the smoother type lowercase for internal canonical representation
153 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
154 if (smoother == "ilu0" || smoother == "paroverilu0") {
155 using Smoother = SeqILU<M, V, V>;
156 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
157 } else if (smoother == "jac") {
158 using Smoother = SeqJac<M, V, V>;
159 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
160 } else if (smoother == "sor") {
161 using Smoother = SeqSOR<M, V, V>;
162 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
163 } else if (smoother == "gs") {
164 using Smoother = SeqGS<M, V, V>;
165 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
166 } else if (smoother == "ssor") {
167 using Smoother = SeqSSOR<M, V, V>;
168 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
169 } else if (smoother == "ilun") {
170 using Smoother = SeqILU<M, V, V>;
171 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
172 } else {
173 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
174 }
175 });
176 F::addCreator("famg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
177 if constexpr (std::is_same_v<typename V::field_type, float>) {
178 OPM_THROW(std::logic_error, "famg requires UMFPack which is not available for floats");
179 return nullptr;
180 } else {
181 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
182 Dune::Amg::Parameters parms;
183 parms.setNoPreSmoothSteps(1);
184 parms.setNoPostSmoothSteps(1);
185 return getRebuildOnUpdateWrapper<Dune::Amg::FastAMG<O, V>>(op, crit, parms);
186 }
187 });
188
189#if HAVE_AMGX
190 // Only add AMGX for scalar matrices
191 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1) {
192 F::addCreator("amgx", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
193 auto prm_copy = prm;
194 prm_copy.put("setup_frequency", Opm::Parameters::Get<Opm::Parameters::CprReuseInterval>());
195 return std::make_shared<Amgx::AmgxPreconditioner<M, V, V>>(op.getmat(), prm_copy);
196 });
197 }
198#endif
199
200#if HAVE_HYPRE
201 // Only add Hypre for scalar matrices
202 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1 &&
203 std::is_same_v<HYPRE_Real, typename V::field_type>) {
204 F::addCreator("hypre", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
205 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, Dune::Amg::SequentialInformation>>(op.getmat(), prm, Dune::Amg::SequentialInformation());
206 });
207 }
208#endif
209 }
210
211 // Add CPRW only for the WellModelMatrixAdapter, as the method requires that the operator
212 // has the addWellPressureEquations() method (and a few more) it can not be combined with
213 // a well-less operator such as Dune::MatrixAdapter. For OPM Flow this corresponds to
214 // requiring --matrix-add-well-contributions=false (which is the default).
215 if constexpr (std::is_same_v<O, WellModelMatrixAdapter<M, V, V>>) {
216 F::addCreator(
217 "cprw",
218 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
219 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
220 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
221 }
222 using Scalar = typename V::field_type;
223 using LevelTransferPolicy
225 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
226 op, prm, weightsCalculator, pressureIndex);
227 });
228 }
229
230 F::addCreator(
231 "cpr",
232 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
233 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
234 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
235 }
236 using Scalar = typename V::field_type;
238 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
239 op, prm, weightsCalculator, pressureIndex);
240 });
241 F::addCreator(
242 "cprt",
243 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
244 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
245 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
246 }
247 using Scalar = typename V::field_type;
249 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
250 op, prm, weightsCalculator, pressureIndex);
251 });
252
253#if HAVE_CUDA
254 // Here we create the *wrapped* GPU preconditioners
255 // meaning they will act as CPU preconditioners on the outside,
256 // but copy data back and forth to the GPU as needed.
257
258 // TODO: Make this use the GPU preconditioner factory once that is up and running.
259 F::addCreator("gpuilu0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
260 const double w = prm.get<double>("relaxation", 1.0);
261 using field_type = typename V::field_type;
262 using GpuILU0 = typename gpuistl::
263 GpuSeqILU0<M, gpuistl::GpuVector<field_type>, gpuistl::GpuVector<field_type>>;
264 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(
265 std::make_shared<GpuILU0>(op.getmat(), w));
266 });
267
268 F::addCreator("gpuilu0float", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
269 const double w = prm.get<double>("relaxation", 1.0);
270 using block_type = typename V::block_type;
271 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
272 using matrix_type_to =
273 typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
274 using GpuILU0 = typename gpuistl::
275 GpuSeqILU0<matrix_type_to, gpuistl::GpuVector<float>, gpuistl::GpuVector<float>>;
278 auto converted = std::make_shared<Converter>(op.getmat());
279 auto adapted = std::make_shared<Adapter>(std::make_shared<GpuILU0>(converted->getConvertedMatrix(), w));
280 converted->setUnderlyingPreconditioner(adapted);
281 return converted;
282 });
283
284 F::addCreator("gpujac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
285 const double w = prm.get<double>("relaxation", 1.0);
286 using field_type = typename V::field_type;
287 using GPUJac =
290
293 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
294 std::make_shared<MatrixOwner>(op.getmat(), w));
295 });
296
297 F::addCreator("opmgpuilu0", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
298 const bool split_matrix = prm.get<bool>("split_matrix", true);
299 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
300 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
301
302 using field_type = typename V::field_type;
306 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
307 // Note: op.getmat() is passed twice, because the ILU0 needs both the CPU and GPU matrix.
308 // The first argument will be converted to a GPU matrix, and the second one is used as a CPU matrix.
309 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme));
310 });
311
312
313 F::addCreator("gpudilu", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
314 const bool split_matrix = prm.get<bool>("split_matrix", true);
315 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
316 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
317 const bool reorder = prm.get<bool>("reorder", true);
318 using field_type = typename V::field_type;
322 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
323 // Note: op.getmat() is passed twice, because the DILU needs both the CPU and GPU matrix.
324 // The first argument will be converted to a GPU matrix, and the second one is used as a CPU matrix.
325 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
326 });
327
328 F::addCreator("gpudilufloat", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
329 const bool split_matrix = prm.get<bool>("split_matrix", true);
330 const bool tune_gpu_kernels = prm.get<bool>("tune_gpu_kernels", true);
331 const int mixed_precision_scheme = prm.get<int>("mixed_precision_scheme", 0);
332 const bool reorder = prm.get<bool>("reorder", true);
333
334 using block_type = typename V::block_type;
335 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
336 using matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
339 gpuistl::GpuVector<float>, GpuDILU, matrix_type_to>;
342
343
344 auto converted = std::make_shared<Converter>(op.getmat());
345 // Note: converted->getConvertedMatrix() is passed twice, because the DILU needs both the CPU and GPU matrix.
346 // The first argument will be converted to a GPU matrix, and the second one is used as a CPU matrix.
347 auto adapted = std::make_shared<Adapter>(std::make_shared<MatrixOwner>(
348 converted->getConvertedMatrix(), converted->getConvertedMatrix(),
349 split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
350 converted->setUnderlyingPreconditioner(adapted);
351 return converted;
352 });
353#endif // HAVE_CUDA
354 }
355};
356
357
358} // namespace Opm
359
360
361#endif // OPM_STANDARDPRECONDITIONERS_SERIAL_HPP
The OpenMP thread parallelized DILU preconditioner.
Definition: DILU.hpp:53
Interface class adding the update() method to the preconditioner interface.
Definition: PreconditionerWithUpdate.hpp:32
Definition: PreconditionerFactory.hpp:64
Definition: PressureBhpTransferPolicy.hpp:99
Definition: PressureTransferPolicy.hpp:55
Hierarchical collection of key/value pairs.
Definition: PropertyTree.hpp:39
Definition: StandardPreconditioners_serial.hpp:36
virtual bool hasPerfectUpdate() const override
Definition: StandardPreconditioners_serial.hpp:40
virtual void pre(X &x, Y &y) override
Definition: StandardPreconditioners_serial.hpp:41
virtual void post(X &x) override
Definition: StandardPreconditioners_serial.hpp:42
virtual void update() override
Definition: StandardPreconditioners_serial.hpp:39
virtual void apply(X &x, const Y &y) override
Definition: StandardPreconditioners_serial.hpp:43
TrivialPreconditioner()
Definition: StandardPreconditioners_serial.hpp:38
virtual Dune::SolverCategory::Category category() const override
Definition: StandardPreconditioners_serial.hpp:44
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:161
Definition: blackoilbioeffectsmodules.hh:43
@ ILU
Do not perform modified ILU.
Definition: PreconditionerFactory.hpp:43
static Criterion criterion(const PropertyTree &prm)
Definition: StandardPreconditioners_mpi.hpp:87
Definition: StandardPreconditioners_mpi.hpp:132