ryujin 2.1.1 revision 6759a3f00bf045f3527c5e7e7dfd18c7d96a6edb
local_index_handling.h
Go to the documentation of this file.
1//
2// SPDX-License-Identifier: Apache-2.0 WITH LLVM-exception
3// Copyright (C) 2020 - 2022 by the ryujin authors
4//
5
6#pragma once
7
8#include <deal.II/base/partitioner.h>
9#include <deal.II/dofs/dof_handler.h>
10#include <deal.II/dofs/dof_renumbering.h>
11#include <deal.II/dofs/dof_tools.h>
12#include <deal.II/lac/affine_constraints.h>
13#include <deal.II/lac/dynamic_sparsity_pattern.h>
14#include <deal.II/lac/sparsity_tools.h>
15
16namespace ryujin
17{
26 template <typename Number>
28 const dealii::Utilities::MPI::Partitioner &partitioner,
29 dealii::AffineConstraints<Number> &affine_constraints)
30 {
31 affine_constraints.close();
32
33 dealii::AffineConstraints<Number> temporary;
34
35 for (auto line : affine_constraints.get_lines()) {
36 /* translate into local index ranges: */
37 line.index = partitioner.global_to_local(line.index);
38 std::transform(line.entries.begin(),
39 line.entries.end(),
40 line.entries.begin(),
41 [&](auto entry) {
42 return std::make_pair(
43 partitioner.global_to_local(entry.first),
44 entry.second);
45 });
46
47 temporary.add_line(line.index);
48 temporary.add_entries(line.index, line.entries);
49 temporary.set_inhomogeneity(line.index, line.inhomogeneity);
50 }
51
52 temporary.close();
53
54 affine_constraints = std::move(temporary);
55 }
56
57
66 template <typename VECTOR>
68 const dealii::Utilities::MPI::Partitioner &partitioner, VECTOR &vector)
69 {
70 std::transform(
71 vector.begin(), vector.end(), vector.begin(), [&](auto index) {
72 return partitioner.global_to_local(index);
73 });
74 }
75
76
83 namespace DoFRenumbering
84 {
89 using dealii::DoFRenumbering::Cuthill_McKee;
90
100 template <int dim>
101 unsigned int export_indices_first(dealii::DoFHandler<dim> &dof_handler,
102 const MPI_Comm &mpi_communicator,
103 const unsigned int n_locally_internal,
104 const std::size_t group_size)
105 {
106 using namespace dealii;
107
108 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
109 const auto n_locally_owned = locally_owned.n_elements();
110
111 /* The locally owned index range has to be contiguous */
112 Assert(locally_owned.is_contiguous() == true,
113 dealii::ExcMessage(
114 "Need a contiguous set of locally owned indices."));
115
116 /* Offset to translate from global to local index range */
117 const auto offset = n_locally_owned != 0 ? *locally_owned.begin() : 0;
118
119 IndexSet locally_relevant;
120 DoFTools::extract_locally_relevant_dofs(dof_handler, locally_relevant);
121
122 /* Create a temporary MPI partitioner: */
123
124 Utilities::MPI::Partitioner partitioner(
125 locally_owned, locally_relevant, mpi_communicator);
126
127 IndexSet export_indices(n_locally_owned);
128 for (const auto &it : partitioner.import_indices()) {
129 export_indices.add_range(it.first, it.second);
130 }
131
132 std::vector<dealii::types::global_dof_index> new_order(n_locally_owned);
133
134 /*
135 * First pass: reorder all strides containing export indices and mark
136 * all other indices with numbers::invalid_dof_index:
137 */
138
139 unsigned int n_export_indices = 0;
140
141 Assert(n_locally_internal <= n_locally_owned, dealii::ExcInternalError());
142
143 for (unsigned int i = 0; i < n_locally_internal; i += group_size) {
144 bool export_index_present = false;
145 for (unsigned int j = 0; j < group_size; ++j) {
146 if (export_indices.is_element(i + j)) {
147 export_index_present = true;
148 break;
149 }
150 }
151
152 if (export_index_present) {
153 Assert(n_export_indices % group_size == 0,
154 dealii::ExcInternalError());
155 for (unsigned int j = 0; j < group_size; ++j) {
156 new_order[i + j] = offset + n_export_indices++;
157 }
158 } else {
159 for (unsigned int j = 0; j < group_size; ++j)
160 new_order[i + j] = dealii::numbers::invalid_dof_index;
161 }
162 }
163
164#if DEBUG
165 unsigned int n_other = 0;
166 for (unsigned int i = n_locally_internal; i < n_locally_owned; ++i)
167 if (export_indices.is_element(i))
168 n_other++;
169
170 Assert(n_other + n_export_indices >= export_indices.n_elements(),
171 dealii::ExcInternalError());
172#endif
173
174 unsigned int running_index = n_export_indices;
175
176 /*
177 * Second pass: append the rest:
178 */
179
180 for (unsigned int i = 0; i < n_locally_internal; i += group_size) {
181 if (new_order[i] == dealii::numbers::invalid_dof_index) {
182 for (unsigned int j = 0; j < group_size; ++j) {
183 Assert(new_order[i + j] == dealii::numbers::invalid_dof_index,
184 dealii::ExcInternalError());
185 new_order[i + j] = offset + running_index++;
186 }
187 }
188 }
189
190 Assert(running_index == n_locally_internal, dealii::ExcInternalError());
191
192 for (unsigned int i = n_locally_internal; i < n_locally_owned; i++) {
193 new_order[i] = offset + running_index++;
194 }
195
196 Assert(running_index == n_locally_owned, dealii::ExcInternalError());
197
198 dof_handler.renumber_dofs(new_order);
199
200 Assert(n_export_indices % group_size == 0, dealii::ExcInternalError());
201 Assert(n_export_indices <= n_locally_internal,
202 dealii::ExcInternalError());
203 return n_export_indices;
204 }
205
206
213 template <int dim>
214 unsigned int
215 inconsistent_strides_last(dealii::DoFHandler<dim> &dof_handler,
216 const dealii::DynamicSparsityPattern &sparsity,
217 const unsigned int n_locally_internal,
218 const std::size_t group_size)
219 {
220 using namespace dealii;
221
222 const IndexSet &locally_owned = dof_handler.locally_owned_dofs();
223 const auto n_locally_owned = locally_owned.n_elements();
224
225 /* The locally owned index range has to be contiguous */
226 Assert(locally_owned.is_contiguous() == true,
227 dealii::ExcMessage(
228 "Need a contiguous set of locally owned indices."));
229
230 /* Offset to translate from global to local index range */
231 const auto offset = n_locally_owned != 0 ? *locally_owned.begin() : 0;
232
233 std::vector<dealii::types::global_dof_index> new_order(n_locally_owned);
234
235 /*
236 * First pass: keep all strides with consistent row length at the
237 * beginning of the locally internal index range and mark all other
238 * indices with numbers::invalid_dof_index:
239 */
240
241 unsigned int n_consistent_range = 0;
242
243 Assert(n_locally_internal <= n_locally_owned, dealii::ExcInternalError());
244
245 for (unsigned int i = 0; i < n_locally_internal; i += group_size) {
246
247 bool stride_is_consistent = true;
248 const auto group_row_length = sparsity.row_length(offset + i);
249 for (unsigned int j = 0; j < group_size; ++j) {
250 if (group_row_length != sparsity.row_length(offset + i + j)) {
251 stride_is_consistent = false;
252 break;
253 }
254 }
255
256 if (stride_is_consistent) {
257 for (unsigned int j = 0; j < group_size; ++j) {
258 new_order[i + j] = offset + n_consistent_range++;
259 }
260 } else {
261 for (unsigned int j = 0; j < group_size; ++j)
262 new_order[i + j] = dealii::numbers::invalid_dof_index;
263 }
264 }
265
266 /*
267 * Second pass: append the rest:
268 */
269
270 unsigned int running_index = n_consistent_range;
271
272 for (unsigned int i = 0; i < n_locally_internal; i += group_size) {
273 if (new_order[i] == dealii::numbers::invalid_dof_index) {
274 for (unsigned int j = 0; j < group_size; ++j) {
275 Assert(new_order[i + j] == dealii::numbers::invalid_dof_index,
276 dealii::ExcInternalError());
277 new_order[i + j] = offset + running_index++;
278 }
279 }
280 }
281
282 Assert(running_index == n_locally_internal, dealii::ExcInternalError());
283
284 for (unsigned int i = n_locally_internal; i < n_locally_owned; i++) {
285 new_order[i] = offset + running_index++;
286 }
287
288 Assert(running_index == n_locally_owned, dealii::ExcInternalError());
289
290 dof_handler.renumber_dofs(new_order);
291
292 Assert(n_consistent_range % group_size == 0, dealii::ExcInternalError());
293 Assert(n_consistent_range < n_locally_internal,
294 dealii::ExcInternalError());
295 return n_consistent_range;
296 }
297
298
313 template <int dim>
314 unsigned int internal_range(dealii::DoFHandler<dim> &dof_handler,
315 const dealii::DynamicSparsityPattern &sparsity,
316 const std::size_t group_size)
317 {
318 using namespace dealii;
319
320 const auto &locally_owned = dof_handler.locally_owned_dofs();
321 const auto n_locally_owned = locally_owned.n_elements();
322
323 /* The locally owned index range has to be contiguous */
324
325 Assert(locally_owned.is_contiguous() == true,
326 dealii::ExcMessage(
327 "Need a contiguous set of locally owned indices."));
328
329 /* Offset to translate from global to local index range */
330 const auto offset = n_locally_owned != 0 ? *locally_owned.begin() : 0;
331
332 using dof_type = dealii::types::global_dof_index;
333 std::vector<dof_type> new_order(n_locally_owned);
334 dof_type current_index = offset;
335
336 /*
337 * Sort degrees of freedom into a map grouped by stencil size. Write
338 * out dof indices into the new_order vector in groups of group_size
339 * and with same stencil size.
340 */
341
342 std::map<unsigned int, std::set<dof_type>> bins;
343
344 for (unsigned int i = 0; i < n_locally_owned; ++i) {
345 const dof_type index = i;
346 const unsigned int row_length = sparsity.row_length(offset + index);
347 bins[row_length].insert(index);
348
349 if (bins[row_length].size() == group_size) {
350 for (const auto &index : bins[row_length])
351 new_order[index] = current_index++;
352 bins.erase(row_length);
353 }
354 }
355
356 unsigned int n_locally_internal = current_index - offset;
357
358 /* Write out the rest. */
359
360 for (const auto &entries : bins) {
361 Assert(entries.second.size() > 0, ExcInternalError());
362 for (const auto &index : entries.second)
363 new_order[index] = current_index++;
364 }
365 Assert(current_index == offset + n_locally_owned, ExcInternalError());
366
367 dof_handler.renumber_dofs(new_order);
368
369 Assert(n_locally_internal % group_size == 0, ExcInternalError());
370 return n_locally_internal;
371 }
372 } // namespace DoFRenumbering
373
374
381 namespace DoFTools
382 {
384 using dealii::DoFTools::extract_locally_relevant_dofs;
385
387 using dealii::DoFTools::make_hanging_node_constraints;
388
390 using dealii::DoFTools::make_periodicity_constraints;
391
393 using dealii::DoFTools::make_sparsity_pattern;
394
402 template <int dim, typename Number, typename SPARSITY>
404 const dealii::DoFHandler<dim> &dof_handler,
405 SPARSITY &dsp,
406 const dealii::AffineConstraints<Number> &affine_constraints,
407 bool keep_constrained)
408 {
409 const unsigned int dofs_per_cell = dof_handler.get_fe().dofs_per_cell;
410 std::vector<dealii::types::global_dof_index> dof_indices(dofs_per_cell);
411
412 for (auto cell : dof_handler.active_cell_iterators()) {
413 /* iterate over locally owned cells and the ghost layer */
414 if (cell->is_artificial())
415 continue;
416
417 /* translate into local index ranges: */
418 cell->get_dof_indices(dof_indices);
419
420 affine_constraints.add_entries_local_to_global(
421 dof_indices, dsp, keep_constrained);
422 }
423 }
424 } // namespace DoFTools
425
426} // namespace ryujin
void make_extended_sparsity_pattern(const dealii::DoFHandler< dim > &dof_handler, SPARSITY &dsp, const dealii::AffineConstraints< Number > &affine_constraints, bool keep_constrained)
void transform_to_local_range(const dealii::Utilities::MPI::Partitioner &partitioner, dealii::AffineConstraints< Number > &affine_constraints)
unsigned int internal_range(dealii::DoFHandler< dim > &dof_handler, const dealii::DynamicSparsityPattern &sparsity, const std::size_t group_size)
unsigned int export_indices_first(dealii::DoFHandler< dim > &dof_handler, const MPI_Comm &mpi_communicator, const unsigned int n_locally_internal, const std::size_t group_size)
unsigned int inconsistent_strides_last(dealii::DoFHandler< dim > &dof_handler, const dealii::DynamicSparsityPattern &sparsity, const unsigned int n_locally_internal, const std::size_t group_size)