DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
SerialAdjReordering.hpp
Go to the documentation of this file.
1#pragma once
2
3#include "Geometric.hpp"
4
5#include "CorrectRCM.hpp"
6
7#include "Metis.hpp"
8
9#include <boost/graph/adjacency_list.hpp>
10#include <boost/graph/graph_utility.hpp>
11#include <boost/graph/minimum_degree_ordering.hpp>
12#include <boost/graph/cuthill_mckee_ordering.hpp>
13#include <boost/graph/properties.hpp>
14#include <boost/graph/bandwidth.hpp>
15
16namespace DNDS::Geom
17{
18
19 /**
20 * @brief Partition a (sub-)graph using Metis.
21 *
22 * Adjacency entries in [mat_begin, mat_end) are absolute indices offset by ind_offset.
23 * This function internally subtracts ind_offset and filters out cross-sub-graph
24 * references (entries outside [ind_offset, ind_offset + n_elem)) before building
25 * the 0-based CSR that Metis requires.
26 *
27 * @param ind_offset Absolute index of the first row in this sub-graph. For a full
28 * graph this is 0 and no filtering occurs.
29 */
31 tLocalMatStruct::const_iterator mat_begin, tLocalMatStruct::const_iterator mat_end, int nPart,
32 index ind_offset = 0,
33 std::string metisType = "KWAY", int metisNcuts = 3, int metisUfactor = 5, int metisSeed = 0)
34 {
35 idx_t nCell = _METIS::indexToIdx(size_t_to_signed<index>(mat_end - mat_begin));
36 idx_t nCon{1}, options[METIS_NOPTIONS];
37 METIS_SetDefaultOptions(options);
38 {
39 options[METIS_OPTION_OBJTYPE] = METIS_OBJTYPE_CUT;
40 options[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM; //? could try shem?
41 options[METIS_OPTION_IPTYPE] = METIS_IPTYPE_GROW;
42 options[METIS_OPTION_RTYPE] = METIS_RTYPE_FM;
43 // options[METIS_OPTION_NO2HOP] = 0; // only available in metis 5.1.0
44 options[METIS_OPTION_NCUTS] = std::max(metisNcuts, 1);
45 options[METIS_OPTION_NITER] = 10;
46 // options[METIS_OPTION_UFACTOR] = 30; // load imbalance factor, fow k-way
47 options[METIS_OPTION_UFACTOR] = metisUfactor;
48 options[METIS_OPTION_MINCONN] = 1;
49 options[METIS_OPTION_CONTIG] = 1; // ! forcing contigious partition now ? necessary?
50 options[METIS_OPTION_SEED] = metisSeed; // ! seeding 0 for determined result
51 options[METIS_OPTION_NUMBERING] = 0;
52 // options[METIS_OPTION_DBGLVL] = METIS_DBG_TIME | METIS_DBG_IPART;
53 // options[METIS_OPTION_DBGLVL] = METIS_DBG_TIME;
54 }
55 auto &cell2cellFaceV = mat_begin;
56
57 //! Build CSR with 0-based local indices: subtract ind_offset and skip
58 //! entries outside [ind_offset, ind_offset + nCell) (cross-sub-graph edges).
59 std::vector<idx_t> adjncy, xadj;
60 xadj.resize(nCell + 1);
61 xadj[0] = 0;
62 for (idx_t iC = 0; iC < nCell; iC++)
63 {
64 idx_t count = 0;
65 for (auto iCOther : cell2cellFaceV[iC])
66 {
67 index iLocal = iCOther - ind_offset;
68 if (iLocal >= 0 && iLocal < nCell)
69 count++;
70 }
71 xadj[iC + 1] = xadj[iC] + count;
72 }
73 adjncy.resize(xadj.back());
74 for (idx_t iC = 0; iC < nCell; iC++)
75 {
76 idx_t pos = xadj[iC];
77 for (auto iCOther : cell2cellFaceV[iC])
78 {
79 index iLocal = iCOther - ind_offset;
80 if (iLocal >= 0 && iLocal < nCell)
81 adjncy[pos++] = _METIS::indexToIdx(iLocal);
82 }
83 }
84
85 idx_t objval;
86 std::vector<idx_t> partOut(nCell);
87
88 int ret{0};
89 if (metisType == "RB")
90 ret = METIS_PartGraphRecursive(
91 &nCell, &nCon, xadj.data(), adjncy.data(), NULL, NULL, NULL,
92 &nPart, NULL, NULL, options, &objval, partOut.data());
93 else if (metisType == "KWAY")
94 ret = METIS_PartGraphKway(
95 &nCell, &nCon, xadj.data(), adjncy.data(), NULL, NULL, NULL,
96 &nPart, NULL, NULL, options, &objval, partOut.data());
97
98 DNDS_assert_info(ret == METIS_OK, fmt::format("Metis return not ok, [{}]", ret));
99
100 return partOut;
101 }
102
103 inline std::pair<std::vector<index>, std::vector<index>> ReorderSerialAdj_Metis(const tLocalMatStruct &mat)
104 {
105 idx_t nCell = _METIS::indexToIdx(size_t_to_signed<index>(mat.size()));
106 idx_t nCon{1}, options[METIS_NOPTIONS];
107 METIS_SetDefaultOptions(options);
108 {
109 options[METIS_OPTION_CTYPE] = METIS_CTYPE_SHEM;
110 options[METIS_OPTION_RTYPE] = METIS_RTYPE_FM;
111 options[METIS_OPTION_IPTYPE] = METIS_IPTYPE_EDGE;
112 options[METIS_OPTION_RTYPE] = METIS_RTYPE_SEP1SIDED;
113 options[METIS_OPTION_NSEPS] = 1;
114 options[METIS_OPTION_NITER] = 10;
115 options[METIS_OPTION_UFACTOR] = 30;
116 options[METIS_OPTION_COMPRESS] = 0; // do not compress
117 // options[METIS_OPTION_CCORDER] = 0; //use default?
118 options[METIS_OPTION_SEED] = 0; // ! seeding 0 for determined result
119 options[METIS_OPTION_PFACTOR] = 0; // not removing large vertices
120 options[METIS_OPTION_NUMBERING] = 0;
121 // options[METIS_OPTION_DBGLVL] = METIS_DBG_TIME | METIS_DBG_IPART;
122 }
123 const std::vector<std::vector<index>> &cell2cellFaceV = mat;
124 std::vector<idx_t> adjncy, xadj, perm, iPerm;
125 xadj.resize(nCell + 1);
126 xadj[0] = 0;
127 for (idx_t iC = 0; iC < nCell; iC++)
128 xadj[iC + 1] = signedIntSafeAdd<idx_t>(xadj[iC], size_t_to_signed<idx_t>(cell2cellFaceV[iC].size())); //! check overflow!
129 adjncy.resize(xadj.back());
130 for (idx_t iC = 0; iC < nCell; iC++)
131 std::copy(cell2cellFaceV[iC].begin(), cell2cellFaceV[iC].end(), adjncy.begin() + xadj[iC]);
132 perm.resize(nCell);
133 iPerm.resize(nCell);
134
135 int ret = METIS_NodeND(&nCell, xadj.data(), adjncy.data(), NULL, options, perm.data(), iPerm.data());
136 DNDS_assert_info(ret == METIS_OK, fmt::format("Metis return not ok, [{}]", ret));
137
138 std::vector<index> localFillOrderingNew2Old, localFillOrderingOld2New;
139
140 localFillOrderingNew2Old.resize(nCell);
141 localFillOrderingOld2New.resize(nCell);
142 for (index i = 0; i < nCell; i++)
143 {
144 localFillOrderingNew2Old[i] = perm[i];
145 localFillOrderingOld2New[i] = iPerm[i];
146 }
147
148 return {localFillOrderingNew2Old, localFillOrderingOld2New};
149 }
150
151 inline std::pair<std::vector<index>, std::vector<index>> ReorderSerialAdj_BoostMMD(const tLocalMatStruct &mat)
152 {
153 std::vector<index> localFillOrderingNew2Old, localFillOrderingOld2New;
154 using namespace boost;
155 using Graph = adjacency_list<vecS, vecS, directedS>;
156 Graph cell2cellG(mat.size());
157 const std::vector<std::vector<index>> &cell2cellFaceV = mat;
158 for (index iCell = 0; iCell < size_t_to_signed<index>(mat.size()); iCell++)
159 for (auto iCOther : cell2cellFaceV[iCell])
160 add_edge(iCell, iCOther, cell2cellG);
161 std::vector<index> supernodeSizes(mat.size(), 1), degree(mat.size(), 0);
162 localFillOrderingNew2Old.resize(mat.size(), 0);
163 localFillOrderingOld2New.resize(mat.size(), 0);
164 boost::property_map<Graph, vertex_index_t>::type id = get(vertex_index, cell2cellG);
165 minimum_degree_ordering(
166 cell2cellG,
167 make_iterator_property_map(degree.data(), id, degree[0]),
168 localFillOrderingOld2New.data(),
169 localFillOrderingNew2Old.data(),
170 make_iterator_property_map(supernodeSizes.data(), id, supernodeSizes[0]),
171 0,
172 id);
173 return {localFillOrderingNew2Old, localFillOrderingOld2New};
174 }
175
176 inline std::pair<std::vector<index>, std::vector<index>> ReorderSerialAdj_BoostRCM(const tLocalMatStruct &mat, index &bandWidthOld, index &bandWidthNew)
177 {
178 std::vector<index> localFillOrderingNew2Old, localFillOrderingOld2New;
179 using namespace boost;
180 using Graph = adjacency_list<vecS, vecS, undirectedS, property<vertex_color_t, default_color_type, property<vertex_degree_t, int>>>;
181 using Vertex = graph_traits<Graph>::vertex_descriptor;
182 // typedef graph_traits<Graph>::vertices_size_type size_type;
183 Graph cell2cellG(mat.size());
184 const std::vector<std::vector<index>> &cell2cellFaceV = mat;
185 bandWidthOld = 0;
186 for (index iCell = 0; iCell < size_t_to_signed<index>(mat.size()); iCell++)
187 for (auto iCOther : cell2cellFaceV[iCell])
188 add_edge(iCell, iCOther, cell2cellG), bandWidthOld = std::max(bandWidthOld, std::abs(iCell - iCOther));
189 localFillOrderingNew2Old.resize(mat.size(), 0);
190 localFillOrderingOld2New.resize(mat.size(), 0);
191 Vertex startVert = vertex(0, cell2cellG);
192 cuthill_mckee_ordering(cell2cellG, startVert, localFillOrderingNew2Old.rbegin(),
193 get(vertex_color, cell2cellG), get(vertex_degree, cell2cellG));
194 std::unordered_set<index> __checkOrder;
195 for (auto v : localFillOrderingNew2Old)
196 DNDS_assert(v < mat.size() && v >= 0), __checkOrder.insert(v);
197 DNDS_assert_info(__checkOrder.size() == localFillOrderingNew2Old.size(), "The output of boost::cuthill_mckee_ordering is invalid!");
198
199 for (index iCell = 0; iCell < size_t_to_signed<index>(mat.size()); iCell++)
200 localFillOrderingOld2New[localFillOrderingNew2Old[iCell]] = iCell;
201 for (auto v : localFillOrderingOld2New)
202 DNDS_assert(v < size_t_to_signed<index>(mat.size()) && v >= 0);
203 bandWidthNew = 0;
204 for (index iCell = 0; iCell < size_t_to_signed<index>(mat.size()); iCell++)
205 for (auto iCOther : cell2cellFaceV[iCell])
206 bandWidthNew = std::max(bandWidthNew, std::abs(localFillOrderingOld2New[iCell] - localFillOrderingOld2New[iCOther]));
207
208 return {localFillOrderingNew2Old, localFillOrderingOld2New};
209 }
210
211 template <class T>
213 {
214 const T *ptr_;
215 T offset_;
216
217 public:
218 using iterator_category = std::input_iterator_tag;
219 using value_type = T;
220 using difference_type = std::ptrdiff_t;
221 using pointer = const T *;
222 using reference = T; // NOT a real reference! It's a temporary.
223
224 OffsetIterator(const T *p, T off) : ptr_(p), offset_(off) {}
225
226 // Dereference returns a transformed VALUE
227 T operator*() const
228 {
229 return *ptr_ + offset_;
230 }
231
232 // Iterator movement
234 {
235 ++ptr_;
236 return *this;
237 }
238
240 {
241 OffsetIterator tmp = *this;
242 ++(*this);
243 return tmp;
244 }
245
246 // Comparison
247 bool operator==(const OffsetIterator &other) const
248 {
249 return ptr_ == other.ptr_;
250 }
251
252 bool operator!=(const OffsetIterator &other) const
253 {
254 return !(*this == other);
255 }
256 };
257 // Proxy wrapper class
258 template <class T>
260 {
261 const T *begin_;
262 const T *end_;
263 T offset_;
264
265 public:
266 // Constructor from a sub-range of a vector
267 OffsetRange(typename std::vector<T>::const_iterator b, typename std::vector<T>::const_iterator e, T offset)
268 : begin_(std::addressof(*b)), end_(std::addressof(*e)), offset_(offset) {}
269
270 // Support for full vector
271 explicit OffsetRange(const std::vector<const T> &vec, T offset)
272 : begin_(vec.data()), end_(vec.data() + vec.size()), offset_(offset) {}
273
274 OffsetIterator<T> begin() const { return {begin_, offset_}; }
275 OffsetIterator<T> end() const { return {end_, offset_}; }
276 ptrdiff_t size() const { return end_ - begin_; }
277 };
278
279 inline std::pair<std::vector<index>, std::vector<index>> ReorderSerialAdj_CorrectRCM(
280 tLocalMatStruct::const_iterator mat_begin, tLocalMatStruct::const_iterator mat_end,
281 index &bandWidthOld, index &bandWidthNew,
282 index offset = 0)
283 {
284 std::vector<index> localFillOrderingNew2Old, localFillOrderingOld2New;
285 bandWidthOld = 0;
286
287 // for (index iCell = 0; iCell < this->NumCell(); ++iCell)
288 // {
289 // cell2cellFaceVLocal[iCell].resize(4);
290 // int x = iCell % 20, y = iCell / 20;
291 // cell2cellFaceVLocal[iCell][0] = mod(x - 1, 20) + y * 20;
292 // cell2cellFaceVLocal[iCell][1] = mod(x + 1, 20) + y * 20;
293 // cell2cellFaceVLocal[iCell][2] = x + mod(y - 1, 20) * 20;
294 // cell2cellFaceVLocal[iCell][3] = x + mod(y + 1, 20) * 20;
295 // }
296
297 index mat_size = mat_end - mat_begin;
298 // std::cout << "HHH " << offset << "," << mat_size << std::endl;
299 for (index iCell = 0; iCell < size_t_to_signed<index>(mat_size); iCell++)
300 for (auto iCOther : mat_begin[iCell])
301 bandWidthOld = std::max(bandWidthOld, std::abs(iCell + offset - iCOther));
302
303 localFillOrderingNew2Old.resize(mat_size, 0);
304 localFillOrderingOld2New.resize(mat_size, 0);
305 auto graphFunctor = [&](index i)
306 {
307 DNDS_assert(i >= 0 && i < mat_size);
308 return OffsetRange<index>(mat_begin[i].cbegin(), mat_begin[i].cend(), -offset); // mind that OffsetRange adds the offset
309 }; // todo: need improvement in CorrectRCM: can pass a temporary functor and store
310 auto graph = CorrectRCM::UndirectedGraphProxy(graphFunctor, size_t_to_signed<int64_t>(mat_size));
311 int ret = graph.CheckAdj();
313 graph,
314 [&](index i) -> index &
315 {
316 return localFillOrderingOld2New.at(i);
317 },
318 0);
319 for (auto &v : localFillOrderingOld2New)
320 v = localFillOrderingOld2New.size() - 1 - v;
321
322 std::unordered_set<index>
323 __checkOrder;
324 for (auto v : localFillOrderingOld2New)
325 DNDS_assert(v < size_t_to_signed<index>(mat_size) && v >= 0), __checkOrder.insert(v);
326 DNDS_check_throw_info(__checkOrder.size() == localFillOrderingOld2New.size(), "The output of CorrectRCM::CuthillMcKeeOrdering is invalid!");
327
328 for (index iCell = 0; iCell < size_t_to_signed<index>(mat_size); iCell++)
329 localFillOrderingNew2Old[localFillOrderingOld2New[iCell]] = iCell;
330 for (auto v : localFillOrderingNew2Old)
331 DNDS_assert(v < size_t_to_signed<index>(mat_size) && v >= 0);
332 bandWidthNew = 0;
333 for (index iCell = 0; iCell < size_t_to_signed<index>(mat_size); iCell++)
334 for (auto iCOther : mat_begin[iCell])
335 bandWidthNew = std::max(bandWidthNew, std::abs(localFillOrderingOld2New[iCell] - localFillOrderingOld2New.at(iCOther - offset)));
336
337 return {localFillOrderingNew2Old, localFillOrderingOld2New};
338 }
339
340 /**
341 * @brief Partition a sub-graph and optionally RCM-reorder within each partition.
342 *
343 * When used for second-level (inner) partitioning, pass full_mat_begin / full_n_elem
344 * so that cross-sub-graph references in the full graph are updated after permutation.
345 * See implementation in SerialAdjReordering.cpp for detailed documentation.
346 */
347 std::vector<index> ReorderSerialAdj_PartitionMetisC(
348 tLocalMatStruct::iterator mat_begin, tLocalMatStruct::iterator mat_end,
349 std::vector<index>::iterator i_new2old_begin, std::vector<index>::iterator i_new2old_end,
350 int nParts,
351 index ind_offset,
352 bool do_rcm,
353 index &bwOldM, index &bwNewM,
354 tLocalMatStruct::iterator full_mat_begin = tLocalMatStruct::iterator{},
355 index full_n_elem = 0);
356
357}
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:113
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:108
#define DNDS_check_throw_info(expr, info)
Same as DNDS_check_throw but attaches a user-supplied info message to the thrown std::runtime_error.
Definition Errors.hpp:96
bool operator==(const OffsetIterator &other) const
bool operator!=(const OffsetIterator &other) const
std::input_iterator_tag iterator_category
OffsetRange(typename std::vector< T >::const_iterator b, typename std::vector< T >::const_iterator e, T offset)
OffsetIterator< T > begin() const
OffsetRange(const std::vector< const T > &vec, T offset)
OffsetIterator< T > end() const
void CuthillMcKeeOrdering(TGraph &&G, TFFiller &&FFiller, TIndex &&root, std::optional< std::reference_wrapper< std::ostream > > os={})
std::pair< std::vector< index >, std::vector< index > > ReorderSerialAdj_BoostRCM(const tLocalMatStruct &mat, index &bandWidthOld, index &bandWidthNew)
std::pair< std::vector< index >, std::vector< index > > ReorderSerialAdj_Metis(const tLocalMatStruct &mat)
std::pair< std::vector< index >, std::vector< index > > ReorderSerialAdj_BoostMMD(const tLocalMatStruct &mat)
std::pair< std::vector< index >, std::vector< index > > ReorderSerialAdj_CorrectRCM(tLocalMatStruct::const_iterator mat_begin, tLocalMatStruct::const_iterator mat_end, index &bandWidthOld, index &bandWidthNew, index offset=0)
std::vector< index > ReorderSerialAdj_PartitionMetisC(tLocalMatStruct::iterator mat_begin, tLocalMatStruct::iterator mat_end, std::vector< index >::iterator i_new2old_begin, std::vector< index >::iterator i_new2old_end, int nParts, index ind_offset, bool do_rcm, index &bwOldM, index &bwNewM, tLocalMatStruct::iterator full_mat_begin, index full_n_elem)
Partition a sub-graph into nParts and optionally apply RCM reordering within each part.
auto PartitionSerialAdj_Metis(tLocalMatStruct::const_iterator mat_begin, tLocalMatStruct::const_iterator mat_end, int nPart, index ind_offset=0, std::string metisType="KWAY", int metisNcuts=3, int metisUfactor=5, int metisSeed=0)
Partition a (sub-)graph using Metis.
std::vector< std::vector< index > > tLocalMatStruct
Definition Geometric.hpp:66
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:107
Eigen::Matrix< real, 5, 1 > v
tVec b(NCells)