135 using namespace Dune;
139 using M =
typename F::Matrix;
140 using V =
typename F::Vector;
142 F::addCreator(
"ilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
143 return createParILU(op, prm, comm, 0);
145 F::addCreator(
"paroverilu0",
146 [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
147 return createParILU(op, prm, comm, prm.get<
int>(
"ilulevel", 0));
149 F::addCreator(
"ilun", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
150 return createParILU(op, prm, comm, prm.get<
int>(
"ilulevel", 0));
152 F::addCreator(
"duneilu", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
153 const int n = prm.get<
int>(
"ilulevel", 0);
154 const double w = prm.get<
double>(
"relaxation", 1.0);
155 const bool resort = prm.get<
bool>(
"resort",
false);
156 return wrapBlockPreconditioner<RebuildOnUpdatePreconditioner<Dune::SeqILU<M, V, V>>>(
157 comm, op.getmat(), n, w, resort);
159 F::addCreator(
"dilu", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
160 DUNE_UNUSED_PARAMETER(prm);
161 return wrapBlockPreconditioner<MultithreadDILU<M, V, V>>(comm, op.getmat());
163 F::addCreator(
"jac", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
164 const int n = prm.get<
int>(
"repeats", 1);
165 const double w = prm.get<
double>(
"relaxation", 1.0);
166 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqJac<M, V, V>>>(comm, op.getmat(), n, w);
168 F::addCreator(
"gs", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
169 const int n = prm.get<
int>(
"repeats", 1);
170 const double w = prm.get<
double>(
"relaxation", 1.0);
171 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqGS<M, V, V>>>(comm, op.getmat(), n, w);
173 F::addCreator(
"sor", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
174 const int n = prm.get<
int>(
"repeats", 1);
175 const double w = prm.get<
double>(
"relaxation", 1.0);
176 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqSOR<M, V, V>>>(comm, op.getmat(), n, w);
178 F::addCreator(
"ssor", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
179 const int n = prm.get<
int>(
"repeats", 1);
180 const double w = prm.get<
double>(
"relaxation", 1.0);
181 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqSSOR<M, V, V>>>(comm, op.getmat(), n, w);
188 if constexpr (std::is_same_v<O, Dune::OverlappingSchwarzOperator<M, V, V, C>> ||
189 std::is_same_v<O, Opm::GhostLastMatrixAdapter<M, V, V, C>>) {
190 F::addCreator(
"amg", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
191 using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<V, V>>;
192 std::string smoother = prm.get<std::string>(
"smoother",
"paroverilu0");
194 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
196 if (smoother ==
"ilu0" || smoother ==
"paroverilu0") {
198 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
199 auto sargs = AMGSmootherArgsHelper<Smoother>::args(prm);
200 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
202 }
else if (smoother ==
"dilu") {
204 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
205 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
207 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
208 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
210 }
else if (smoother ==
"jac") {
211 using SeqSmoother = SeqJac<M, V, V>;
212 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
213 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
215 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
216 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
218 }
else if (smoother ==
"gs") {
219 using SeqSmoother = SeqGS<M, V, V>;
220 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
221 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
223 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
224 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
226 }
else if (smoother ==
"sor") {
227 using SeqSmoother = SeqSOR<M, V, V>;
228 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
229 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
231 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
232 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
234 }
else if (smoother ==
"ssor") {
235 using SeqSmoother = SeqSSOR<M, V, V>;
236 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
237 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
239 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
240 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
242 }
else if (smoother ==
"ilun") {
243 using SeqSmoother = SeqILU<M, V, V>;
244 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
245 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
247 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
248 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
251 OPM_THROW(std::invalid_argument,
"Properties: No smoother with name " + smoother +
".");
255 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1
256 && std::is_same_v<HYPRE_Real, typename V::field_type>) {
258 "hypre", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
259 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, C>>(op.getmat(), prm, comm);
269 const std::function<V()> weightsCalculator,
270 std::size_t pressureIndex,
272 assert(weightsCalculator);
273 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
274 OPM_THROW(std::logic_error,
275 "Pressure index out of bounds. It needs to specified for CPR");
277 using Scalar =
typename V::field_type;
279 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
280 op, prm, weightsCalculator, pressureIndex, comm);
282 F::addCreator(
"cprt",
285 const std::function<V()> weightsCalculator,
286 std::size_t pressureIndex,
288 assert(weightsCalculator);
289 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
290 OPM_THROW(std::logic_error,
291 "Pressure index out of bounds. It needs to specified for CPR");
293 using Scalar =
typename V::field_type;
295 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
296 op, prm, weightsCalculator, pressureIndex, comm);
304 if constexpr (std::is_same_v<O, WellModelGhostLastMatrixAdapter<M, V, V, true>>) {
305 F::addCreator(
"cprw",
308 const std::function<V()> weightsCalculator,
309 std::size_t pressureIndex,
311 assert(weightsCalculator);
312 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
313 OPM_THROW(std::logic_error,
314 "Pressure index out of bounds. It needs to specified for CPR");
316 using Scalar =
typename V::field_type;
318 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
319 op, prm, weightsCalculator, pressureIndex, comm);
329 F::addCreator(
"gpuilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
330 const double w = prm.get<
double>(
"relaxation", 1.0);
331 using field_type =
typename V::field_type;
332 using GpuILU0 =
typename gpuistl::
334 auto gpuILU0 = std::make_shared<GpuILU0>(op.getmat(), w);
336 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(gpuILU0);
337 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
341 F::addCreator(
"gpujac", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
342 const double w = prm.get<
double>(
"relaxation", 1.0);
343 using field_type =
typename V::field_type;
350 auto gpuJac = std::make_shared<MatrixOwner>(op.getmat(), w);
352 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(gpuJac);
353 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
357 F::addCreator(
"gpudilu", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
358 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
359 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
360 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
361 const bool reorder = prm.get<
bool>(
"reorder",
true);
362 using field_type =
typename V::field_type;
369 auto gpuDILU = std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder);
371 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(gpuDILU);
372 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
376 F::addCreator(
"opmgpuilu0", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
377 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
378 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
379 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
380 using field_type =
typename V::field_type;
388 auto gpuilu0 = std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme);
390 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(gpuilu0);
391 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
399 createParILU(
const Operator& op,
const PropertyTree& prm,
const Comm& comm,
const int ilulevel)
402 using M =
typename F::Matrix;
403 using V =
typename F::Vector;
405 const double w = prm.
get<
double>(
"relaxation", 1.0);
406 const bool redblack = prm.
get<
bool>(
"redblack",
false);
407 const bool reorder_spheres = prm.
get<
bool>(
"reorder_spheres",
false);
411 assert(num_interior <= op.getmat().N());
412 return std::make_shared<ParallelOverlappingILU0<M, V, V, Comm>>(
413 op.getmat(), comm, w,
MILU_VARIANT::ILU, num_interior, redblack, reorder_spheres);
415 return std::make_shared<ParallelOverlappingILU0<M, V, V, Comm>>(
426 std::size_t interior_count = 0;
427 std::size_t highest_interior_index = 0;
428 const auto& is = comm.indexSet();
429 for (
const auto& ind : is) {
430 if (Comm::OwnerSet::contains(ind.local().attribute())) {
432 highest_interior_index = std::max(highest_interior_index, ind.local().local());
435 if (highest_interior_index + 1 == interior_count) {
436 return interior_count;