Dune::Copasi
Loading...
Searching...
No Matches
umfpack.hh
Go to the documentation of this file.
1#ifndef DUNE_COPASI_SOLVER_ISTL_UMFPACK_HH
2#define DUNE_COPASI_SOLVER_ISTL_UMFPACK_HH
3
4#if !HAVE_SUITESPARSE_UMFPACK
5#error Trying to use umfpack wrapper without including umfpack
6#endif
7
9
10#include <dune/istl/solver.hh>
11
12#include <umfpack.h>
13
14#include <algorithm>
15#include <ranges>
16#include <version>
17
18#if __cpp_lib_execution >= 201603L
19#include <execution>
20#endif
21
22namespace Dune::Copasi::ISTL {
23
24template<Concept::AssembledLinearOperator O, class Alloc = std::allocator<void>>
25class UMFPackWapper final : public InverseOperator<typename O::domain_type, typename O::range_type>
26{
27 using Matrix = typename O::matrix_type;
28 using X = typename O::domain_type;
29 using Y = typename O::range_type;
30 using MI = PDELab::MultiIndex<std::size_t, blockLevel<Matrix>()>;
31
32 static const inline std::map<std::string, double> umf_strategy_map = {
33 { "AUTO", UMFPACK_STRATEGY_AUTO },
34 { "UNSYMMETRIC", UMFPACK_STRATEGY_UNSYMMETRIC },
35 { "UMFPACK_STRATEGY_SYMMETRIC", UMFPACK_STRATEGY_SYMMETRIC },
36 };
37
38 static const inline std::map<std::string, double> umf_ordering_map = {
39 { "CHOLMOD", UMFPACK_ORDERING_CHOLMOD },
40 { "AMD", UMFPACK_ORDERING_AMD },
41 { "METIS", UMFPACK_ORDERING_METIS },
42 { "BEST", UMFPACK_ORDERING_BEST },
43 { "NONE", UMFPACK_ORDERING_NONE },
44#if defined(UMFPACK_VER) && (UMFPACK >= UMFPACK_VER_CODE(6, 0))
45 { "METIS_GUARD", UMFPACK_ORDERING_METIS_GUARD },
46#endif
47 };
48
49 // matrix entries as a (val, row, col) triplet
50 struct Entry
51 {
52 double val;
53 MI row;
54 MI col;
55 };
56
57 // rebind allocators to correct types
58 using Int64Alloc = typename std::allocator_traits<Alloc>::template rebind_alloc<int64_t>;
59 using DoubleAlloc = typename std::allocator_traits<Alloc>::template rebind_alloc<double>;
60 using EntryAlloc = typename std::allocator_traits<Alloc>::template rebind_alloc<Entry>;
61
62 struct SymbolicDeleter
63 {
64 void operator()(void* ptr) const { umfpack_dl_free_symbolic(&ptr); }
65 };
66 struct NumericDeleter
67 {
68 void operator()(void* ptr) const { umfpack_dl_free_numeric(&ptr); }
69 };
70
71public:
72 // note allocator only allocates data for what UMFPack allows to be externally owned. Otherwise,
73 // umprfack still allocates some memory internally
74 UMFPackWapper(const std::shared_ptr<const O>& op,
75 const ParameterTree& config,
76 const Alloc& alloc = Alloc())
77 : _alloc(alloc)
78 , _umf_offsets(Int64Alloc(_alloc))
79 , _umf_rows(Int64Alloc(_alloc))
80 , _umf_iworkspace(Int64Alloc(_alloc))
81 , _umf_dworkspace(DoubleAlloc(_alloc))
82 , _umf_values(DoubleAlloc(_alloc))
83 , _umf_x(DoubleAlloc(_alloc))
84 , _umf_b(DoubleAlloc(_alloc))
85 , _report_symbolic{config.get("report_control", false)}
86 , _report_numeric{config.get("report_numeric", false)}
87 , _report_info{config.get("report_info", true)}
88 {
89 // set up umfpack configuration options
91 auto set_config = [&](auto id, auto key) {
92 _umf_control[id] = config.get(key, _umf_control[id]);
93 };
94 set_config(UMFPACK_PRL, "verbosity");
95 set_config(UMFPACK_DENSE_ROW, "dense_row");
96 set_config(UMFPACK_DENSE_COL, "dense_col");
97 set_config(UMFPACK_PIVOT_TOLERANCE, "pivot_tolerance");
98 set_config(UMFPACK_BLOCK_SIZE, "block_size");
99 set_config(UMFPACK_DENSE_ROW, "dense_row");
100 set_config(UMFPACK_ALLOC_INIT, "alloc_init");
101 set_config(UMFPACK_IRSTEP, "irstep");
102 set_config(UMFPACK_FIXQ, "fixq");
103 set_config(UMFPACK_AMD_DENSE, "amd_dense");
104 set_config(UMFPACK_SYM_PIVOT_TOLERANCE, "sym_pivot_tolerance");
105 set_config(UMFPACK_SCALE, "scale");
106 set_config(UMFPACK_FRONT_ALLOC_INIT, "front_alloc_init");
107 set_config(UMFPACK_DROPTOL, "droptol");
108 set_config(UMFPACK_AGGRESSIVE, "aggressive");
109 set_config(UMFPACK_SINGLETONS, "singletons");
110 _umf_control[UMFPACK_STRATEGY] = umf_strategy_map.at(config.get("strategy", "AUTO"));
111 _umf_control[UMFPACK_ORDERING] = umf_ordering_map.at(config.get("ordering", "CHOLMOD"));
112 // report about umfpack settings
113 if (config.get("report_control", false))
115 // store matrix into a format that unfpack understand
116 storeFlatMatrix(op->getmat());
117 }
118
119private:
120 template<class M>
121 void forEachMatrixEntry(M&& mat, const auto& call_back, MI row = MI(), MI col = MI())
122 {
123 using MD = std::decay_t<M>;
124 if constexpr (IsNumber<MD>::value) {
125 std::invoke(call_back, mat, row, col);
126 } else { // assume row mayor access!
127 for (auto row_it = mat.begin(); row_it != mat.end(); ++row_it) {
128 auto this_row = push_back(row, row_it.index());
129 for (auto col_it = row_it->begin(); col_it != row_it->end(); ++col_it) {
130 forEachMatrixEntry(*col_it, call_back, this_row, push_back(col, col_it.index()));
131 }
132 }
133 }
134 }
135
136 void storeFlatMatrix(const Matrix& mat)
137 {
138 EntryAlloc entryalloc(_alloc);
139 std::vector<Entry, EntryAlloc> entries(entryalloc);
140 forEachMatrixEntry(mat, [&](double val, MI row, MI col) {
141 entries.emplace_back(Entry{ val, row, col });
142 });
143
144 if (blockLevel<Matrix>() > 1) {
145
146 auto row_compare = [](auto lhs, auto rhs) {
147 return std::ranges::lexicographical_compare(lhs.row, rhs.row);
148 };
149 // transform row multi-index into a flat index (only needed if matrix is nested)
150 // note that sort does not need to order the col entry
151 std::sort(
152#if __cpp_lib_execution >= 201603L
153 std::execution::par,
154#endif
155 entries.begin(),
156 entries.end(),
158
159 MI last_row;
160 std::size_t flat_row = -1;
161 for (auto& [val, row, col] : entries) {
162 // row is new [[pre: rows are sorted]]: increase flat row count
163 if (row != std::exchange(last_row, row))
164 ++flat_row;
165 // store flat row in first index of row multi-index
166 row.resize(1);
167 row[0] = flat_row;
168 }
169 // number of flat rows
170 _umf_n = flat_row + 1;
171 } else {
172 _umf_n = mat.N();
173 }
174
175 auto col_row_compare = [](auto lhs, auto rhs) {
176 return std::ranges::lexicographical_compare(std::array{ lhs.col, lhs.row },
177 std::array{ rhs.col, rhs.row },
178 std::ranges::lexicographical_compare);
179 };
180
181 // order indices by column (note that sort does need to be order on row entry too)
182 std::sort(
183#if __cpp_lib_execution >= 201603L
184 std::execution::par,
185#endif
186 entries.begin(),
187 entries.end(),
189
190 // initialize umfpack data
191 _umf_values.clear();
192 _umf_offsets.clear();
193 _umf_rows.clear();
194 _umf_rows.reserve(entries.size());
195 _umf_values.reserve(entries.size());
196
197 // transform col multi-index into a flat offset on the row flat indices
198 MI last_col;
199 for (std::size_t i = 0; i != entries.size(); ++i) {
200 auto [val, row, col] = entries[i];
201 // col is new [[pre: cols are sorted]]: store new column offset
202 if (col != std::exchange(last_col, col))
203 _umf_offsets.push_back(i);
204 // store value and flat row index
205 _umf_values.push_back(val);
206 // [[pre: rows are sorted within column]]
207 _umf_rows.push_back(static_cast<int64_t>(row[0]));
208 }
209 _umf_m = _umf_offsets.size();
210 _umf_offsets.push_back(entries.size());
211 }
212
213 void factorize()
214 {
215 assert(not _umf_offsets.empty());
216 assert(not _umf_rows.empty());
217 assert(not _umf_values.empty());
218
219 // pre-order matrix columns
220 void* symbolic_ptr = nullptr;
222 _umf_m,
223 _umf_offsets.data(),
224 _umf_rows.data(),
225 _umf_values.data(),
228 _umf_info);
229 std::unique_ptr<void, SymbolicDeleter> _umf_symbolic_ptr{ std::exchange(symbolic_ptr, nullptr) };
232
233 // factorize matrix
234 void* numeric_ptr = nullptr;
236 _umf_rows.data(),
237 _umf_values.data(),
238 _umf_symbolic_ptr.get(),
241 _umf_info);
242 _umf_numeric_ptr = std::unique_ptr<void, NumericDeleter>{ std::exchange(numeric_ptr, nullptr) };
243 // report results
244 if (_report_numeric)
246 }
247
248public:
249 void apply(X& x, Y& b, InverseOperatorResult& res) override
250 {
251 // make sure matrix and buffers are setup
253 factorize();
254 _umf_iworkspace.resize(_umf_n);
255 _umf_dworkspace.resize(_umf_n * (_umf_control[UMFPACK_IRSTEP] > 0. ? 5 : 1));
256
257 // assign b to working array
258 _umf_b.resize(_umf_n);
259 std::size_t flat_i = 0;
260 PDELab::forEachContainerEntry(b, [&](auto val) { _umf_b[flat_i++] = val; });
261
262 // reserve data for result
263 _umf_x.resize(_umf_m);
264
265 // actual solve
267 _umf_offsets.data(),
268 _umf_rows.data(),
269 _umf_values.data(),
270 _umf_x.data(),
271 _umf_b.data(),
272 _umf_numeric_ptr.get(),
274 _umf_info,
275 _umf_iworkspace.data(),
276 _umf_dworkspace.data());
277 // report results
278 if (_report_info)
280
281 // assign x from working array
282 flat_i = 0;
283 PDELab::forEachContainerEntry(x, [&](auto& val) { val = _umf_x[flat_i++]; });
284
285 // store result info
286 res.iterations = 1;
287 res.converged = (status == 0);
288 }
289
290 void apply(X& x, Y& b, [[maybe_unused]] double reduction, InverseOperatorResult& res) override
291 {
292 this->apply(x, b, res);
293 }
294
295 SolverCategory::Category category() const override
296 {
297 return SolverCategory::Category::sequential;
298 }
299
300 static auto factory()
301 {
302 return [](const std::shared_ptr<LinearOperator<X, Y>>& op,
303 const ParameterTree& config,
304 const Alloc& alloc) -> std::shared_ptr<InverseOperator<X, Y>> {
305 // make an allocator for this preconditioner
306 using OpAlloc = typename std::allocator_traits<Alloc>::template rebind_alloc<UMFPackWapper>;
308
309 // cast operator to something the preconditioner expects
310 auto aop = std::dynamic_pointer_cast<O>(op);
311 if (not aop)
312 throw format_exception(InvalidStateException{}, "Linear operator does not hold a matrix!");
313 // construct preconditioner instance
314 return std::allocate_shared<UMFPackWapper>(opalloc, aop, config, alloc);
315 };
316 }
317
320 std::vector<int64_t, Int64Alloc> _umf_offsets;
321 std::vector<int64_t, Int64Alloc> _umf_rows;
322 std::vector<int64_t, Int64Alloc> _umf_iworkspace;
323 std::vector<double, DoubleAlloc> _umf_dworkspace;
324 std::vector<double, DoubleAlloc> _umf_values;
325 std::vector<double, DoubleAlloc> _umf_x;
326 std::vector<double, DoubleAlloc> _umf_b;
327 std::unique_ptr<void, NumericDeleter> _umf_numeric_ptr;
333};
334
335} // namespace Dune::Copasi::ISTL
336
337#endif // DUNE_COPASI_SOLVER_ISTL_UMFPACK_HH
Definition umfpack.hh:26
double _umf_info[UMFPACK_INFO]
Definition umfpack.hh:329
void apply(X &x, Y &b, InverseOperatorResult &res) override
Definition umfpack.hh:249
double _umf_control[UMFPACK_CONTROL]
Definition umfpack.hh:328
std::vector< int64_t, Int64Alloc > _umf_offsets
Definition umfpack.hh:320
std::unique_ptr< void, NumericDeleter > _umf_numeric_ptr
Definition umfpack.hh:327
int _umf_m
Definition umfpack.hh:318
std::vector< int64_t, Int64Alloc > _umf_iworkspace
Definition umfpack.hh:322
bool _report_numeric
Definition umfpack.hh:331
bool _report_info
Definition umfpack.hh:332
UMFPackWapper(const std::shared_ptr< const O > &op, const ParameterTree &config, const Alloc &alloc=Alloc())
Definition umfpack.hh:74
std::vector< double, DoubleAlloc > _umf_values
Definition umfpack.hh:324
std::vector< double, DoubleAlloc > _umf_b
Definition umfpack.hh:326
SolverCategory::Category category() const override
Definition umfpack.hh:295
static auto factory()
Definition umfpack.hh:300
Alloc _alloc
Definition umfpack.hh:319
int _umf_n
Definition umfpack.hh:318
void apply(X &x, Y &b, double reduction, InverseOperatorResult &res) override
Definition umfpack.hh:290
std::vector< double, DoubleAlloc > _umf_x
Definition umfpack.hh:325
std::vector< int64_t, Int64Alloc > _umf_rows
Definition umfpack.hh:321
bool _report_symbolic
Definition umfpack.hh:330
std::vector< double, DoubleAlloc > _umf_dworkspace
Definition umfpack.hh:323
Definition block_jacobi.hh:14
auto format_exception(Exception &&e, fmt::format_string< Args... > format, Args &&... args)
Definition exceptions.hh:23
constexpr bool is_bitflags_v
Alias for Bitflag indicator.
Definition bit_flags.hh:24