138 using namespace Dune;
142 using M =
typename F::Matrix;
143 using V =
typename F::Vector;
145 F::addCreator(
"ilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
146 return createParILU(op, prm, comm, 0);
148 F::addCreator(
"paroverilu0",
149 [](
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(
"ilun", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
153 return createParILU(op, prm, comm, prm.get<
int>(
"ilulevel", 0));
155 F::addCreator(
"duneilu", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
156 const int n = prm.get<
int>(
"ilulevel", 0);
157 const double w = prm.get<
double>(
"relaxation", 1.0);
158 const bool resort = prm.get<
bool>(
"resort",
false);
159 return wrapBlockPreconditioner<RebuildOnUpdatePreconditioner<Dune::SeqILU<M, V, V>>>(
160 comm, std::cref(op.getmat()), n, w, resort);
162 F::addCreator(
"dilu", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
163 DUNE_UNUSED_PARAMETER(prm);
164 return wrapBlockPreconditioner<MultithreadDILU<M, V, V>>(comm, op.getmat());
166 F::addCreator(
"jac", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
167 const int n = prm.get<
int>(
"repeats", 1);
168 const double w = prm.get<
double>(
"relaxation", 1.0);
169 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqJac<M, V, V>>>(comm, op.getmat(), n, w);
171 F::addCreator(
"gs", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
172 const int n = prm.get<
int>(
"repeats", 1);
173 const double w = prm.get<
double>(
"relaxation", 1.0);
174 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqGS<M, V, V>>>(comm, op.getmat(), n, w);
176 F::addCreator(
"sor", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
177 const int n = prm.get<
int>(
"repeats", 1);
178 const double w = prm.get<
double>(
"relaxation", 1.0);
179 return wrapBlockPreconditioner<DummyUpdatePreconditioner<SeqSOR<M, V, V>>>(comm, op.getmat(), n, w);
181 F::addCreator(
"ssor", [](
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<SeqSSOR<M, V, V>>>(comm, op.getmat(), n, w);
191 if constexpr (std::is_same_v<O, Dune::OverlappingSchwarzOperator<M, V, V, C>> ||
192 std::is_same_v<O, Opm::GhostLastMatrixAdapter<M, V, V, C>>) {
193 F::addCreator(
"amg", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
194 using PrecPtr = std::shared_ptr<Dune::PreconditionerWithUpdate<V, V>>;
195 std::string smoother = prm.get<std::string>(
"smoother",
"paroverilu0");
197 std::ranges::transform(smoother, smoother.begin(), ::tolower);
199 if (smoother ==
"ilu0" || smoother ==
"paroverilu0") {
201 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
202 auto sargs = AMGSmootherArgsHelper<Smoother>::args(prm);
203 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
205 }
else if (smoother ==
"dilu") {
207 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
208 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
210 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
211 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
213 }
else if (smoother ==
"jac") {
214 using SeqSmoother = SeqJac<M, V, V>;
215 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
216 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
218 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
219 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
221 }
else if (smoother ==
"gs") {
222 using SeqSmoother = SeqGS<M, V, V>;
223 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
224 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
226 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
227 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
229 }
else if (smoother ==
"sor") {
230 using SeqSmoother = SeqSOR<M, V, V>;
231 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
232 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
234 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
235 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
237 }
else if (smoother ==
"ssor") {
238 using SeqSmoother = SeqSSOR<M, V, V>;
239 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
240 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
242 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
243 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
245 }
else if (smoother ==
"ilun") {
246 using SeqSmoother = SeqILU<M, V, V>;
247 using Smoother = Dune::BlockPreconditioner<V, V, C, SeqSmoother>;
248 using SmootherArgs =
typename Dune::Amg::SmootherTraits<Smoother>::Arguments;
250 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
251 PrecPtr prec = std::make_shared<Dune::Amg::AMGCPR<O, V, Smoother, C>>(op, crit, sargs, comm);
254 OPM_THROW(std::invalid_argument,
"Properties: No smoother with name " + smoother +
".");
258 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1
259 && std::is_same_v<HYPRE_Real, typename V::field_type>) {
261 "hypre", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
262 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, C>>(op.getmat(), prm, comm);
272 const std::function<V()> weightsCalculator,
273 std::size_t pressureIndex,
275 assert(weightsCalculator);
276 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
277 OPM_THROW(std::logic_error,
278 "Pressure index out of bounds. It needs to specified for CPR");
280 using Scalar =
typename V::field_type;
282 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
283 op, prm, weightsCalculator, pressureIndex, comm);
285 F::addCreator(
"cprt",
288 const std::function<V()> weightsCalculator,
289 std::size_t pressureIndex,
291 assert(weightsCalculator);
292 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
293 OPM_THROW(std::logic_error,
294 "Pressure index out of bounds. It needs to specified for CPR");
296 using Scalar =
typename V::field_type;
298 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
299 op, prm, weightsCalculator, pressureIndex, comm);
307 if constexpr (std::is_same_v<O, WellModelGhostLastMatrixAdapter<M, V, V, true>>) {
308 F::addCreator(
"cprw",
311 const std::function<V()> weightsCalculator,
312 std::size_t pressureIndex,
314 assert(weightsCalculator);
315 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
316 OPM_THROW(std::logic_error,
317 "Pressure index out of bounds. It needs to specified for CPR");
319 using Scalar =
typename V::field_type;
321 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy, Comm>>(
322 op, prm, weightsCalculator, pressureIndex, comm);
332 F::addCreator(
"gpuilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
333 const double w = prm.get<
double>(
"relaxation", 1.0);
334 using field_type =
typename V::field_type;
335 using GpuILU0 =
typename gpuistl::
337 auto gpuILU0 = std::make_shared<GpuILU0>(op.getmat(), w);
339 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(gpuILU0);
340 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
344 F::addCreator(
"gpujac", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
345 const double w = prm.get<
double>(
"relaxation", 1.0);
346 using field_type =
typename V::field_type;
353 auto gpuJac = std::make_shared<MatrixOwner>(op.getmat(), w);
355 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(gpuJac);
356 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
360 F::addCreator(
"gpudilu", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
361 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
362 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
363 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
364 const bool reorder = prm.get<
bool>(
"reorder",
true);
365 using field_type =
typename V::field_type;
372 auto gpuDILU = std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder);
374 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(gpuDILU);
375 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
379 F::addCreator(
"opmgpuilu0", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t,
const C& comm) {
380 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
381 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
382 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
383 using field_type =
typename V::field_type;
391 auto gpuilu0 = std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme);
393 auto adapted = std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(gpuilu0);
394 auto wrapped = std::make_shared<gpuistl::GpuBlockPreconditioner<V, V, Comm>>(adapted, comm);
402 createParILU(
const Operator& op,
const PropertyTree& prm,
const Comm& comm,
const int ilulevel)
405 using M =
typename F::Matrix;
406 using V =
typename F::Vector;
408 const double w = prm.
get<
double>(
"relaxation", 1.0);
409 const bool redblack = prm.
get<
bool>(
"redblack",
false);
410 const bool reorder_spheres = prm.
get<
bool>(
"reorder_spheres",
false);
414 assert(num_interior <= op.getmat().N());
415 return std::make_shared<ParallelOverlappingILU0<M, V, V, Comm>>(
416 op.getmat(), comm, w,
MILU_VARIANT::ILU, num_interior, redblack, reorder_spheres);
418 return std::make_shared<ParallelOverlappingILU0<M, V, V, Comm>>(
429 std::size_t interior_count = 0;
430 std::size_t highest_interior_index = 0;
431 const auto& is = comm.indexSet();
432 for (
const auto& ind : is) {
433 if (Comm::OwnerSet::contains(ind.local().attribute())) {
435 highest_interior_index = std::max(highest_interior_index, ind.local().local());
438 if (highest_interior_index + 1 == interior_count) {
439 return interior_count;