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