36struct StandardPreconditioners<Operator, Dune::Amg::SequentialInformation, typename std::enable_if_t<!Opm::is_gpu_operator_v<Operator>>>
42 using C = Dune::Amg::SequentialInformation;
44 using M =
typename F::Matrix;
45 using V =
typename F::Vector;
47 F::addCreator(
"ilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
48 const double w = prm.get<
double>(
"relaxation", 1.0);
49 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
52 F::addCreator(
"duneilu", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
53 const double w = prm.get<
double>(
"relaxation", 1.0);
54 const int n = prm.get<
int>(
"ilulevel", 0);
55 const bool resort = prm.get<
bool>(
"resort",
false);
56 return getRebuildOnUpdateWrapper<Dune::SeqILU<M, V, V>>(op.getmat(), n, w, resort);
58 F::addCreator(
"paroverilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
59 const double w = prm.get<
double>(
"relaxation", 1.0);
60 const int n = prm.get<
int>(
"ilulevel", 0);
61 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
64 F::addCreator(
"ilun", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
65 const int n = prm.get<
int>(
"ilulevel", 0);
66 const double w = prm.get<
double>(
"relaxation", 1.0);
67 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
70 F::addCreator(
"dilu", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
71 DUNE_UNUSED_PARAMETER(prm);
72 return std::make_shared<MultithreadDILU<M, V, V>>(op.getmat());
74 F::addCreator(
"jac", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
75 const int n = prm.get<
int>(
"repeats", 1);
76 const double w = prm.get<
double>(
"relaxation", 1.0);
77 return getDummyUpdateWrapper<SeqJac<M, V, V>>(op.getmat(), n, w);
79 F::addCreator(
"gs", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
80 const int n = prm.get<
int>(
"repeats", 1);
81 const double w = prm.get<
double>(
"relaxation", 1.0);
82 return getDummyUpdateWrapper<SeqGS<M, V, V>>(op.getmat(), n, w);
84 F::addCreator(
"sor", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
85 const int n = prm.get<
int>(
"repeats", 1);
86 const double w = prm.get<
double>(
"relaxation", 1.0);
87 return getDummyUpdateWrapper<SeqSOR<M, V, V>>(op.getmat(), n, w);
89 F::addCreator(
"ssor", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
90 const int n = prm.get<
int>(
"repeats", 1);
91 const double w = prm.get<
double>(
"relaxation", 1.0);
92 return getDummyUpdateWrapper<SeqSSOR<M, V, V>>(op.getmat(), n, w);
97 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
98 F::addCreator(
"amg", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
99 std::string smoother = prm.get<std::string>(
"smoother",
"paroverilu0");
101 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
102 if (smoother ==
"ilu0" || smoother ==
"paroverilu0") {
103 using Smoother = SeqILU<M, V, V>;
105 }
else if (smoother ==
"jac") {
106 using Smoother = SeqJac<M, V, V>;
108 }
else if (smoother ==
"gs") {
109 using Smoother = SeqGS<M, V, V>;
111 }
else if (smoother ==
"dilu") {
114 }
else if (smoother ==
"sor") {
115 using Smoother = SeqSOR<M, V, V>;
117 }
else if (smoother ==
"ssor") {
118 using Smoother = SeqSSOR<M, V, V>;
120 }
else if (smoother ==
"ilun") {
121 using Smoother = SeqILU<M, V, V>;
124 OPM_THROW(std::invalid_argument,
"Properties: No smoother with name " + smoother +
".");
127 F::addCreator(
"kamg", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
128 std::string smoother = prm.get<std::string>(
"smoother",
"paroverilu0");
130 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
131 if (smoother ==
"ilu0" || smoother ==
"paroverilu0") {
132 using Smoother = SeqILU<M, V, V>;
134 }
else if (smoother ==
"jac") {
135 using Smoother = SeqJac<M, V, V>;
137 }
else if (smoother ==
"sor") {
138 using Smoother = SeqSOR<M, V, V>;
140 }
else if (smoother ==
"gs") {
141 using Smoother = SeqGS<M, V, V>;
143 }
else if (smoother ==
"ssor") {
144 using Smoother = SeqSSOR<M, V, V>;
146 }
else if (smoother ==
"ilun") {
147 using Smoother = SeqILU<M, V, V>;
150 OPM_THROW(std::invalid_argument,
"Properties: No smoother with name " + smoother +
".");
153 F::addCreator(
"famg", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
154 if constexpr (std::is_same_v<typename V::field_type, float>) {
155 OPM_THROW(std::logic_error,
"famg requires UMFPack which is not available for floats");
158 auto crit = AMGHelper<O, C, M, V>::criterion(prm);
159 Dune::Amg::Parameters parms;
160 parms.setNoPreSmoothSteps(1);
161 parms.setNoPostSmoothSteps(1);
162 return getRebuildOnUpdateWrapper<Dune::Amg::FastAMG<O, V>>(op, crit, parms);
168 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1) {
169 F::addCreator(
"amgx", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
172 return std::make_shared<Amgx::AmgxPreconditioner<M, V, V>>(op.getmat(), prm_copy);
179 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1 &&
180 std::is_same_v<HYPRE_Real, typename V::field_type>) {
181 F::addCreator(
"hypre", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
182 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, Dune::Amg::SequentialInformation>>(op.getmat(), prm, Dune::Amg::SequentialInformation());
192 if constexpr (std::is_same_v<O, WellModelMatrixAdapter<M, V, V>>) {
195 [](
const O& op,
const P& prm,
const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
196 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
197 OPM_THROW(std::logic_error,
"Pressure index out of bounds. It needs to specified for CPR");
199 using Scalar =
typename V::field_type;
200 using LevelTransferPolicy
202 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
203 op, prm, weightsCalculator, pressureIndex);
209 [](
const O& op,
const P& prm,
const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
210 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
211 OPM_THROW(std::logic_error,
"Pressure index out of bounds. It needs to specified for CPR");
213 using Scalar =
typename V::field_type;
215 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
216 op, prm, weightsCalculator, pressureIndex);
220 [](
const O& op,
const P& prm,
const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
221 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
222 OPM_THROW(std::logic_error,
"Pressure index out of bounds. It needs to specified for CPR");
224 using Scalar =
typename V::field_type;
226 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
227 op, prm, weightsCalculator, pressureIndex);
236 F::addCreator(
"gpuilu0", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
237 const double w = prm.get<
double>(
"relaxation", 1.0);
238 using field_type =
typename V::field_type;
239 using GpuILU0 =
typename gpuistl::
241 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(
242 std::make_shared<GpuILU0>(op.getmat(), w));
245 F::addCreator(
"gpuilu0float", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
246 const double w = prm.get<
double>(
"relaxation", 1.0);
247 using block_type =
typename V::block_type;
248 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
249 using matrix_type_to =
250 typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
251 using GpuILU0 =
typename gpuistl::
255 auto converted = std::make_shared<Converter>(op.getmat());
256 auto adapted = std::make_shared<Adapter>(std::make_shared<GpuILU0>(converted->getConvertedMatrix(), w));
257 converted->setUnderlyingPreconditioner(adapted);
261 F::addCreator(
"gpujac", [](
const O& op,
const P& prm,
const std::function<V()>&, std::size_t) {
262 const double w = prm.get<
double>(
"relaxation", 1.0);
263 using field_type =
typename V::field_type;
270 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
271 std::make_shared<MatrixOwner>(op.getmat(), w));
274 F::addCreator(
"opmgpuilu0", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t) {
275 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
276 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
277 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
279 using field_type =
typename V::field_type;
283 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
286 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme));
290 F::addCreator(
"gpudilu", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t) {
291 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
292 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
293 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
294 const bool reorder = prm.get<
bool>(
"reorder",
true);
295 using field_type =
typename V::field_type;
299 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
302 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
305 F::addCreator(
"gpudilufloat", [](
const O& op, [[maybe_unused]]
const P& prm,
const std::function<V()>&, std::size_t) {
306 const bool split_matrix = prm.get<
bool>(
"split_matrix",
true);
307 const bool tune_gpu_kernels = prm.get<
bool>(
"tune_gpu_kernels",
true);
308 const int mixed_precision_scheme = prm.get<
int>(
"mixed_precision_scheme", 0);
309 const bool reorder = prm.get<
bool>(
"reorder",
true);
311 using block_type =
typename V::block_type;
312 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
313 using matrix_type_to =
typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
321 auto converted = std::make_shared<Converter>(op.getmat());
324 auto adapted = std::make_shared<Adapter>(std::make_shared<MatrixOwner>(
325 converted->getConvertedMatrix(), converted->getConvertedMatrix(),
326 split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
327 converted->setUnderlyingPreconditioner(adapted);