PreconditionerFactory_impl.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#include <opm/common/ErrorMacros.hpp>
21#include <opm/common/TimingMacros.hpp>
22
24
25
37#include <opm/simulators/linalg/ilufirstelement.hh>
38#include <opm/simulators/linalg/matrixblock.hh>
39
40#include <dune/common/unused.hh>
41#include <dune/istl/owneroverlapcopy.hh>
42#include <dune/istl/paamg/amg.hh>
43#include <dune/istl/paamg/fastamg.hh>
44#include <dune/istl/paamg/kamg.hh>
45#include <dune/istl/preconditioners.hh>
46
47#include <config.h>
48#if HAVE_CUDA
55
56#endif
57
58
59namespace Opm
60{
61
62template <class Smoother>
64 static auto args(const PropertyTree& prm)
65 {
66 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
67 SmootherArgs smootherArgs;
68 smootherArgs.iterations = prm.get<int>("iterations", 1);
69 // smootherArgs.overlap=SmootherArgs::vertex;
70 // smootherArgs.overlap=SmootherArgs::none;
71 // smootherArgs.overlap=SmootherArgs::aggregate;
72 smootherArgs.relaxationFactor = prm.get<double>("relaxation", 1.0);
73 return smootherArgs;
74 }
75};
76
77template <class M, class V, class C>
79 static auto args(const PropertyTree& prm)
80 {
82 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
83 SmootherArgs smootherArgs;
84 smootherArgs.iterations = prm.get<int>("iterations", 1);
85 const int iluwitdh = prm.get<int>("iluwidth", 0);
86 smootherArgs.setN(iluwitdh);
87 const MILU_VARIANT milu = convertString2Milu(prm.get<std::string>("milutype", std::string("ilu")));
88 smootherArgs.setMilu(milu);
89 // smootherArgs.overlap=SmootherArgs::vertex;
90 // smootherArgs.overlap=SmootherArgs::none;
91 // smootherArgs.overlap=SmootherArgs::aggregate;
92 smootherArgs.relaxationFactor = prm.get<double>("relaxation", 1.0);
93 return smootherArgs;
94 }
95};
96
97
98// trailing return type with decltype used for detecting existence of setUseFixedOrder member function by overloading the setUseFixedOrder function
99template <typename C>
100auto setUseFixedOrder(C criterion, bool booleanValue) -> decltype(criterion.setUseFixedOrder(booleanValue))
101{
102 return criterion.setUseFixedOrder(booleanValue); // Set flag to ensure that the matrices in the AMG hierarchy are constructed with deterministic indices.
103}
104template <typename C>
106{
107 // do nothing, since the function setUseFixedOrder does not exist yet
108}
109
110template <class Operator, class Comm, class Matrix, class Vector>
113{
114 Criterion criterion(15, prm.get<int>("coarsenTarget", 1200));
115 criterion.setDefaultValuesIsotropic(2);
116 criterion.setAlpha(prm.get<double>("alpha", 0.33));
117 criterion.setBeta(prm.get<double>("beta", 1e-5));
118 criterion.setMaxLevel(prm.get<int>("maxlevel", 15));
119 criterion.setSkipIsolated(prm.get<bool>("skip_isolated", false));
120 criterion.setNoPreSmoothSteps(prm.get<int>("pre_smooth", 1));
121 criterion.setNoPostSmoothSteps(prm.get<int>("post_smooth", 1));
122 criterion.setDebugLevel(prm.get<int>("verbosity", 0));
123 // As the default we request to accumulate data to 1 process always as our matrix
124 // graph might be unsymmetric and hence not supported by the PTScotch/ParMetis
125 // calls in DUNE. Accumulating to 1 skips PTScotch/ParMetis
126 criterion.setAccumulate(static_cast<Dune::Amg::AccumulationMode>(prm.get<int>("accumulate", 1)));
127 criterion.setProlongationDampingFactor(prm.get<double>("prolongationdamping", 1.6));
128 criterion.setMaxDistance(prm.get<int>("maxdistance", 2));
129 criterion.setMaxConnectivity(prm.get<int>("maxconnectivity", 15));
130 criterion.setMaxAggregateSize(prm.get<int>("maxaggsize", 6));
131 criterion.setMinAggregateSize(prm.get<int>("minaggsize", 4));
132 setUseFixedOrder(criterion, true); // If possible, set flag to ensure that the matrices in the AMG hierarchy are constructed with deterministic indices.
133 return criterion;
134}
135
136template <class Operator, class Comm, class Matrix, class Vector>
137template <class Smoother>
140 const PropertyTree& prm,
141 bool useKamg)
142{
143 auto crit = criterion(prm);
145 if (useKamg) {
147 return std::make_shared<Type>(
148 op, crit, sargs, prm.get<std::size_t>("max_krylov", 1), prm.get<double>("min_reduction", 1e-1));
149 } else {
151 return std::make_shared<Type>(op, crit, sargs);
152 }
153}
154
155template <class Operator, class Comm>
157 static void add()
158 {
159 using namespace Dune;
160 using O = Operator;
161 using C = Comm;
163 using M = typename F::Matrix;
164 using V = typename F::Vector;
165 using P = PropertyTree;
166 F::addCreator("ILU0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
167 return createParILU(op, prm, comm, 0);
168 });
169 F::addCreator("ParOverILU0",
170 [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
171 return createParILU(op, prm, comm, prm.get<int>("ilulevel", 0));
172 });
173 F::addCreator("ILUn", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
174 return createParILU(op, prm, comm, prm.get<int>("ilulevel", 0));
175 });
176 F::addCreator("DuneILU", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
177 const int n = prm.get<int>("ilulevel", 0);
178 const double w = prm.get<double>("relaxation", 1.0);
179 const bool resort = prm.get<bool>("resort", false);
180 return wrapBlockPreconditioner<RebuildOnUpdatePreconditioner<Dune::SeqILU<M, V, V>>>(
181 comm, op.getmat(), n, w, resort);
182 });
183 F::addCreator("DILU", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
184 DUNE_UNUSED_PARAMETER(prm);
185 return wrapBlockPreconditioner<MultithreadDILU<M, V, V>>(comm, op.getmat());
186 });
187 F::addCreator("Jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
188 const int n = prm.get<int>("repeats", 1);
189 const double w = prm.get<double>("relaxation", 1.0);
190 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqJac<M, V, V>>>(comm, op.getmat(), n, w);
191 });
192 F::addCreator("GS", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
193 const int n = prm.get<int>("repeats", 1);
194 const double w = prm.get<double>("relaxation", 1.0);
195 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqGS<M, V, V>>>(comm, op.getmat(), n, w);
196 });
197 F::addCreator("SOR", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
198 const int n = prm.get<int>("repeats", 1);
199 const double w = prm.get<double>("relaxation", 1.0);
200 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqSOR<M, V, V>>>(comm, op.getmat(), n, w);
201 });
202 F::addCreator("SSOR", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
203 const int n = prm.get<int>("repeats", 1);
204 const double w = prm.get<double>("relaxation", 1.0);
205 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqSSOR<M, V, V>>>(comm, op.getmat(), n, w);
206 });
207
208 // Only add AMG preconditioners to the factory if the operator
209 // is the overlapping schwarz operator. This could be extended
210 // later, but at this point no other operators are compatible
211 // with the AMG hierarchy construction.
212 if constexpr (std::is_same_v<O, Dune::OverlappingSchwarzOperator<M, V, V, C>>) {
213 F::addCreator("amg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
214 using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<V, V>>;
215 const std::string smoother = prm.get<std::string>("smoother", "ParOverILU0");
216 // TODO: merge this with ILUn, and possibly simplify the factory to only work with ILU?
217 if (smoother == "ILU0" || smoother == "ParOverILU0") {
219 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
221 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
222 return prec;
223 } else if (smoother == "DILU") {
224 using SeqSmoother = Dune::MultithreadDILU<M, V, V>;
225 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
226 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
227 SmootherArgs sargs;
228 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
229 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
230 return prec;
231 } else if (smoother == "Jac") {
232 using SeqSmoother = SeqJac<M, V, V>;
233 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
234 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
235 SmootherArgs sargs;
236 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
237 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
238 return prec;
239 } else if (smoother == "GS") {
240 using SeqSmoother = SeqGS<M, V, V>;
241 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
242 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
243 SmootherArgs sargs;
244 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
245 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
246 return prec;
247 } else if (smoother == "SOR") {
248 using SeqSmoother = SeqSOR<M, V, V>;
249 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
250 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
251 SmootherArgs sargs;
252 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
253 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
254 return prec;
255 } else if (smoother == "SSOR") {
256 using SeqSmoother = SeqSSOR<M, V, V>;
257 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
258 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
259 SmootherArgs sargs;
260 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
261 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
262 return prec;
263 } else if (smoother == "ILUn") {
264 using SeqSmoother = SeqILU<M, V, V>;
265 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
266 using SmootherArgs = typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
267 SmootherArgs sargs;
268 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
269 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
270 return prec;
271 } else {
272 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
273 }
274 });
275 }
276
277 F::addCreator("cpr",
278 [](const O& op,
279 const P& prm,
280 const std::function<V()> weightsCalculator,
281 std::size_t pressureIndex,
282 const C& comm) {
283 assert(weightsCalculator);
284 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
285 OPM_THROW(std::logic_error,
286 "Pressure index out of bounds. It needs to specified for CPR");
287 }
288 using LevelTransferPolicy = Opm::PressureTransferPolicy<O, Comm, false>;
289 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
290 op, prm, weightsCalculator, pressureIndex, comm);
291 });
292 F::addCreator("cprt",
293 [](const O& op,
294 const P& prm,
295 const std::function<V()> weightsCalculator,
296 std::size_t pressureIndex,
297 const C& comm) {
298 assert(weightsCalculator);
299 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
300 OPM_THROW(std::logic_error,
301 "Pressure index out of bounds. It needs to specified for CPR");
302 }
303 using LevelTransferPolicy = Opm::PressureTransferPolicy<O, Comm, true>;
304 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
305 op, prm, weightsCalculator, pressureIndex, comm);
306 });
307
308 if constexpr (std::is_same_v<O, WellModelGhostLastMatrixAdapter<M, V, V, true>>) {
309 F::addCreator("cprw",
310 [](const O& op,
311 const P& prm,
312 const std::function<V()> weightsCalculator,
313 std::size_t pressureIndex,
314 const C& comm) {
315 assert(weightsCalculator);
316 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
317 OPM_THROW(std::logic_error,
318 "Pressure index out of bounds. It needs to specified for CPR");
319 }
320 using LevelTransferPolicy = Opm::PressureBhpTransferPolicy<O, Comm, false>;
321 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
322 op, prm, weightsCalculator, pressureIndex, comm);
323 });
324 }
325
326#if HAVE_CUDA
327 F::addCreator("CUILU0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
328 const double w = prm.get<double>("relaxation", 1.0);
329 using field_type = typename V::field_type;
330 using CuILU0 = typename Opm::cuistl::
331 CuSeqILU0<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
332 auto cuILU0 = std::make_shared<CuILU0>(op.getmat(), w);
333
334 auto adapted = std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CuILU0>>(cuILU0);
335 auto wrapped = std::make_shared<Opm::cuistl::CuBlockPreconditioner<V, V, Comm>>(adapted, comm);
336 return wrapped;
337 });
338
339 F::addCreator("CUJac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
340 const double w = prm.get<double>("relaxation", 1.0);
341 using field_type = typename V::field_type;
342 using CuJac =
343 typename Opm::cuistl::CuJac<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
344 auto cuJac = std::make_shared<CuJac>(op.getmat(), w);
345
346 auto adapted = std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CuJac>>(cuJac);
347 auto wrapped = std::make_shared<Opm::cuistl::CuBlockPreconditioner<V, V, Comm>>(adapted, comm);
348 return wrapped;
349 });
350
351 F::addCreator("CUDILU", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t, const C& comm) {
352 using field_type = typename V::field_type;
353 using CuDILU = typename Opm::cuistl::CuDILU<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
354 auto cuDILU = std::make_shared<CuDILU>(op.getmat());
355
356 auto adapted = std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CuDILU>>(cuDILU);
357 auto wrapped = std::make_shared<Opm::cuistl::CuBlockPreconditioner<V, V, Comm>>(adapted, comm);
358 return wrapped;
359 });
360#endif
361 }
362
363
365 createParILU(const Operator& op, const PropertyTree& prm, const Comm& comm, const int ilulevel)
366 {
368 using M = typename F::Matrix;
369 using V = typename F::Vector;
370
371 const double w = prm.get<double>("relaxation", 1.0);
372 const bool redblack = prm.get<bool>("redblack", false);
373 const bool reorder_spheres = prm.get<bool>("reorder_spheres", false);
374 // Already a parallel preconditioner. Need to pass comm, but no need to wrap it in a BlockPreconditioner.
375 if (ilulevel == 0) {
376 const std::size_t num_interior = interiorIfGhostLast(comm);
377 return std::make_shared<Opm::ParallelOverlappingILU0<M, V, V, Comm>>(
378 op.getmat(), comm, w, Opm::MILU_VARIANT::ILU, num_interior, redblack, reorder_spheres);
379 } else {
380 return std::make_shared<Opm::ParallelOverlappingILU0<M, V, V, Comm>>(
381 op.getmat(), comm, ilulevel, w, Opm::MILU_VARIANT::ILU, redblack, reorder_spheres);
382 }
383 }
384
389 static std::size_t interiorIfGhostLast(const Comm& comm)
390 {
391 std::size_t interior_count = 0;
392 std::size_t highest_interior_index = 0;
393 const auto& is = comm.indexSet();
394 for (const auto& ind : is) {
395 if (Comm::OwnerSet::contains(ind.local().attribute())) {
396 ++interior_count;
397 highest_interior_index = std::max(highest_interior_index, ind.local().local());
398 }
399 }
400 if (highest_interior_index + 1 == interior_count) {
401 return interior_count;
402 } else {
403 return is.size();
404 }
405 }
406};
407
408template <class Operator>
409struct StandardPreconditioners<Operator, Dune::Amg::SequentialInformation> {
410 static void add()
411 {
412 using namespace Dune;
413 using O = Operator;
414 using C = Dune::Amg::SequentialInformation;
416 using M = typename F::Matrix;
417 using V = typename F::Vector;
418 using P = PropertyTree;
419 F::addCreator("ILU0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
420 const double w = prm.get<double>("relaxation", 1.0);
421 return std::make_shared<Opm::ParallelOverlappingILU0<M, V, V, C>>(
422 op.getmat(), 0, w, Opm::MILU_VARIANT::ILU);
423 });
424 F::addCreator("DuneILU", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
425 const double w = prm.get<double>("relaxation", 1.0);
426 const int n = prm.get<int>("ilulevel", 0);
427 const bool resort = prm.get<bool>("resort", false);
428 return getRebuildOnUpdateWrapper<Dune::SeqILU<M, V, V>>(op.getmat(), n, w, resort);
429 });
430 F::addCreator("ParOverILU0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
431 const double w = prm.get<double>("relaxation", 1.0);
432 const int n = prm.get<int>("ilulevel", 0);
433 return std::make_shared<Opm::ParallelOverlappingILU0<M, V, V, C>>(
434 op.getmat(), n, w, Opm::MILU_VARIANT::ILU);
435 });
436 F::addCreator("ILUn", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
437 const int n = prm.get<int>("ilulevel", 0);
438 const double w = prm.get<double>("relaxation", 1.0);
439 return std::make_shared<Opm::ParallelOverlappingILU0<M, V, V, C>>(
440 op.getmat(), n, w, Opm::MILU_VARIANT::ILU);
441 });
442 F::addCreator("DILU", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
443 DUNE_UNUSED_PARAMETER(prm);
444 return std::make_shared<MultithreadDILU<M, V, V>>(op.getmat());
445 });
446 F::addCreator("Jac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
447 const int n = prm.get<int>("repeats", 1);
448 const double w = prm.get<double>("relaxation", 1.0);
449 return getDummyUpdateWrapper<SeqJac<M, V, V>>(op.getmat(), n, w);
450 });
451 F::addCreator("GS", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
452 const int n = prm.get<int>("repeats", 1);
453 const double w = prm.get<double>("relaxation", 1.0);
454 return getDummyUpdateWrapper<SeqGS<M, V, V>>(op.getmat(), n, w);
455 });
456 F::addCreator("SOR", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
457 const int n = prm.get<int>("repeats", 1);
458 const double w = prm.get<double>("relaxation", 1.0);
459 return getDummyUpdateWrapper<SeqSOR<M, V, V>>(op.getmat(), n, w);
460 });
461 F::addCreator("SSOR", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
462 const int n = prm.get<int>("repeats", 1);
463 const double w = prm.get<double>("relaxation", 1.0);
464 return getDummyUpdateWrapper<SeqSSOR<M, V, V>>(op.getmat(), n, w);
465 });
466
467 // Only add AMG preconditioners to the factory if the operator
468 // is an actual matrix operator.
469 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
470 F::addCreator("amg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
471 const std::string smoother = prm.get<std::string>("smoother", "ParOverILU0");
472 if (smoother == "ILU0" || smoother == "ParOverILU0") {
473 using Smoother = SeqILU<M, V, V>;
474 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
475 } else if (smoother == "Jac") {
476 using Smoother = SeqJac<M, V, V>;
477 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
478 } else if (smoother == "GS") {
479 using Smoother = SeqGS<M, V, V>;
480 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
481 } else if (smoother == "DILU") {
482 using Smoother = MultithreadDILU<M, V, V>;
483 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
484 } else if (smoother == "SOR") {
485 using Smoother = SeqSOR<M, V, V>;
486 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
487 } else if (smoother == "SSOR") {
488 using Smoother = SeqSSOR<M, V, V>;
489 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
490 } else if (smoother == "ILUn") {
491 using Smoother = SeqILU<M, V, V>;
492 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm);
493 } else {
494 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
495 }
496 });
497 F::addCreator("kamg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
498 const std::string smoother = prm.get<std::string>("smoother", "ParOverILU0");
499 if (smoother == "ILU0" || smoother == "ParOverILU0") {
500 using Smoother = SeqILU<M, V, V>;
501 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
502 } else if (smoother == "Jac") {
503 using Smoother = SeqJac<M, V, V>;
504 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
505 } else if (smoother == "SOR") {
506 using Smoother = SeqSOR<M, V, V>;
507 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
508 } else if (smoother == "GS") {
509 using Smoother = SeqGS<M, V, V>;
510 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
511 } else if (smoother == "SSOR") {
512 using Smoother = SeqSSOR<M, V, V>;
513 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
514 } else if (smoother == "ILUn") {
515 using Smoother = SeqILU<M, V, V>;
516 return AMGHelper<O, C, M, V>::template makeAmgPreconditioner<Smoother>(op, prm, true);
517 } else {
518 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
519 }
520 });
521 F::addCreator("famg", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
522 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
523 Dune::Amg::Parameters parms;
524 parms.setNoPreSmoothSteps(1);
525 parms.setNoPostSmoothSteps(1);
526 return getRebuildOnUpdateWrapper<Dune::Amg::FastAMG<O, V>>(op, crit, parms);
527 });
528 }
529 if constexpr (std::is_same_v<O, WellModelMatrixAdapter<M, V, V, false>>) {
530 F::addCreator(
531 "cprw",
532 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
533 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
534 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
535 }
536 using LevelTransferPolicy
538 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
539 op, prm, weightsCalculator, pressureIndex);
540 });
541 }
542
543 F::addCreator(
544 "cpr",
545 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
546 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
547 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
548 }
550 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
551 op, prm, weightsCalculator, pressureIndex);
552 });
553 F::addCreator(
554 "cprt",
555 [](const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
556 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
557 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
558 }
560 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
561 op, prm, weightsCalculator, pressureIndex);
562 });
563
564#if HAVE_CUDA
565 F::addCreator("CUILU0", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
566 const double w = prm.get<double>("relaxation", 1.0);
567 using field_type = typename V::field_type;
568 using CuILU0 = typename Opm::cuistl::
569 CuSeqILU0<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
570 return std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CuILU0>>(
571 std::make_shared<CuILU0>(op.getmat(), w));
572 });
573
574 F::addCreator("CUILU0Float", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
575 const double w = prm.get<double>("relaxation", 1.0);
576 using block_type = typename V::block_type;
577 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
578 using matrix_type_to =
579 typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
580 using CuILU0 = typename Opm::cuistl::
581 CuSeqILU0<matrix_type_to, Opm::cuistl::CuVector<float>, Opm::cuistl::CuVector<float>>;
584 auto converted = std::make_shared<Converter>(op.getmat());
585 auto adapted = std::make_shared<Adapter>(std::make_shared<CuILU0>(converted->getConvertedMatrix(), w));
586 converted->setUnderlyingPreconditioner(adapted);
587 return converted;
588 });
589
590 F::addCreator("CUJac", [](const O& op, const P& prm, const std::function<V()>&, std::size_t) {
591 const double w = prm.get<double>("relaxation", 1.0);
592 using field_type = typename V::field_type;
593 using CUJac =
594 typename Opm::cuistl::CuJac<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
595 return std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CUJac>>(
596 std::make_shared<CUJac>(op.getmat(), w));
597 });
598
599 F::addCreator("CUDILU", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
600 using field_type = typename V::field_type;
601 using CUDILU = typename Opm::cuistl::CuDILU<M, Opm::cuistl::CuVector<field_type>, Opm::cuistl::CuVector<field_type>>;
602 return std::make_shared<Opm::cuistl::PreconditionerAdapter<V, V, CUDILU>>(std::make_shared<CUDILU>(op.getmat()));
603 });
604
605 F::addCreator("CUDILUFloat", [](const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
606 using block_type = typename V::block_type;
607 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
608 using matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
609 using CuDILU = typename Opm::cuistl::CuDILU<matrix_type_to, Opm::cuistl::CuVector<float>, Opm::cuistl::CuVector<float>>;
612 auto converted = std::make_shared<Converter>(op.getmat());
613 auto adapted = std::make_shared<Adapter>(std::make_shared<CuDILU>(converted->getConvertedMatrix()));
614 converted->setUnderlyingPreconditioner(adapted);
615 return converted;
616 });
617#endif
618 }
619};
620
621template <class Operator, class Comm>
622PreconditionerFactory<Operator, Comm>::PreconditionerFactory()
623{
624}
625
626
627template <class Operator, class Comm>
628PreconditionerFactory<Operator, Comm>&
629PreconditionerFactory<Operator, Comm>::instance()
630{
631 static PreconditionerFactory singleton;
632 return singleton;
633}
634
635template <class Operator, class Comm>
637PreconditionerFactory<Operator, Comm>::doCreate(const Operator& op,
638 const PropertyTree& prm,
639 const std::function<Vector()> weightsCalculator,
640 std::size_t pressureIndex)
641{
642 if (!defAdded_) {
644 defAdded_ = true;
645 }
646 const std::string& type = prm.get<std::string>("type", "ParOverILU0");
647 auto it = creators_.find(type);
648 if (it == creators_.end()) {
649 std::ostringstream msg;
650 msg << "Preconditioner type " << type << " is not registered in the factory. Available types are: ";
651 for (const auto& prec : creators_) {
652 msg << prec.first << ' ';
653 }
654 msg << std::endl;
655 OPM_THROW(std::invalid_argument, msg.str());
656 }
657 return it->second(op, prm, weightsCalculator, pressureIndex);
658}
659
660template <class Operator, class Comm>
662PreconditionerFactory<Operator, Comm>::doCreate(const Operator& op,
663 const PropertyTree& prm,
664 const std::function<Vector()> weightsCalculator,
665 std::size_t pressureIndex,
666 const Comm& comm)
667{
668 if (!defAdded_) {
670 defAdded_ = true;
671 }
672 const std::string& type = prm.get<std::string>("type", "ParOverILU0");
673 auto it = parallel_creators_.find(type);
674 if (it == parallel_creators_.end()) {
675 std::ostringstream msg;
676 msg << "Parallel preconditioner type " << type << " is not registered in the factory. Available types are: ";
677 for (const auto& prec : parallel_creators_) {
678 msg << prec.first << ' ';
679 }
680 msg << std::endl;
681 OPM_THROW(std::invalid_argument, msg.str());
682 }
683 return it->second(op, prm, weightsCalculator, pressureIndex, comm);
684}
685
686template <class Operator, class Comm>
687void
688PreconditionerFactory<Operator, Comm>::doAddCreator(const std::string& type, Creator c)
689{
690 creators_[type] = c;
691}
692
693template <class Operator, class Comm>
694void
695PreconditionerFactory<Operator, Comm>::doAddCreator(const std::string& type, ParCreator c)
696{
697 parallel_creators_[type] = c;
698}
699
700template <class Operator, class Comm>
703 const PropertyTree& prm,
704 const std::function<Vector()>& weightsCalculator,
705 std::size_t pressureIndex)
706{
707 return instance().doCreate(op, prm, weightsCalculator, pressureIndex);
708}
709
710template <class Operator, class Comm>
713 const PropertyTree& prm,
714 const std::function<Vector()>& weightsCalculator,
715 const Comm& comm,
716 std::size_t pressureIndex)
717{
718 return instance().doCreate(op, prm, weightsCalculator, pressureIndex, comm);
719}
720
721
722template <class Operator, class Comm>
725 const PropertyTree& prm,
726 const Comm& comm,
727 std::size_t pressureIndex)
728{
729 return instance().doCreate(op, prm, std::function<Vector()>(), pressureIndex, comm);
730}
731
732template <class Operator, class Comm>
733void
735{
736 instance().doAddCreator(type, creator);
737}
738
739template <class Operator, class Comm>
740void
742{
743 instance().doAddCreator(type, creator);
744}
745
746using CommSeq = Dune::Amg::SequentialInformation;
747
748template <int Dim>
749using OpFSeq = Dune::MatrixAdapter<Dune::BCRSMatrix<Dune::FieldMatrix<double, Dim, Dim>>,
750 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
751 Dune::BlockVector<Dune::FieldVector<double, Dim>>>;
752template <int Dim>
753using OpBSeq = Dune::MatrixAdapter<Dune::BCRSMatrix<Opm::MatrixBlock<double, Dim, Dim>>,
754 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
755 Dune::BlockVector<Dune::FieldVector<double, Dim>>>;
756
757template <int Dim, bool overlap>
759 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
760 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
761 overlap>;
762
763template <int Dim, bool overlap>
765 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
766 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
767 overlap>;
768
769#if HAVE_MPI
770using CommPar = Dune::OwnerOverlapCopyCommunication<int, int>;
771
772template <int Dim>
773using OpFPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<Dune::FieldMatrix<double, Dim, Dim>>,
774 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
775 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
776 CommPar>;
777
778template <int Dim>
779using OpBPar = Dune::OverlappingSchwarzOperator<Dune::BCRSMatrix<MatrixBlock<double, Dim, Dim>>,
780 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
781 Dune::BlockVector<Dune::FieldVector<double, Dim>>,
782 CommPar>;
783
784#define INSTANCE_PF_PAR(Dim) \
785 template class PreconditionerFactory<OpBSeq<Dim>, CommPar>; \
786 template class PreconditionerFactory<OpFPar<Dim>, CommPar>; \
787 template class PreconditionerFactory<OpBPar<Dim>, CommPar>; \
788 template class PreconditionerFactory<OpW<Dim, false>, CommPar>; \
789 template class PreconditionerFactory<OpWG<Dim, true>, CommPar>; \
790 template class PreconditionerFactory<OpBPar<Dim>, CommSeq>;
791#endif
792
793#define INSTANCE_PF_SEQ(Dim) \
794 template class PreconditionerFactory<OpFSeq<Dim>, CommSeq>; \
795 template class PreconditionerFactory<OpBSeq<Dim>, CommSeq>; \
796 template class PreconditionerFactory<OpW<Dim, false>, CommSeq>; \
797 template class PreconditionerFactory<OpWG<Dim, true>, CommSeq>;
798
799#if HAVE_MPI
800#define INSTANCE_PF(Dim) \
801 INSTANCE_PF_PAR(Dim) \
802 INSTANCE_PF_SEQ(Dim)
803#else
804#define INSTANCE_PF(Dim) INSTANCE_PF_SEQ(Dim)
805#endif
806} // namespace Opm
Dune::OwnerOverlapCopyCommunication< int, int > Comm
Definition: FlexibleSolver_impl.hpp:270
The AMG preconditioner.
Parallel algebraic multigrid based on agglomeration.
Definition: amgcpr.hh:88
Definition: PreconditionerWithUpdate.hpp:40
The OpenMP thread parallelized DILU preconditioner.
Definition: DILU.hpp:53
A two-step version of an overlapping Schwarz preconditioner using one step ILU0 as.
Definition: ParallelOverlappingILU0.hpp:131
Definition: PreconditionerFactory.hpp:64
static PrecPtr create(const Operator &op, const PropertyTree &prm, const std::function< Vector()> &weightsCalculator={}, std::size_t pressureIndex=std::numeric_limits< std::size_t >::max())
Definition: PreconditionerFactory_impl.hpp:702
std::function< PrecPtr(const Operator &, const PropertyTree &, const std::function< Vector()> &, std::size_t, const Comm &)> ParCreator
Definition: PreconditionerFactory.hpp:77
typename Operator::domain_type Vector
Definition: PreconditionerFactory.hpp:68
std::function< PrecPtr(const Operator &, const PropertyTree &, const std::function< Vector()> &, std::size_t)> Creator
The type of creator functions passed to addCreator().
Definition: PreconditionerFactory.hpp:75
static void addCreator(const std::string &type, Creator creator)
Definition: PreconditionerFactory_impl.hpp:734
std::shared_ptr< Dune::PreconditionerWithUpdate< Vector, Vector > > PrecPtr
The type of pointer returned by create().
Definition: PreconditionerFactory.hpp:71
Definition: PressureBhpTransferPolicy.hpp:93
Definition: PressureTransferPolicy.hpp:53
Definition: PropertyTree.hpp:37
T get(const std::string &key) const
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:224
Adapter to combine a matrix and another linear operator into a combined linear operator.
Definition: WellOperators.hpp:134
DILU preconditioner on the GPU.
Definition: CuDILU.hpp:48
Jacobi preconditioner on the GPU.
Definition: CuJac.hpp:47
Makes a CUDA preconditioner available to a CPU simulator.
Definition: PreconditionerAdapter.hpp:43
Converts the field type (eg. double to float) to benchmark single precision preconditioners.
Definition: PreconditionerConvertFieldTypeAdapter.hpp:86
Definition: SupportsFaceTag.hpp:27
Definition: BlackoilPhases.hpp:27
MILU_VARIANT
Definition: MILU.hpp:34
@ ILU
Do not perform modified ILU.
Dune::OverlappingSchwarzOperator< Dune::BCRSMatrix< MatrixBlock< double, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > >, CommPar > OpBPar
Definition: PreconditionerFactory_impl.hpp:782
Dune::MatrixAdapter< Dune::BCRSMatrix< Opm::MatrixBlock< double, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > > > OpBSeq
Definition: PreconditionerFactory_impl.hpp:755
Dune::OverlappingSchwarzOperator< Dune::BCRSMatrix< Dune::FieldMatrix< double, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > >, CommPar > OpFPar
Definition: PreconditionerFactory_impl.hpp:776
auto setUseFixedOrder(C criterion, bool booleanValue) -> decltype(criterion.setUseFixedOrder(booleanValue))
Definition: PreconditionerFactory_impl.hpp:100
Dune::Amg::SequentialInformation CommSeq
Definition: PreconditionerFactory_impl.hpp:746
Dune::MatrixAdapter< Dune::BCRSMatrix< Dune::FieldMatrix< double, Dim, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > >, Dune::BlockVector< Dune::FieldVector< double, Dim > > > OpFSeq
Definition: PreconditionerFactory_impl.hpp:751
MILU_VARIANT convertString2Milu(const std::string &milu)
Dune::OwnerOverlapCopyCommunication< int, int > CommPar
Definition: PreconditionerFactory_impl.hpp:770
Definition: PreconditionerFactory.hpp:43
Dune::Amg::CoarsenCriterion< CriterionBase > Criterion
Definition: PreconditionerFactory.hpp:47
static Criterion criterion(const PropertyTree &prm)
Definition: PreconditionerFactory_impl.hpp:112
std::shared_ptr< Dune::PreconditionerWithUpdate< Vector, Vector > > PrecPtr
Definition: PreconditionerFactory.hpp:44
static PrecPtr makeAmgPreconditioner(const Operator &op, const PropertyTree &prm, bool useKamg=false)
Definition: PreconditionerFactory_impl.hpp:139
static auto args(const PropertyTree &prm)
Definition: PreconditionerFactory_impl.hpp:79
Definition: PreconditionerFactory_impl.hpp:63
static auto args(const PropertyTree &prm)
Definition: PreconditionerFactory_impl.hpp:64
static void add()
Definition: PreconditionerFactory_impl.hpp:410
Definition: PreconditionerFactory_impl.hpp:156
static std::size_t interiorIfGhostLast(const Comm &comm)
Definition: PreconditionerFactory_impl.hpp:389
static PreconditionerFactory< Operator, Comm >::PrecPtr createParILU(const Operator &op, const PropertyTree &prm, const Comm &comm, const int ilulevel)
Definition: PreconditionerFactory_impl.hpp:365
static void add()
Definition: PreconditionerFactory_impl.hpp:157