/*
Copyright 2009, 2010 SINTEF ICT, Applied Mathematics.
Copyright 2019 SINTEF Digital, Mathematics and Cybernetics.
This file is part of the Open Porous Media project (OPM).
OPM is free software: you can redistribute it and/or modify
it under the terms of the GNU General Public License as published by
the Free Software Foundation, either version 3 of the License, or
(at your option) any later version.
OPM is distributed in the hope that it will be useful,
but WITHOUT ANY WARRANTY; without even the implied warranty of
MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
GNU General Public License for more details.
You should have received a copy of the GNU General Public License
along with OPM. If not, see .
*/
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
#include
// Include all cuistl/GPU preconditioners inside of this headerfile
#include
namespace Opm {
template
struct AMGSmootherArgsHelper
{
static auto args(const PropertyTree& prm)
{
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs smootherArgs;
smootherArgs.iterations = prm.get("iterations", 1);
// smootherArgs.overlap=SmootherArgs::vertex;
// smootherArgs.overlap=SmootherArgs::none;
// smootherArgs.overlap=SmootherArgs::aggregate;
smootherArgs.relaxationFactor = prm.get("relaxation", 1.0);
return smootherArgs;
}
};
template
struct AMGSmootherArgsHelper>
{
static auto args(const PropertyTree& prm)
{
using Smoother = ParallelOverlappingILU0;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs smootherArgs;
smootherArgs.iterations = prm.get("iterations", 1);
const int iluwitdh = prm.get("iluwidth", 0);
smootherArgs.setN(iluwitdh);
const MILU_VARIANT milu = convertString2Milu(prm.get("milutype", std::string("ilu")));
smootherArgs.setMilu(milu);
// smootherArgs.overlap=SmootherArgs::vertex;
// smootherArgs.overlap=SmootherArgs::none;
// smootherArgs.overlap=SmootherArgs::aggregate;
smootherArgs.relaxationFactor = prm.get("relaxation", 1.0);
return smootherArgs;
}
};
// trailing return type with decltype used for detecting existence of setUseFixedOrder member function by overloading the setUseFixedOrder function
template
auto setUseFixedOrder(C criterion, bool booleanValue) -> decltype(criterion.setUseFixedOrder(booleanValue))
{
return criterion.setUseFixedOrder(booleanValue); // Set flag to ensure that the matrices in the AMG hierarchy are constructed with deterministic indices.
}
template
void setUseFixedOrder(C, ...)
{
// do nothing, since the function setUseFixedOrder does not exist yet
}
template
typename AMGHelper::Criterion
AMGHelper::criterion(const PropertyTree& prm)
{
Criterion criterion(15, prm.get("coarsenTarget", 1200));
criterion.setDefaultValuesIsotropic(2);
criterion.setAlpha(prm.get("alpha", 0.33));
criterion.setBeta(prm.get("beta", 1e-5));
criterion.setMaxLevel(prm.get("maxlevel", 15));
criterion.setSkipIsolated(prm.get("skip_isolated", false));
criterion.setNoPreSmoothSteps(prm.get("pre_smooth", 1));
criterion.setNoPostSmoothSteps(prm.get("post_smooth", 1));
criterion.setDebugLevel(prm.get("verbosity", 0));
// As the default we request to accumulate data to 1 process always as our matrix
// graph might be unsymmetric and hence not supported by the PTScotch/ParMetis
// calls in DUNE. Accumulating to 1 skips PTScotch/ParMetis
criterion.setAccumulate(static_cast(prm.get("accumulate", 1)));
criterion.setProlongationDampingFactor(prm.get("prolongationdamping", 1.6));
criterion.setMaxDistance(prm.get("maxdistance", 2));
criterion.setMaxConnectivity(prm.get("maxconnectivity", 15));
criterion.setMaxAggregateSize(prm.get("maxaggsize", 6));
criterion.setMinAggregateSize(prm.get("minaggsize", 4));
setUseFixedOrder(criterion, true); // If possible, set flag to ensure that the matrices in the AMG hierarchy are constructed with deterministic indices.
return criterion;
}
template
template
typename AMGHelper::PrecPtr
AMGHelper::makeAmgPreconditioner(const Operator& op,
const PropertyTree& prm,
bool useKamg)
{
auto crit = criterion(prm);
auto sargs = AMGSmootherArgsHelper::args(prm);
if (useKamg) {
using Type = Dune::DummyUpdatePreconditioner>;
return std::make_shared(
op, crit, sargs, prm.get("max_krylov", 1), prm.get("min_reduction", 1e-1));
} else {
using Type = Dune::Amg::AMGCPR;
return std::make_shared(op, crit, sargs);
}
}
template
struct StandardPreconditioners {
static void add()
{
using namespace Dune;
using O = Operator;
using C = Comm;
using F = PreconditionerFactory;
using M = typename F::Matrix;
using V = typename F::Vector;
using P = PropertyTree;
F::addCreator("ILU0", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
return createParILU(op, prm, comm, 0);
});
F::addCreator("ParOverILU0",
[](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
return createParILU(op, prm, comm, prm.get("ilulevel", 0));
});
F::addCreator("ILUn", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
return createParILU(op, prm, comm, prm.get("ilulevel", 0));
});
F::addCreator("DuneILU", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const int n = prm.get("ilulevel", 0);
const double w = prm.get("relaxation", 1.0);
const bool resort = prm.get("resort", false);
return wrapBlockPreconditioner>>(
comm, op.getmat(), n, w, resort);
});
F::addCreator("DILU", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
DUNE_UNUSED_PARAMETER(prm);
return wrapBlockPreconditioner>(comm, op.getmat());
});
F::addCreator("Jac", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return wrapBlockPreconditioner>>(comm, op.getmat(), n, w);
});
F::addCreator("GS", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return wrapBlockPreconditioner>>(comm, op.getmat(), n, w);
});
F::addCreator("SOR", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return wrapBlockPreconditioner>>(comm, op.getmat(), n, w);
});
F::addCreator("SSOR", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return wrapBlockPreconditioner>>(comm, op.getmat(), n, w);
});
// Only add AMG preconditioners to the factory if the operator
// is the overlapping schwarz operator. This could be extended
// later, but at this point no other operators are compatible
// with the AMG hierarchy construction.
if constexpr (std::is_same_v> ||
std::is_same_v>) {
F::addCreator("amg", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
using PrecPtr = std::shared_ptr>;
const std::string smoother = prm.get("smoother", "ParOverILU0");
// TODO: merge this with ILUn, and possibly simplify the factory to only work with ILU?
if (smoother == "ILU0" || smoother == "ParOverILU0") {
using Smoother = ParallelOverlappingILU0;
auto crit = AMGHelper::criterion(prm);
auto sargs = AMGSmootherArgsHelper::args(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else if (smoother == "DILU") {
using SeqSmoother = Dune::MultithreadDILU;
using Smoother = Dune::BlockPreconditioner;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs sargs;
auto crit = AMGHelper::criterion(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else if (smoother == "Jac") {
using SeqSmoother = SeqJac;
using Smoother = Dune::BlockPreconditioner;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs sargs;
auto crit = AMGHelper::criterion(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else if (smoother == "GS") {
using SeqSmoother = SeqGS;
using Smoother = Dune::BlockPreconditioner;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs sargs;
auto crit = AMGHelper::criterion(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else if (smoother == "SOR") {
using SeqSmoother = SeqSOR;
using Smoother = Dune::BlockPreconditioner;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs sargs;
auto crit = AMGHelper::criterion(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else if (smoother == "SSOR") {
using SeqSmoother = SeqSSOR;
using Smoother = Dune::BlockPreconditioner;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs sargs;
auto crit = AMGHelper::criterion(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else if (smoother == "ILUn") {
using SeqSmoother = SeqILU;
using Smoother = Dune::BlockPreconditioner;
using SmootherArgs = typename Dune::Amg::SmootherTraits::Arguments;
SmootherArgs sargs;
auto crit = AMGHelper::criterion(prm);
PrecPtr prec = std::make_shared>(op, crit, sargs, comm);
return prec;
} else {
OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
}
});
}
F::addCreator("cpr",
[](const O& op,
const P& prm,
const std::function weightsCalculator,
std::size_t pressureIndex,
const C& comm) {
assert(weightsCalculator);
if (pressureIndex == std::numeric_limits::max()) {
OPM_THROW(std::logic_error,
"Pressure index out of bounds. It needs to specified for CPR");
}
using Scalar = typename V::field_type;
using LevelTransferPolicy = PressureTransferPolicy;
return std::make_shared>(
op, prm, weightsCalculator, pressureIndex, comm);
});
F::addCreator("cprt",
[](const O& op,
const P& prm,
const std::function weightsCalculator,
std::size_t pressureIndex,
const C& comm) {
assert(weightsCalculator);
if (pressureIndex == std::numeric_limits::max()) {
OPM_THROW(std::logic_error,
"Pressure index out of bounds. It needs to specified for CPR");
}
using Scalar = typename V::field_type;
using LevelTransferPolicy = PressureTransferPolicy;
return std::make_shared>(
op, prm, weightsCalculator, pressureIndex, comm);
});
if constexpr (std::is_same_v>) {
F::addCreator("cprw",
[](const O& op,
const P& prm,
const std::function weightsCalculator,
std::size_t pressureIndex,
const C& comm) {
assert(weightsCalculator);
if (pressureIndex == std::numeric_limits::max()) {
OPM_THROW(std::logic_error,
"Pressure index out of bounds. It needs to specified for CPR");
}
using Scalar = typename V::field_type;
using LevelTransferPolicy = PressureBhpTransferPolicy;
return std::make_shared>(
op, prm, weightsCalculator, pressureIndex, comm);
});
}
#if HAVE_CUDA
F::addCreator("GPUILU0", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const double w = prm.get("relaxation", 1.0);
using field_type = typename V::field_type;
using GpuILU0 = typename gpuistl::
GpuSeqILU0, gpuistl::GpuVector>;
auto gpuILU0 = std::make_shared(op.getmat(), w);
auto adapted = std::make_shared>(gpuILU0);
auto wrapped = std::make_shared>(adapted, comm);
return wrapped;
});
F::addCreator("GPUJAC", [](const O& op, const P& prm, const std::function&, std::size_t, const C& comm) {
const double w = prm.get("relaxation", 1.0);
using field_type = typename V::field_type;
using GpuJac =
typename gpuistl::GpuJac, gpuistl::GpuVector>;
auto gpuJac = std::make_shared(op.getmat(), w);
auto adapted = std::make_shared>(gpuJac);
auto wrapped = std::make_shared>(adapted, comm);
return wrapped;
});
F::addCreator("GPUDILU", [](const O& op, [[maybe_unused]] const P& prm, const std::function&, std::size_t, const C& comm) {
const bool split_matrix = prm.get("split_matrix", true);
const bool tune_gpu_kernels = prm.get("tune_gpu_kernels", true);
using field_type = typename V::field_type;
using GpuDILU = typename gpuistl::GpuDILU, gpuistl::GpuVector>;
auto gpuDILU = std::make_shared(op.getmat(), split_matrix, tune_gpu_kernels);
auto adapted = std::make_shared>(gpuDILU);
auto wrapped = std::make_shared>(adapted, comm);
return wrapped;
});
F::addCreator("OPMCUILU0", [](const O& op, [[maybe_unused]] const P& prm, const std::function&, std::size_t, const C& comm) {
const bool split_matrix = prm.get("split_matrix", true);
const bool tune_gpu_kernels = prm.get("tune_gpu_kernels", true);
using field_type = typename V::field_type;
using OpmCuILU0 = typename gpuistl::OpmCuILU0, gpuistl::GpuVector>;
auto cuilu0 = std::make_shared(op.getmat(), split_matrix, tune_gpu_kernels);
auto adapted = std::make_shared>(cuilu0);
auto wrapped = std::make_shared>(adapted, comm);
return wrapped;
});
#endif
}
static typename PreconditionerFactory::PrecPtr
createParILU(const Operator& op, const PropertyTree& prm, const Comm& comm, const int ilulevel)
{
using F = PreconditionerFactory;
using M = typename F::Matrix;
using V = typename F::Vector;
const double w = prm.get("relaxation", 1.0);
const bool redblack = prm.get("redblack", false);
const bool reorder_spheres = prm.get("reorder_spheres", false);
// Already a parallel preconditioner. Need to pass comm, but no need to wrap it in a BlockPreconditioner.
if (ilulevel == 0) {
const std::size_t num_interior = interiorIfGhostLast(comm);
return std::make_shared>(
op.getmat(), comm, w, MILU_VARIANT::ILU, num_interior, redblack, reorder_spheres);
} else {
return std::make_shared>(
op.getmat(), comm, ilulevel, w, MILU_VARIANT::ILU, redblack, reorder_spheres);
}
}
/// Helper method to determine if the local partitioning has the
/// K interior cells from [0, K-1] and ghost cells from [K, N-1].
/// Returns K if true, otherwise returns N. This is motivated by
/// usage in the ParallelOverlappingILU0 preconditioner.
static std::size_t interiorIfGhostLast(const Comm& comm)
{
std::size_t interior_count = 0;
std::size_t highest_interior_index = 0;
const auto& is = comm.indexSet();
for (const auto& ind : is) {
if (Comm::OwnerSet::contains(ind.local().attribute())) {
++interior_count;
highest_interior_index = std::max(highest_interior_index, ind.local().local());
}
}
if (highest_interior_index + 1 == interior_count) {
return interior_count;
} else {
return is.size();
}
}
};
template
struct StandardPreconditioners {
static void add()
{
using namespace Dune;
using O = Operator;
using C = Dune::Amg::SequentialInformation;
using F = PreconditionerFactory;
using M = typename F::Matrix;
using V = typename F::Vector;
using P = PropertyTree;
F::addCreator("ILU0", [](const O& op, const P& prm, const std::function&, std::size_t) {
const double w = prm.get("relaxation", 1.0);
return std::make_shared>(
op.getmat(), 0, w, MILU_VARIANT::ILU);
});
F::addCreator("DuneILU", [](const O& op, const P& prm, const std::function&, std::size_t) {
const double w = prm.get("relaxation", 1.0);
const int n = prm.get("ilulevel", 0);
const bool resort = prm.get("resort", false);
return getRebuildOnUpdateWrapper>(op.getmat(), n, w, resort);
});
F::addCreator("ParOverILU0", [](const O& op, const P& prm, const std::function&, std::size_t) {
const double w = prm.get("relaxation", 1.0);
const int n = prm.get("ilulevel", 0);
return std::make_shared>(
op.getmat(), n, w, MILU_VARIANT::ILU);
});
F::addCreator("ILUn", [](const O& op, const P& prm, const std::function&, std::size_t) {
const int n = prm.get("ilulevel", 0);
const double w = prm.get("relaxation", 1.0);
return std::make_shared>(
op.getmat(), n, w, MILU_VARIANT::ILU);
});
F::addCreator("DILU", [](const O& op, const P& prm, const std::function&, std::size_t) {
DUNE_UNUSED_PARAMETER(prm);
return std::make_shared>(op.getmat());
});
F::addCreator("Jac", [](const O& op, const P& prm, const std::function&, std::size_t) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return getDummyUpdateWrapper>(op.getmat(), n, w);
});
F::addCreator("GS", [](const O& op, const P& prm, const std::function&, std::size_t) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return getDummyUpdateWrapper>(op.getmat(), n, w);
});
F::addCreator("SOR", [](const O& op, const P& prm, const std::function&, std::size_t) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return getDummyUpdateWrapper>(op.getmat(), n, w);
});
F::addCreator("SSOR", [](const O& op, const P& prm, const std::function&, std::size_t) {
const int n = prm.get("repeats", 1);
const double w = prm.get("relaxation", 1.0);
return getDummyUpdateWrapper>(op.getmat(), n, w);
});
// Only add AMG preconditioners to the factory if the operator
// is an actual matrix operator.
if constexpr (std::is_same_v>) {
F::addCreator("amg", [](const O& op, const P& prm, const std::function&, std::size_t) {
const std::string smoother = prm.get("smoother", "ParOverILU0");
if (smoother == "ILU0" || smoother == "ParOverILU0") {
using Smoother = SeqILU;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else if (smoother == "Jac") {
using Smoother = SeqJac;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else if (smoother == "GS") {
using Smoother = SeqGS;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else if (smoother == "DILU") {
using Smoother = MultithreadDILU;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else if (smoother == "SOR") {
using Smoother = SeqSOR;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else if (smoother == "SSOR") {
using Smoother = SeqSSOR;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else if (smoother == "ILUn") {
using Smoother = SeqILU;
return AMGHelper::template makeAmgPreconditioner(op, prm);
} else {
OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
}
});
F::addCreator("kamg", [](const O& op, const P& prm, const std::function&, std::size_t) {
const std::string smoother = prm.get("smoother", "ParOverILU0");
if (smoother == "ILU0" || smoother == "ParOverILU0") {
using Smoother = SeqILU;
return AMGHelper::template makeAmgPreconditioner(op, prm, true);
} else if (smoother == "Jac") {
using Smoother = SeqJac;
return AMGHelper::template makeAmgPreconditioner