DNDSR 0.1.0.dev1+gcd065ad
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
SerialAdjReordering.cpp
Go to the documentation of this file.
2
3namespace DNDS::Geom
4{
5 /**
6 * @brief Partition a sub-graph into nParts and optionally apply RCM reordering within each part.
7 *
8 * This function operates on the sub-range [mat_begin, mat_end), which represents
9 * n_elem rows of an adjacency graph. Adjacency entries use absolute indices
10 * (offset by ind_offset), i.e., cell i's neighbors are stored as (ind_offset + local_index).
11 *
12 * When called for second-level (inner) partitioning, the sub-range's adjacency rows
13 * may contain cross-sub-graph references pointing outside [ind_offset, ind_offset + n_elem).
14 * These arise because the first-level partitioning reindexed the full graph, and boundary
15 * cells retain edges to cells in other first-level partitions. Two problems follow:
16 *
17 * 1. Metis and the partition-adjacency builder require a self-contained local graph;
18 * cross-sub-graph edges must be excluded from these operations.
19 *
20 * 2. When cells within this sub-range are permuted, all references to them -- including
21 * those in other sub-ranges' rows -- must be updated to maintain the bidirectional
22 * graph invariant. The optional full_mat_begin / full_n_elem parameters provide
23 * access to the full graph for this cross-sub-graph fixup.
24 *
25 * The RCM section trims cross-sub-partition edges before reordering (each sub-partition
26 * only needs internal connectivity for bandwidth reduction). The RCM reindex uses
27 * (ind_offset + local_offset) as the base for correct absolute index arithmetic.
28 *
29 * @param full_mat_begin Iterator to the start of the full graph (for cross-sub-graph
30 * reference updates). Pass default (empty) when this sub-range
31 * IS the full graph (first-level call).
32 * @param full_n_elem Total number of rows in the full graph. 0 means no cross-sub-graph
33 * fixup is needed.
34 */
36 tLocalMatStruct::iterator mat_begin, tLocalMatStruct::iterator mat_end,
37 std::vector<index>::iterator i_new2old_begin, std::vector<index>::iterator i_new2old_end,
38 int nParts,
39 index ind_offset,
40 bool do_rcm,
41 index &bwOldM, index &bwNewM,
42 tLocalMatStruct::iterator full_mat_begin,
43 index full_n_elem)
44 {
45 index n_elem = mat_end - mat_begin;
46 DNDS_check_throw(i_new2old_end - i_new2old_begin == n_elem);
47
48 nParts = std::max(1, nParts);
49 std::vector<index> part_start(nParts + 1, 0);
50
51 if (nParts > 1)
52 {
53 //! PartitionSerialAdj_Metis now handles the filtering and rebasing internally
54 //! via its ind_offset parameter.
55 auto partition_local = PartitionSerialAdj_Metis(mat_begin, mat_end, nParts, ind_offset);
56
57 //! Build partition-level adjacency (which partitions neighbor which).
58 //! Only intra-sub-graph edges matter for partition adjacency.
59 std::vector<std::set<int>> localPartsAdjV(nParts);
60 for (index iC = 0; iC < n_elem; iC++)
61 {
62 auto iP = partition_local.at(iC);
63 for (auto iCOther : mat_begin[iC])
64 {
65 index iLocal = iCOther - ind_offset;
66 if (iLocal >= 0 && iLocal < n_elem)
67 {
68 auto iPOther = partition_local.at(iLocal);
69 if (iPOther != iP)
70 localPartsAdjV.at(iP).insert(iPOther);
71 }
72 }
73 }
74 tLocalMatStruct localPartsAdj;
75 localPartsAdj.reserve(nParts);
76 for (auto &row : localPartsAdjV)
77 {
78 localPartsAdj.emplace_back();
79 localPartsAdj.back().reserve(row.size());
80 for (auto &v : row)
81 localPartsAdj.back().emplace_back(v);
82 }
83
85 auto [partNew2Old_, partOld2New_] = ReorderSerialAdj_CorrectRCM(
86 localPartsAdj.begin(), localPartsAdj.end(), bwOld, bwNew, 0);
87 for (auto &iP : partition_local)
88 {
89 DNDS_assert(partOld2New_[iP] < nParts);
90 iP = partOld2New_[iP];
91 }
92
93 std::vector<index> n_elem_part(nParts, 0);
94 for (index iC = 0; iC < n_elem; iC++)
95 {
96 auto iP = partition_local.at(iC);
97 n_elem_part[iP]++;
98 }
99 for (int i = 0; i < nParts; i++)
100 part_start[i + 1] = part_start[i] + n_elem_part[i];
101
102 //! Permute rows according to partition assignment; build i_old2new mapping.
103 tLocalMatStruct new_mat(n_elem);
104 std::vector<index> new_i_new2old(n_elem);
105 std::vector<index> i_old2new(n_elem);
106 n_elem_part.assign(n_elem, 0);
107 for (index iC = 0; iC < n_elem; iC++)
108 {
109 auto iP = partition_local[iC];
110 index iNew = part_start[iP] + n_elem_part[iP];
111 DNDS_assert(iNew < n_elem);
112 new_i_new2old[iNew] = i_new2old_begin[iC];
113 new_mat[iNew] = std::move(mat_begin[iC]);
114 i_old2new[iC] = iNew;
115 n_elem_part[iP]++;
116 }
117
118 //! Remap in-range references within the permuted sub-range rows.
119 //! Cross-sub-graph references (outside [ind_offset, ind_offset + n_elem)) are
120 //! left unchanged here; they point to cells not touched by this permutation.
121 for (auto &row : new_mat)
122 for (auto &iCOther : row)
123 {
124 index iLocal = iCOther - ind_offset;
125 if (iLocal >= 0 && iLocal < n_elem)
126 iCOther = ind_offset + i_old2new[iLocal];
127 }
128
129 for (index iC = 0; iC < n_elem; iC++)
130 {
131 i_new2old_begin[iC] = new_i_new2old[iC]; // composing the oldest mapping
132 mat_begin[iC] = std::move(new_mat[iC]);
133 }
134
135 //! Update cross-sub-graph references: rows OUTSIDE [mat_begin, mat_end) may
136 //! contain edges pointing into this sub-range with pre-permutation indices.
137 //! Walk the full graph and apply i_old2new to any such references.
138 if (full_n_elem > 0)
139 {
140 index sub_begin = mat_begin - full_mat_begin;
141 index sub_end = mat_end - full_mat_begin;
142 for (index iC = 0; iC < full_n_elem; iC++)
143 {
144 if (iC >= sub_begin && iC < sub_end)
145 continue; // already remapped above
146 for (auto &iCOther : full_mat_begin[iC])
147 {
148 index iLocal = iCOther - ind_offset;
149 if (iLocal >= 0 && iLocal < n_elem)
150 iCOther = ind_offset + i_old2new[iLocal];
151 }
152 }
153 }
154 }
155 else
156 part_start[1] = n_elem;
157
158 bwOldM = 0, bwNewM = 0;
159 if (!do_rcm)
160 return part_start;
161
162 //! Per-partition RCM reordering. Cross-sub-partition edges are trimmed before RCM
163 //! since bandwidth reduction only considers intra-partition connectivity.
164 //! The reindex step uses (ind_offset + local_offset) as the absolute base.
165 for (int iPart = 0; iPart < nParts; iPart++)
166 {
167 index local_offset = part_start[iPart];
168 index local_nelem = part_start[iPart + 1] - part_start[iPart];
169 //! Trim: remove edges outside [ind_offset + local_offset, ind_offset + local_offset + local_nelem)
170 for (index iC = 0; iC < local_nelem; iC++)
171 {
172 auto &row = mat_begin[iC + local_offset];
173 auto last = std::remove_if(row.begin(), row.end(), [&](index v)
174 { return (v < local_offset + ind_offset) ||
175 (v >= local_offset + ind_offset + local_nelem); });
176 row.erase(last, row.end());
177 }
179 auto [cNew2Old_, cOld2New_] = ReorderSerialAdj_CorrectRCM(
180 mat_begin + local_offset,
181 mat_begin + local_offset + local_nelem,
182 bwOld, bwNew,
183 ind_offset + local_offset);
184 bwNewM = std::max(bwNewM, bwNew);
185 bwOldM = std::max(bwOldM, bwOld);
186
187 tLocalMatStruct new_mat(local_nelem);
188 std::vector<index> new_i_new2old(local_nelem);
189 for (index iC = 0; iC < local_nelem; iC++)
190 {
191 index iNew = cOld2New_[iC];
192 new_i_new2old[iNew] = i_new2old_begin[iC + local_offset];
193 new_mat[iNew] = std::move(mat_begin[iC + local_offset]);
194 }
195 for (auto &row : new_mat)
196 for (auto &iCOther : row)
197 iCOther = ind_offset + local_offset + cOld2New_.at(iCOther - ind_offset - local_offset);
198
199 for (index iC = 0; iC < local_nelem; iC++)
200 {
201 i_new2old_begin[iC + local_offset] = new_i_new2old[iC]; // composing the oldest mapping
202 mat_begin[iC + local_offset] = std::move(new_mat[iC]);
203 }
204 }
205
206 return part_start;
207 }
208}
#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(expr)
Runtime check active in both debug and release builds. Throws std::runtime_error if expr evaluates to...
Definition Errors.hpp:89
index bwNew
Definition Mesh.cpp:2254
index bwOld
Definition Mesh.cpp:2253
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