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>
38template < class X, class Y>
45 virtual void pre ([[maybe_unused]] X& x, [[maybe_unused]] Y& y) override {};
46 virtual void post ([[maybe_unused]] X& x) override {};
47 virtual void apply ([[maybe_unused]] X& x, [[maybe_unused]] const Y& y) override {};
48 virtual Dune::SolverCategory::Category category() const override { return Dune::SolverCategory::sequential; };
52template < class Operator>
59 using C = Dune::Amg::SequentialInformation;
61 using M = typename F::Matrix;
62 using V = typename F::Vector;
64 F::addCreator( "ilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
65 const double w = prm.get< double>( "relaxation", 1.0);
66 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
69 F::addCreator( "duneilu", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
70 const double w = prm.get< double>( "relaxation", 1.0);
71 const int n = prm.get< int>( "ilulevel", 0);
72 const bool resort = prm.get< bool>( "resort", false);
73 return getRebuildOnUpdateWrapper<Dune::SeqILU<M, V, V>>(std::cref(op.getmat()), n, w, resort);
75 F::addCreator( "paroverilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
76 const double w = prm.get< double>( "relaxation", 1.0);
77 const int n = prm.get< int>( "ilulevel", 0);
78 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
81 F::addCreator( "ilun", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
82 const int n = prm.get< int>( "ilulevel", 0);
83 const double w = prm.get< double>( "relaxation", 1.0);
84 return std::make_shared<ParallelOverlappingILU0<M, V, V, C>>(
87 F::addCreator( "dilu", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
88 DUNE_UNUSED_PARAMETER(prm);
89 return std::make_shared<MultithreadDILU<M, V, V>>(op.getmat());
91 F::addCreator( "mixed-ilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
92 DUNE_UNUSED_PARAMETER(prm);
93 DUNE_UNUSED_PARAMETER(op);
94 return std::make_shared<TrivialPreconditioner<V,V>>();
96 F::addCreator( "mixed-dilu", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
97 DUNE_UNUSED_PARAMETER(prm);
98 DUNE_UNUSED_PARAMETER(op);
99 return std::make_shared<TrivialPreconditioner<V,V>>();
101 F::addCreator( "jac", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
102 const int n = prm.get< int>( "repeats", 1);
103 const double w = prm.get< double>( "relaxation", 1.0);
104 return getDummyUpdateWrapper<SeqJac<M, V, V>>(op.getmat(), n, w);
106 F::addCreator( "gs", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
107 const int n = prm.get< int>( "repeats", 1);
108 const double w = prm.get< double>( "relaxation", 1.0);
109 return getDummyUpdateWrapper<SeqGS<M, V, V>>(op.getmat(), n, w);
111 F::addCreator( "sor", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
112 const int n = prm.get< int>( "repeats", 1);
113 const double w = prm.get< double>( "relaxation", 1.0);
114 return getDummyUpdateWrapper<SeqSOR<M, V, V>>(op.getmat(), n, w);
116 F::addCreator( "ssor", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
117 const int n = prm.get< int>( "repeats", 1);
118 const double w = prm.get< double>( "relaxation", 1.0);
119 return getDummyUpdateWrapper<SeqSSOR<M, V, V>>(op.getmat(), n, w);
124 if constexpr (std::is_same_v<O, Dune::MatrixAdapter<M, V, V>>) {
125 F::addCreator( "amg", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
126 std::string smoother = prm.get<std::string>( "smoother", "paroverilu0");
128 std::ranges::transform(smoother, smoother.begin(), ::tolower);
129 if (smoother == "ilu0" || smoother == "paroverilu0") {
130 using Smoother = SeqILU<M, V, V>;
132 } else if (smoother == "jac") {
133 using Smoother = SeqJac<M, V, V>;
135 } else if (smoother == "gs") {
136 using Smoother = SeqGS<M, V, V>;
138 } else if (smoother == "dilu") {
141 } else if (smoother == "sor") {
142 using Smoother = SeqSOR<M, V, V>;
144 } else if (smoother == "ssor") {
145 using Smoother = SeqSSOR<M, V, V>;
147 } else if (smoother == "ilun") {
148 using Smoother = SeqILU<M, V, V>;
151 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
154 F::addCreator( "kamg", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
155 std::string smoother = prm.get<std::string>( "smoother", "paroverilu0");
157 std::ranges::transform(smoother, smoother.begin(), ::tolower);
158 if (smoother == "ilu0" || smoother == "paroverilu0") {
159 using Smoother = SeqILU<M, V, V>;
161 } else if (smoother == "jac") {
162 using Smoother = SeqJac<M, V, V>;
164 } else if (smoother == "sor") {
165 using Smoother = SeqSOR<M, V, V>;
167 } else if (smoother == "gs") {
168 using Smoother = SeqGS<M, V, V>;
170 } else if (smoother == "ssor") {
171 using Smoother = SeqSSOR<M, V, V>;
173 } else if (smoother == "ilun") {
174 using Smoother = SeqILU<M, V, V>;
177 OPM_THROW(std::invalid_argument, "Properties: No smoother with name " + smoother + ".");
180 F::addCreator( "famg", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
181 if constexpr (std::is_same_v<typename V::field_type, float>) {
182 OPM_THROW(std::logic_error, "famg requires UMFPack which is not available for floats");
186 Dune::Amg::Parameters parms;
187 parms.setNoPreSmoothSteps(1);
188 parms.setNoPostSmoothSteps(1);
189 return getRebuildOnUpdateWrapper<Dune::Amg::FastAMG<O, V>>(op, crit, parms);
195 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1) {
196 F::addCreator( "amgx", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
198 prm_copy.put( "setup_frequency", Opm::Parameters::Get<Opm::Parameters::CprReuseInterval>());
199 return std::make_shared<Amgx::AmgxPreconditioner<M, V, V>>(op.getmat(), prm_copy);
206 if constexpr (M::block_type::rows == 1 && M::block_type::cols == 1 &&
207 std::is_same_v<HYPRE_Real, typename V::field_type>) {
208 F::addCreator( "hypre", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
209 return std::make_shared<Hypre::HyprePreconditioner<M, V, V, Dune::Amg::SequentialInformation>>(op.getmat(), prm, Dune::Amg::SequentialInformation());
219 if constexpr (std::is_same_v<O, WellModelMatrixAdapter<M, V, V>>) {
222 []( const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
223 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
224 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
226 using Scalar = typename V::field_type;
227 using LevelTransferPolicy
229 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
230 op, prm, weightsCalculator, pressureIndex);
236 []( const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
237 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
238 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
240 using Scalar = typename V::field_type;
242 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
243 op, prm, weightsCalculator, pressureIndex);
247 []( const O& op, const P& prm, const std::function<V()>& weightsCalculator, std::size_t pressureIndex) {
248 if (pressureIndex == std::numeric_limits<std::size_t>::max()) {
249 OPM_THROW(std::logic_error, "Pressure index out of bounds. It needs to specified for CPR");
251 using Scalar = typename V::field_type;
253 return std::make_shared<OwningTwoLevelPreconditioner<O, V, LevelTransferPolicy>>(
254 op, prm, weightsCalculator, pressureIndex);
263 F::addCreator( "gpuilu0", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
264 const double w = prm.get< double>( "relaxation", 1.0);
265 using field_type = typename V::field_type;
266 using GpuILU0 = typename gpuistl::
268 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, GpuILU0>>(
269 std::make_shared<GpuILU0>(op.getmat(), w));
272 F::addCreator( "gpuilu0float", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
273 const double w = prm.get< double>( "relaxation", 1.0);
274 using block_type = typename V::block_type;
275 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
276 using matrix_type_to =
277 typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
278 using GpuILU0 = typename gpuistl::
282 auto converted = std::make_shared<Converter>(op.getmat());
283 auto adapted = std::make_shared<Adapter>(std::make_shared<GpuILU0>(converted->getConvertedMatrix(), w));
284 converted->setUnderlyingPreconditioner(adapted);
288 F::addCreator( "gpujac", []( const O& op, const P& prm, const std::function<V()>&, std::size_t) {
289 const double w = prm.get< double>( "relaxation", 1.0);
290 using field_type = typename V::field_type;
297 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
298 std::make_shared<MatrixOwner>(op.getmat(), w));
301 F::addCreator( "opmgpuilu0", []( const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
302 const bool split_matrix = prm.get< bool>( "split_matrix", true);
303 const bool tune_gpu_kernels = prm.get< bool>( "tune_gpu_kernels", true);
304 const int mixed_precision_scheme = prm.get< int>( "mixed_precision_scheme", 0);
306 using field_type = typename V::field_type;
310 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
313 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme));
317 F::addCreator( "gpudilu", []( const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
318 const bool split_matrix = prm.get< bool>( "split_matrix", true);
319 const bool tune_gpu_kernels = prm.get< bool>( "tune_gpu_kernels", true);
320 const int mixed_precision_scheme = prm.get< int>( "mixed_precision_scheme", 0);
321 const bool reorder = prm.get< bool>( "reorder", true);
322 using field_type = typename V::field_type;
326 return std::make_shared<gpuistl::PreconditionerAdapter<V, V, MatrixOwner>>(
329 std::make_shared<MatrixOwner>(op.getmat(), op.getmat(), split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
332 F::addCreator( "gpudilufloat", []( const O& op, [[maybe_unused]] const P& prm, const std::function<V()>&, std::size_t) {
333 const bool split_matrix = prm.get< bool>( "split_matrix", true);
334 const bool tune_gpu_kernels = prm.get< bool>( "tune_gpu_kernels", true);
335 const int mixed_precision_scheme = prm.get< int>( "mixed_precision_scheme", 0);
336 const bool reorder = prm.get< bool>( "reorder", true);
338 using block_type = typename V::block_type;
339 using VTo = Dune::BlockVector<Dune::FieldVector<float, block_type::dimension>>;
340 using matrix_type_to = typename Dune::BCRSMatrix<Dune::FieldMatrix<float, block_type::dimension, block_type::dimension>>;
348 auto converted = std::make_shared<Converter>(op.getmat());
351 auto adapted = std::make_shared<Adapter>(std::make_shared<MatrixOwner>(
352 converted->getConvertedMatrix(), converted->getConvertedMatrix(),
353 split_matrix, tune_gpu_kernels, mixed_precision_scheme, reorder));
354 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:34
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:40
virtual bool hasPerfectUpdate() const override Definition: StandardPreconditioners_serial.hpp:44
virtual void pre(X &x, Y &y) override Definition: StandardPreconditioners_serial.hpp:45
virtual void post(X &x) override Definition: StandardPreconditioners_serial.hpp:46
virtual void update() override Definition: StandardPreconditioners_serial.hpp:43
virtual void apply(X &x, const Y &y) override Definition: StandardPreconditioners_serial.hpp:47
TrivialPreconditioner() Definition: StandardPreconditioners_serial.hpp:42
virtual Dune::SolverCategory::Category category() const override Definition: StandardPreconditioners_serial.hpp:48
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:45
@ ILU Do not perform modified ILU.
Definition: PreconditionerFactory.hpp:43
static Criterion criterion(const PropertyTree &prm) Definition: StandardPreconditioners_mpi.hpp:90
Definition: StandardPreconditioners_mpi.hpp:135
|