1#ifndef DUNE_COPASI_SOLVER_ISTL_BLOCK_JACOBI_HH
2#define DUNE_COPASI_SOLVER_ISTL_BLOCK_JACOBI_HH
7#include <dune/pdelab/common/algebra.hh>
9#include <dune/istl/preconditioner.hh>
10#include <dune/istl/solver.hh>
16template<Concept::AssembledLinearOperator O,
class Alloc = std::allocator<
void>>
17class BlockJacobi final :
public Preconditioner<typename O::domain_type, typename O::range_type>
19 using M =
typename O::matrix_type;
20 using X =
typename O::domain_type;
21 using Y =
typename O::range_type;
25 const ParameterTree& config,
26 const Alloc& alloc = Alloc())
30 , _parallel{ _config.get(
"parallel", false) }
31 , _iterations{ _config.get(
"iterations", std::size_t{ 1 }) }
32 , _weight{ _config.get(
"relaxation", 1.) }
38 void pre(X& x, Y& b)
override
41 this->
pre(PDELab::Execution::par, x, b);
43 this->
pre(PDELab::Execution::seq, x, b);
46 void pre(
auto policy, X& , Y& )
48 const ParameterTree& sub_solver = _config.sub(
"diagonal_solver");
49 if constexpr (IsNumber<typename M::block_type>{}) {
50 if (sub_solver.getSubKeys().size() != 0)
51 spdlog::warn(
"'diagona_solver' contains sub-sections that will be ignored");
52 if (sub_solver.get(
"type",
"DenseInverse") !=
"DenseInverse")
53 spdlog::warn(
"Block jacobi at this level is scalar, meaning that the diagonal solver can only be 'DenseInverse'");
56 const M& mat = _op->getmat();
58 _diag_inv.resize(mat.M());
59 PDELab::forEach(policy, mat, [&](
const auto& row,
auto i) {
60 auto row_it = row.begin();
61 for (; row_it != row.end() and row_it.index() != i; ++row_it) {
63 if (row_it == row.end())
64 throw format_exception(MathError{},
"BlockJacobi can only be done on matrices with diagonal entries");
65 if constexpr (IsNumber<typename M::block_type>{}) {
66 _diag_inv[i] = 1. / (*row_it);
68 auto i_str = fmt::format(
"{}", i);
69 const ParameterTree& sub_solver_i =
70 (sub_solver.hasSub(i_str)) ? sub_solver.sub(i_str) : sub_solver;
72 using MatrixAdapter = Dune::
73 MatrixAdapter<typename M::block_type, typename X::value_type, typename Y::value_type>;
74 using MatrixAdapterAlloc =
75 typename std::allocator_traits<Alloc>::template rebind_alloc<MatrixAdapter>;
76 MatrixAdapterAlloc maalloc(_alloc);
77 auto diag_op = std::allocate_shared<MatrixAdapter>(maalloc, *row_it);
83 void apply(X& v,
const Y& rhs)
override
86 this->
apply(PDELab::Execution::par, v, rhs);
88 this->
apply(PDELab::Execution::seq, v, rhs);
91 void apply(
auto policy, X& v,
const Y& rhs)
95 std::optional<X> x_tmp;
96 const M& mat = _op->getmat();
98 using Xi =
typename X::value_type;
99 [[maybe_unused]] std::optional<Xi> vi;
100 X& x = (v.size() > 1) ? x_tmp.emplace(v) : v;
102 for (std::size_t it = 0; it != _iterations; ++it) {
103 if constexpr (IsNumber<typename M::block_type>{}) {
104 PDELab::forEach(policy, mat, [&b, &x, v, weight = _weight, diag_inv = _diag_inv.data()](
const auto& row,
auto i)
mutable {
106 for (auto col = row.begin(); col != row.end(); ++col)
107 b[i] -= (*col) * v[col.index()];
109 x[i] += weight * diag_inv[i] * b[i];
112 PDELab::forEach(policy, mat, [&b, &x, v, vi, weight = _weight, diag_inv = _diag_inv.data()](
const auto& row,
auto i)
mutable {
114 InverseOperatorResult res;
115 for (auto col = row.begin(); col != row.end(); ++col)
116 col->mmv(v[col.index()], b[i]);
118 auto& vir = vi ? *vi : vi.emplace();
121 diag_inv[i]->apply(vir, b[i], res);
122 PDELab::axpy(x[i], weight, vir);
132 SolverCategory::Category
category()
const override {
return SolverCategory::sequential; }
136 return [](
const std::shared_ptr<LinearOperator<X, Y>>& op,
137 const ParameterTree& config,
138 const Alloc& alloc) -> std::shared_ptr<Preconditioner<X, Y>> {
140 using PrecAlloc =
typename std::allocator_traits<Alloc>::template rebind_alloc<BlockJacobi>;
141 PrecAlloc palloc(alloc);
144 auto aop = std::dynamic_pointer_cast<O>(op);
146 throw format_exception(InvalidStateException{},
"Linear operator does not hold a matrix!");
148 return std::allocate_shared<BlockJacobi>(palloc, aop, config);
153 using DiagInverse = std::conditional_t<
154 IsNumber<typename M::block_type>::value,
155 typename M::block_type,
156 std::shared_ptr<InverseOperator<typename X::value_type, typename Y::value_type>>>;
158 using DiagInverseAlloc =
159 typename std::allocator_traits<Alloc>::template rebind_alloc<DiagInverse>;
162 std::shared_ptr<const O> _op;
163 std::vector<DiagInverse, DiagInverseAlloc> _diag_inv;
164 ParameterTree _config;
166 std::size_t _iterations;
Definition: block_jacobi.hh:18
void apply(auto policy, X &v, const Y &rhs)
Definition: block_jacobi.hh:91
void pre(X &x, Y &b) override
Definition: block_jacobi.hh:38
void apply(X &v, const Y &rhs) override
Definition: block_jacobi.hh:83
BlockJacobi(const std::shared_ptr< const O > &op, const ParameterTree &config, const Alloc &alloc=Alloc())
Definition: block_jacobi.hh:24
void pre(auto policy, X &, Y &)
Definition: block_jacobi.hh:46
~BlockJacobi() override=default
void post(X &) override
Definition: block_jacobi.hh:130
static auto factory()
Definition: block_jacobi.hh:134
SolverCategory::Category category() const override
Definition: block_jacobi.hh:132
Definition: block_jacobi.hh:14
std::shared_ptr< InverseOperator< typename Op::domain_type, typename Op::range_type > > makeInverseOperator(const std::shared_ptr< Op > &op, const ParameterTree &config, const Alloc &alloc=Alloc())
Definition: inverse.hh:65
auto format_exception(Exception &&e, fmt::format_string< Args... > format, Args &&... args)
Definition: exceptions.hh:23