Go to the documentation of this file.
21#ifndef OPM_STANDARDPRECONDITIONERS_SERIAL_HPP
22#define OPM_STANDARDPRECONDITIONERS_SERIAL_HPP
26#include <opm/simulators/linalg/gpuistl_hip/PreconditionerCPUMatrixToGPUMatrix.hpp>
34template < class X, class Y>
41 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& y) override {};
42 virtual void post ([[maybe_unused]] X& x) override {};
43 virtual void apply ([[maybe_unused]] X& x, [[maybe_unused]] const Y& y) override {};
44 virtual Dune::SolverCategory::Category category() const override { return Dune::SolverCategory::sequential; };
48template < class Operator>
55 using C = Dune::Amg::SequentialInformation;
57 using M = typename F::Matrix;
58 using V = typename F::Vector;
60 F::addCreator( "ilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
61 const double w = prm.get< double>( "relaxation", 1.0);
62 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
65 F::addCreator( "duneilu", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
66 const double w = prm.get< double>( "relaxation", 1.0);
67 const int n = prm.get< int>( "ilulevel", 0);
68 const bool resort = prm.get< bool>( "resort", false);
69 return getRebuildOnUpdateWrapper<Dune::SeqILU<M, V, V>>(op.getmat(), n, w, resort);
71 F::addCreator( "paroverilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
72 const double w = prm.get< double>( "relaxation", 1.0);
73 const int n = prm.get< int>( "ilulevel", 0);
74 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
77 F::addCreator( "ilun", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
78 const int n = prm.get< int>( "ilulevel", 0);
79 const double w = prm.get< double>( "relaxation", 1.0);
80 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
83 F::addCreator( "dilu", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
84 DUNE_UNUSED_PARAMETER(prm);
85 return std::make_shared<MultithreadDILU<M, V, V>>(op.getmat());
87 F::addCreator( "mixed-ilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
88 DUNE_UNUSED_PARAMETER(prm);
89 DUNE_UNUSED_PARAMETER(op);
90 return std::make_shared<TrivialPreconditioner<V,V>>();
92 F::addCreator( "mixed-dilu", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
93 DUNE_UNUSED_PARAMETER(prm);
94 DUNE_UNUSED_PARAMETER(op);
95 return std::make_shared<TrivialPreconditioner<V,V>>();
97 F::addCreator( "jac", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
98 const int n = prm.get< int>( "repeats", 1);
99 const double w = prm.get< double>( "relaxation", 1.0);
100 return getDummyUpdateWrapper<SeqJac<M, V, V>>(op.getmat(), n, w);
102 F::addCreator( "gs", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
103 const int n = prm.get< int>( "repeats", 1);
104 const double w = prm.get< double>( "relaxation", 1.0);
105 return getDummyUpdateWrapper<SeqGS<M, V, V>>(op.getmat(), n, w);
107 F::addCreator( "sor", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
108 const int n = prm.get< int>( "repeats", 1);
109 const double w = prm.get< double>( "relaxation", 1.0);
110 return getDummyUpdateWrapper<SeqSOR<M, V, V>>(op.getmat(), n, w);
112 F::addCreator( "ssor", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
113 const int n = prm.get< int>( "repeats", 1);
114 const double w = prm.get< double>( "relaxation", 1.0);
115 return getDummyUpdateWrapper<SeqSSOR<M, V, V>>(op.getmat(), n, w);
120 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
121 F::addCreator( "amg", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
122 std::string smoother = prm.get<std::string>( "smoother", "paroverilu0");
124 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
125 if (smoother == "ilu0" || smoother == "paroverilu0") {
126 using Smoother = SeqILU<M, V, V>;
128 } else if (smoother == "jac") {
129 using Smoother = SeqJac<M, V, V>;
131 } else if (smoother == "gs") {
132 using Smoother = SeqGS<M, V, V>;
134 } else if (smoother == "dilu") {
137 } else if (smoother == "sor") {
138 using Smoother = SeqSOR<M, V, V>;
140 } else if (smoother == "ssor") {
141 using Smoother = SeqSSOR<M, V, V>;
143 } else if (smoother == "ilun") {
144 using Smoother = SeqILU<M, V, V>;
147 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
150 F::addCreator( "kamg", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
151 std::string smoother = prm.get<std::string>( "smoother", "paroverilu0");
153 std::transform(smoother.begin(), smoother.end(), smoother.begin(), ::tolower);
154 if (smoother == "ilu0" || smoother == "paroverilu0") {
155 using Smoother = SeqILU<M, V, V>;
157 } else if (smoother == "jac") {
158 using Smoother = SeqJac<M, V, V>;
160 } else if (smoother == "sor") {
161 using Smoother = SeqSOR<M, V, V>;
163 } else if (smoother == "gs") {
164 using Smoother = SeqGS<M, V, V>;
166 } else if (smoother == "ssor") {
167 using Smoother = SeqSSOR<M, V, V>;
169 } else if (smoother == "ilun") {
170 using Smoother = SeqILU<M, V, V>;
173 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
176 F::addCreator( "famg", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
177 if constexpr (std::is_same_v<typename V::field_type, float>) {
178 OPM_THROW(std::logic_error, "famg requires UMFPack which is not available for floats");
182 Dune::Amg::Parameters parms;
183 parms.setNoPreSmoothSteps(1);
184 parms.setNoPostSmoothSteps(1);
185 return getRebuildOnUpdateWrapper<Dune::Amg::FastAMG<O, V>>(op, crit, parms);
191 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1) {
192 F::addCreator( "amgx", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
194 prm_copy.put( "setup_frequency", Opm::Parameters::Get<Opm::Parameters::CprReuseInterval>());
195 return std::make_shared<Amgx::AmgxPreconditioner<M, V, V>>(op.getmat(), prm_copy);
202 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1 &&
203 std::is_same_v<HYPRE_Real, typename V::field_type>) {
204 F::addCreator( "hypre", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
205 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, Dune::Amg::SequentialInformation>>(op.getmat(), prm, Dune::Amg::SequentialInformation());
215 if constexpr (std::is_same_v<O, WellModelMatrixAdapter<M, V, V>>) {
218 []( const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
219 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
220 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
222 using Scalar = typename V::field_type;
223 using LevelTransferPolicy
225 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
226 op, prm, weightsCalculator, pressureIndex);
232 []( const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
233 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
234 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
236 using Scalar = typename V::field_type;
238 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
239 op, prm, weightsCalculator, pressureIndex);
243 []( const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
244 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
245 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
247 using Scalar = typename V::field_type;
249 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
250 op, prm, weightsCalculator, pressureIndex);
259 F::addCreator( "gpuilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
260 const double w = prm.get< double>( "relaxation", 1.0);
261 using field_type = typename V::field_type;
262 using GpuILU0 = typename gpuistl::
264 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(
265 std::make_shared<GpuILU0>(op.getmat(), w));
268 F::addCreator( "gpuilu0float", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
269 const double w = prm.get< double>( "relaxation", 1.0);
270 using block_type = typename V::block_type;
271 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
272 using matrix_type_to =
273 typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
274 using GpuILU0 = typename gpuistl::
278 auto converted = std::make_shared<Converter>(op.getmat());
279 auto adapted = std::make_shared<Adapter>(std::make_shared<GpuILU0>(converted->getConvertedMatrix(), w));
280 converted->setUnderlyingPreconditioner(adapted);
284 F::addCreator( "gpujac", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
285 const double w = prm.get< double>( "relaxation", 1.0);
286 using field_type = typename V::field_type;
293 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
294 std::make_shared<MatrixOwner>(op.getmat(), w));
297 F::addCreator( "opmgpuilu0", []( const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
298 const bool split_matrix = prm.get< bool>( "split_matrix", true);
299 const bool tune_gpu_kernels = prm.get< bool>( "tune_gpu_kernels", true);
300 const int mixed_precision_scheme = prm.get< int>( "mixed_precision_scheme", 0);
302 using field_type = typename V::field_type;
306 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
309 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme));
313 F::addCreator( "gpudilu", []( const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
314 const bool split_matrix = prm.get< bool>( "split_matrix", true);
315 const bool tune_gpu_kernels = prm.get< bool>( "tune_gpu_kernels", true);
316 const int mixed_precision_scheme = prm.get< int>( "mixed_precision_scheme", 0);
317 const bool reorder = prm.get< bool>( "reorder", true);
318 using field_type = typename V::field_type;
322 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
325 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
328 F::addCreator( "gpudilufloat", []( const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
329 const bool split_matrix = prm.get< bool>( "split_matrix", true);
330 const bool tune_gpu_kernels = prm.get< bool>( "tune_gpu_kernels", true);
331 const int mixed_precision_scheme = prm.get< int>( "mixed_precision_scheme", 0);
332 const bool reorder = prm.get< bool>( "reorder", true);
334 using block_type = typename V::block_type;
335 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
336 using matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
344 auto converted = std::make_shared<Converter>(op.getmat());
347 auto adapted = std::make_shared<Adapter>(std::make_shared<MatrixOwner>(
348 converted->getConvertedMatrix(), converted->getConvertedMatrix(),
349 split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
350 converted->setUnderlyingPreconditioner(adapted);
The OpenMP thread parallelized DILU preconditioner. Definition: DILU.hpp:53
Interface class adding the update() method to the preconditioner interface. Definition: PreconditionerWithUpdate.hpp:32
Definition: PreconditionerFactory.hpp:64
Definition: PressureBhpTransferPolicy.hpp:99
Definition: PressureTransferPolicy.hpp:55
Hierarchical collection of key/value pairs. Definition: PropertyTree.hpp:39
Definition: StandardPreconditioners_serial.hpp:36
virtual bool hasPerfectUpdate() const override Definition: StandardPreconditioners_serial.hpp:40
virtual void pre(X &x, Y &y) override Definition: StandardPreconditioners_serial.hpp:41
virtual void post(X &x) override Definition: StandardPreconditioners_serial.hpp:42
virtual void update() override Definition: StandardPreconditioners_serial.hpp:39
virtual void apply(X &x, const Y &y) override Definition: StandardPreconditioners_serial.hpp:43
TrivialPreconditioner() Definition: StandardPreconditioners_serial.hpp:38
virtual Dune::SolverCategory::Category category() const override Definition: StandardPreconditioners_serial.hpp:44
DILU preconditioner on the GPU. Definition: GpuDILU.hpp:53
Jacobi preconditioner on the GPU. Definition: GpuJac.hpp:47
ILU0 preconditioner on the GPU. Definition: OpmGpuILU0.hpp:51
Makes a CUDA preconditioner available to a CPU simulator. Definition: PreconditionerAdapter.hpp:43
Convert a CPU matrix to a GPU matrix and use a CUDA preconditioner on the GPU. Definition: PreconditionerCPUMatrixToGPUMatrix.hpp:42
Converts the field type (eg. double to float) to benchmark single precision preconditioners. Definition: PreconditionerConvertFieldTypeAdapter.hpp:86
Definition: fvbaseprimaryvariables.hh:161
Definition: blackoilbioeffectsmodules.hh:43
@ ILU Do not perform modified ILU.
Definition: PreconditionerFactory.hpp:43
static Criterion criterion(const PropertyTree &prm) Definition: StandardPreconditioners_mpi.hpp:87
Definition: StandardPreconditioners_mpi.hpp:132
|