71 template <
int dim,
typename Number =
double>
83 double downstream_size = 0;
85 for (
unsigned int d = 0; d < dim; ++d)
87 downstream_size += (rhs[d] - lhs[d]) * weight;
90 if (downstream_size < 0)
92 else if (downstream_size > 0)
96 for (
unsigned int d = 0; d < dim; ++d)
100 return lhs[d] < rhs[d];
120 template <
int dim,
int spacedim>
121 std::vector<unsigned char>
125 std::vector<unsigned char> local_component_association(
134 local_component_association[i] =
145 const unsigned int first_comp =
150 local_component_association[i] = first_comp;
157 for (
unsigned int c = first_comp; c < fe.
n_components(); ++c)
158 if (component_mask[c] ==
true)
160 local_component_association[i] = c;
165 Assert(std::find(local_component_association.begin(),
166 local_component_association.end(),
167 static_cast<unsigned char>(-1)) ==
168 local_component_association.end(),
171 return local_component_association;
191 template <
int dim,
int spacedim>
196 std::vector<unsigned char> &dofs_by_component,
199 const ::hp::FECollection<dim, spacedim> &fe_collection =
203 const auto &locally_owned_dofs =
209 locally_owned_dofs.n_elements());
218 std::vector<std::vector<unsigned char>> local_component_association(
219 fe_collection.size());
220 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
223 local_component_association[f] =
228 std::vector<types::global_dof_index> indices;
230 const auto runner = [&](
const auto &task) {
236 indices.resize(cell->get_fe().n_dofs_per_cell());
237 cell->get_dof_indices(indices);
244 for (
const auto &cell :
246 if (cell->is_locally_owned_on_level())
248 indices.resize(cell->get_fe().n_dofs_per_cell());
249 cell->get_mg_dof_indices(indices);
256 runner([&](
const auto &cell) {
258 for (
unsigned int i = 0; i < cell->get_fe().n_dofs_per_cell(); ++i)
259 if (locally_owned_dofs.is_element(indices[i]))
260 dofs_by_component[locally_owned_dofs.index_within_set(indices[i])] =
261 local_component_association[fe_index][i];
273 template <
int dim,
int spacedim>
276 std::vector<unsigned char> &dofs_by_block)
278 const ::hp::FECollection<dim, spacedim> &fe_collection =
292 std::vector<std::vector<unsigned char>> local_block_association(
293 fe_collection.size());
294 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
298 static_cast<unsigned char>(-1));
302 Assert(std::find(local_block_association[f].begin(),
303 local_block_association[f].end(),
304 static_cast<unsigned char>(-1)) ==
305 local_block_association[f].end(),
310 std::vector<types::global_dof_index> indices;
312 if (cell->is_locally_owned())
315 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
316 indices.resize(dofs_per_cell);
317 cell->get_dof_indices(indices);
318 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
321 indices[i])] = local_block_association[fe_index][i];
328 template <
int dim,
int spacedim,
typename Number>
333 const unsigned int component)
342 for (
unsigned int i = 0; i < fe_collection.size(); ++i)
344 Assert(fe_collection[i].is_primitive() ==
true,
351 const bool consider_components =
355 if (consider_components ==
false)
359 std::vector<unsigned char> component_dofs(
366 for (
unsigned int i = 0; i < dof_data.
size(); ++i)
367 if (component_dofs[i] ==
static_cast<unsigned char>(component))
372 std::vector<unsigned char> touch_count(dof_handler.
n_dofs(), 0);
374 std::vector<types::global_dof_index> dof_indices;
375 dof_indices.reserve(fe_collection.max_dofs_per_cell());
379 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
380 dof_indices.resize(dofs_per_cell);
381 cell->get_dof_indices(dof_indices);
383 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
386 if (!consider_components ||
387 (cell->get_fe().system_to_component_index(i).first == component))
390 dof_data(dof_indices[i]) += cell_data(cell->active_cell_index());
393 ++touch_count[dof_indices[i]];
403 Assert(consider_components || (touch_count[i] != 0),
405 if (touch_count[i] != 0)
406 dof_data(i) /= touch_count[i];
412 template <
int dim,
int spacedim>
420 "The given component mask is not sized correctly to represent the "
421 "components of the given finite element."));
439 std::vector<types::global_dof_index> selected_dofs;
442 if (component_mask[dofs_by_component[i]] ==
true)
447 result.
add_indices(selected_dofs.begin(), selected_dofs.end());
453 template <
int dim,
int spacedim>
465 template <
int dim,
int spacedim>
466 std::vector<IndexSet>
473 "The given component mask is not sized correctly to represent the "
474 "components of the given finite element."));
483 std::vector<IndexSet> index_per_comp(n_comps,
IndexSet(dof.
n_dofs()));
487 const auto &comp_i = dofs_by_component[i];
488 if (component_mask[comp_i])
489 index_per_comp[comp_i].add_index(
490 locally_owned_dofs.nth_index_in_set(i));
492 for (
const auto &c : index_per_comp)
494 return index_per_comp;
499 template <
int dim,
int spacedim>
504 std::vector<bool> &selected_dofs)
511 "The given component mask is not sized correctly to represent the "
512 "components of the given finite element."));
521 std::fill_n(selected_dofs.begin(), dof.
n_dofs(
level),
false);
528 std::fill_n(selected_dofs.begin(), dof.
n_dofs(
level),
true);
533 std::fill_n(selected_dofs.begin(), dof.
n_dofs(
level),
false);
537 std::vector<unsigned char> local_component_association =
541 local_selected_dofs[i] = component_mask[local_component_association[i]];
547 cell->get_mg_dof_indices(indices);
549 selected_dofs[indices[i]] = local_selected_dofs[i];
555 template <
int dim,
int spacedim>
560 std::vector<bool> &selected_dofs)
571 template <
int dim,
int spacedim>
575 std::vector<bool> &selected_dofs,
576 const std::set<types::boundary_id> &boundary_ids)
582 "This function can not be used with distributed triangulations. "
583 "See the documentation for more information."));
589 selected_dofs.clear();
590 selected_dofs.resize(dof_handler.
n_dofs(),
false);
598 template <
int dim,
int spacedim>
603 const std::set<types::boundary_id> &boundary_ids)
612 template <
int dim,
int spacedim>
616 const std::set<types::boundary_id> &boundary_ids)
620 ExcMessage(
"Component mask has invalid size."));
629 const bool check_boundary_id = (boundary_ids.size() != 0);
633 const bool check_vector_component =
639 std::vector<types::global_dof_index> face_dof_indices;
640 face_dof_indices.reserve(
652 if (cell->is_artificial() ==
false)
653 for (
const unsigned int face : cell->face_indices())
654 if (cell->at_boundary(face))
655 if (!check_boundary_id ||
656 (boundary_ids.find(cell->face(face)->boundary_id()) !=
661 const auto reference_cell = cell->reference_cell();
663 const unsigned int n_vertices_per_cell =
664 reference_cell.n_vertices();
665 const unsigned int n_lines_per_cell = reference_cell.n_lines();
666 const unsigned int n_vertices_per_face =
667 reference_cell.face_reference_cell(face).n_vertices();
668 const unsigned int n_lines_per_face =
669 reference_cell.face_reference_cell(face).n_lines();
672 face_dof_indices.resize(dofs_per_face);
673 cell->face(face)->get_dof_indices(face_dof_indices,
674 cell->active_fe_index());
677 if (!check_vector_component)
678 selected_dofs.
add_index(face_dof_indices[i]);
693 (dim == 3 ? (i < n_vertices_per_face *
696 (i < n_vertices_per_face *
700 (i - n_vertices_per_face *
702 n_vertices_per_cell *
705 n_vertices_per_face *
709 n_vertices_per_cell *
719 selected_dofs.
add_index(face_dof_indices[i]);
723 const unsigned int first_nonzero_comp =
729 if (component_mask[first_nonzero_comp] ==
true)
730 selected_dofs.
add_index(face_dof_indices[i]);
735 return selected_dofs;
740 template <
int dim,
int spacedim>
745 std::vector<bool> &selected_dofs,
746 const std::set<types::boundary_id> &boundary_ids)
750 ExcMessage(
"This component mask has the wrong size."));
757 const bool check_boundary_id = (boundary_ids.size() != 0);
761 const bool check_vector_component =
765 selected_dofs.clear();
766 selected_dofs.resize(dof_handler.
n_dofs(),
false);
767 std::vector<types::global_dof_index> cell_dof_indices;
768 cell_dof_indices.reserve(
778 for (
const unsigned int face : cell->face_indices())
779 if (cell->at_boundary(face))
780 if (!check_boundary_id ||
781 (boundary_ids.find(cell->face(face)->boundary_id()) !=
787 cell_dof_indices.resize(dofs_per_cell);
788 cell->get_dof_indices(cell_dof_indices);
793 if (!check_vector_component)
794 selected_dofs[cell_dof_indices[i]] =
true;
801 selected_dofs[cell_dof_indices[i]] =
806 const unsigned int first_nonzero_comp =
812 selected_dofs[cell_dof_indices[i]] =
813 (component_mask[first_nonzero_comp] ==
true);
822 template <
int dim,
int spacedim,
typename number>
831 const std::function<
bool(
836 ->
bool {
return cell->is_locally_owned() && predicate(cell); };
838 std::vector<types::global_dof_index> local_dof_indices;
839 local_dof_indices.reserve(
843 std::set<types::global_dof_index> predicate_dofs;
846 if (!cell->is_artificial() && predicate(cell))
848 local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
849 cell->get_dof_indices(local_dof_indices);
850 predicate_dofs.insert(local_dof_indices.begin(),
851 local_dof_indices.end());
855 std::set<types::global_dof_index> dofs_with_support_on_halo_cells;
857 const std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>
861 active_cell_iterator>::const_iterator it =
863 it != halo_cells.end();
875 const unsigned int dofs_per_cell = (*it)->get_fe().n_dofs_per_cell();
876 local_dof_indices.resize(dofs_per_cell);
877 (*it)->get_dof_indices(local_dof_indices);
878 dofs_with_support_on_halo_cells.insert(local_dof_indices.begin(),
879 local_dof_indices.end());
891 std::set<types::global_dof_index> extra_halo;
892 for (
const auto dof : dofs_with_support_on_halo_cells)
898 const unsigned int line_size = line_ptr->size();
899 for (
unsigned int j = 0; j < line_size; ++j)
900 extra_halo.insert((*line_ptr)[j].first);
903 dofs_with_support_on_halo_cells.insert(extra_halo.begin(),
910 support_set.
add_indices(predicate_dofs.begin(), predicate_dofs.end());
914 halo_set.
add_indices(dofs_with_support_on_halo_cells.begin(),
915 dofs_with_support_on_halo_cells.end());
930 template <
int spacedim>
939 template <
int spacedim>
943 const unsigned int dim = 2;
951 for (
const auto &cell : dof_handler.active_cell_iterators())
952 if (!cell->is_artificial())
954 for (
const unsigned int face : cell->face_indices())
955 if (cell->face(face)->has_children())
958 line = cell->face(face);
962 selected_dofs.add_index(
963 line->child(0)->vertex_dof_index(1, dof));
965 for (
unsigned int child = 0; child < 2; ++child)
967 if (cell->neighbor_child_on_subface(face, child)
972 selected_dofs.add_index(
973 line->child(child)->dof_index(dof));
978 selected_dofs.compress();
979 return selected_dofs;
983 template <
int spacedim>
987 const unsigned int dim = 3;
994 for (
const auto &cell : dof_handler.active_cell_iterators())
995 if (!cell->is_artificial())
996 for (auto f : cell->face_indices())
1000 if (cell->face(f)->has_children())
1002 for (
unsigned int child = 0; child < 4; ++child)
1003 if (!cell->neighbor_child_on_subface(f, child)
1007 std::vector<types::global_dof_index> ldi(
1009 face->child(child)->get_dof_indices(ldi);
1010 selected_dofs.add_indices(ldi.begin(), ldi.end());
1015 for (
unsigned int vertex = 0; vertex < 4; ++vertex)
1018 unconstrained_dofs.add_index(
1019 face->vertex_dof_index(vertex, dof));
1022 selected_dofs.subtract_set(unconstrained_dofs);
1023 return selected_dofs;
1030 template <
int dim,
int spacedim>
1034 return internal::extract_hanging_node_dofs(dof_handler);
1039 template <
int dim,
int spacedim>
1043 std::vector<bool> &selected_dofs)
1049 std::fill_n(selected_dofs.begin(), dof_handler.
n_dofs(),
false);
1051 std::vector<types::global_dof_index> local_dof_indices;
1052 local_dof_indices.reserve(
1058 if (cell->subdomain_id() == subdomain_id)
1060 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1061 local_dof_indices.resize(dofs_per_cell);
1062 cell->get_dof_indices(local_dof_indices);
1063 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1064 selected_dofs[local_dof_indices[i]] =
true;
1070 template <
int dim,
int spacedim>
1080 std::vector<types::global_dof_index> dof_indices;
1081 std::set<types::global_dof_index> global_dof_indices;
1086 dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1087 cell->get_dof_indices(dof_indices);
1091 global_dof_indices.insert(dof_index);
1094 dof_set.
add_indices(global_dof_indices.begin(), global_dof_indices.end());
1103 template <
int dim,
int spacedim>
1113 template <
int dim,
int spacedim>
1117 const unsigned int level)
1125 std::vector<types::global_dof_index> dof_indices;
1126 std::set<types::global_dof_index> global_dof_indices;
1128 const auto filtered_iterators_range =
1131 for (
const auto &cell : filtered_iterators_range)
1133 dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1134 cell->get_mg_dof_indices(dof_indices);
1138 global_dof_indices.insert(dof_index);
1141 dof_set.
add_indices(global_dof_indices.begin(), global_dof_indices.end());
1150 template <
int dim,
int spacedim>
1155 const unsigned int level)
1162 template <
int dim,
int spacedim>
1177 std::vector<types::global_dof_index> dof_indices;
1178 std::vector<types::global_dof_index> dofs_on_ghosts;
1181 if (cell->is_ghost())
1183 dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1184 cell->get_dof_indices(dof_indices);
1185 for (
const auto dof_index : dof_indices)
1187 dofs_on_ghosts.push_back(dof_index);
1191 std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1192 dof_set.
add_indices(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1200 template <
int dim,
int spacedim>
1210 template <
int dim,
int spacedim>
1214 const unsigned int level)
1227 std::vector<types::global_dof_index> dof_indices;
1228 std::vector<types::global_dof_index> dofs_on_ghosts;
1239 dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1240 cell->get_mg_dof_indices(dof_indices);
1241 for (
const auto dof_index : dof_indices)
1243 dofs_on_ghosts.push_back(dof_index);
1247 std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1248 dof_set.
add_indices(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1256 template <
int dim,
int spacedim>
1260 const unsigned int level,
1269 template <
int dim,
int spacedim>
1273 const unsigned int mg_level,
1274 std::vector<std::vector<bool>> &constant_modes)
1276 const auto &locally_owned_dofs =
1283 if (locally_owned_dofs.n_elements() == 0)
1285 constant_modes = std::vector<std::vector<bool>>(0);
1293 std::vector<unsigned char> dofs_by_component(
1294 locally_owned_dofs.n_elements());
1299 unsigned int n_selected_dofs = 0;
1300 for (
unsigned int i = 0; i < n_components; ++i)
1301 if (component_mask[i] ==
true)
1303 std::count(dofs_by_component.begin(), dofs_by_component.end(), i);
1306 std::vector<unsigned int> component_numbering(
1308 for (
unsigned int i = 0, count = 0; i < locally_owned_dofs.n_elements();
1310 if (component_mask[dofs_by_component[i]])
1311 component_numbering[i] = count++;
1318 const ::hp::FECollection<dim, spacedim> &fe_collection =
1320 std::vector<Table<2, bool>> element_constant_modes;
1321 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
1322 constant_mode_to_component_translation(n_components);
1324 unsigned int n_constant_modes = 0;
1325 int first_non_empty_constant_mode = -1;
1326 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
1328 std::pair<Table<2, bool>, std::vector<unsigned int>> data =
1329 fe_collection[f].get_constant_modes();
1333 if (first_non_empty_constant_mode < 0 && data.first.n_rows() > 0)
1335 first_non_empty_constant_mode = f;
1339 for (
unsigned int i = 0; i < data.second.size(); ++i)
1340 if (component_mask[data.second[i]])
1341 constant_mode_to_component_translation[data.second[i]]
1342 .emplace_back(n_constant_modes++, i);
1348 element_constant_modes.push_back(data.first);
1349 Assert(element_constant_modes.back().n_rows() == 0 ||
1350 element_constant_modes.back().n_rows() ==
1351 element_constant_modes[first_non_empty_constant_mode]
1359 constant_modes.clear();
1360 constant_modes.resize(n_constant_modes,
1361 std::vector<bool>(n_selected_dofs,
false));
1365 std::vector<types::global_dof_index> dof_indices;
1367 const auto runner = [&](
const auto &task) {
1373 dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1374 cell->get_dof_indices(dof_indices);
1381 for (
const auto &cell :
1383 if (cell->is_locally_owned_on_level())
1385 dof_indices.resize(cell->get_fe().n_dofs_per_cell());
1386 cell->get_mg_dof_indices(dof_indices);
1393 runner([&](
const auto &cell) {
1394 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
1395 if (locally_owned_dofs.is_element(dof_indices[i]))
1397 const unsigned int loc_index =
1398 locally_owned_dofs.index_within_set(dof_indices[i]);
1399 const unsigned int comp = dofs_by_component[loc_index];
1400 if (component_mask[comp])
1401 for (auto &indices :
1402 constant_mode_to_component_translation[comp])
1403 constant_modes[indices
1404 .first][component_numbering[loc_index]] =
1405 element_constant_modes[cell->active_fe_index()](
1414 template <
int dim,
int spacedim>
1418 std::vector<std::vector<bool>> &constant_modes)
1428 template <
int dim,
int spacedim>
1433 std::vector<std::vector<bool>> &constant_modes)
1443 template <
int dim,
int spacedim>
1444 std::map<
typename DoFHandler<dim - 1, spacedim>::active_cell_iterator,
1445 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1456 std::map<
typename DoFHandler<dim - 1, spacedim>::active_cell_iterator,
1457 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1459 c1_to_c0_cell_and_face;
1462 if (c1_to_c0_face.empty())
1463 return c1_to_c0_cell_and_face;
1466 std::map<typename Triangulation<dim, spacedim>::face_iterator,
1467 typename DoFHandler<dim - 1, spacedim>::active_cell_iterator>
1471 for (
const auto &[c1_cell, c0_cell] : c1_to_c0_face)
1472 if (!c1_cell->has_children())
1473 c0_to_c1[c0_cell] = c1_cell->as_dof_handler_iterator(c1_dh);
1477 for (
const auto &cell :
1479 for (
const auto f : cell->face_indices())
1480 if (cell->face(f)->at_boundary())
1482 const auto &it = c0_to_c1.find(cell->face(f));
1483 if (it != c0_to_c1.end())
1485 const auto &c1_cell = it->second;
1486 c1_to_c0_cell_and_face[c1_cell] = {cell, f};
1492 return c1_to_c0_cell_and_face;
1497 template <
int dim,
int spacedim>
1500 std::vector<unsigned int> &active_fe_indices)
1507 active_fe_indices.assign(indices.begin(), indices.end());
1510 template <
int dim,
int spacedim>
1511 std::vector<IndexSet>
1515 ExcMessage(
"The given DoFHandler has no DoFs."));
1523 "For parallel::distributed::Triangulation objects and "
1524 "associated DoF handler objects, asking for any information "
1525 "related to a subdomain other than the locally owned one does "
1526 "not make sense."));
1530 std::vector<::types::subdomain_id> subdomain_association(
1533 subdomain_association);
1547 const unsigned int n_subdomains =
1551 unsigned int max_subdomain_id = 0;
1554 std::max(max_subdomain_id, cell->subdomain_id());
1555 return max_subdomain_id + 1;
1561 Assert(n_subdomains > *std::max_element(subdomain_association.begin(),
1562 subdomain_association.end()),
1565 std::vector<::IndexSet> index_sets(
1574 index < subdomain_association.size();
1578 if (subdomain_association[index] != this_subdomain)
1580 index_sets[this_subdomain].add_range(i_min, index);
1582 this_subdomain = subdomain_association[index];
1587 if (i_min == subdomain_association.size() - 1)
1589 index_sets[this_subdomain].add_index(i_min);
1595 index_sets[this_subdomain].add_range(i_min,
1596 subdomain_association.size());
1599 for (
unsigned int i = 0; i < n_subdomains; ++i)
1600 index_sets[i].compress();
1605 template <
int dim,
int spacedim>
1606 std::vector<IndexSet>
1616 "For parallel::distributed::Triangulation objects and "
1617 "associated DoF handler objects, asking for any information "
1618 "related to a subdomain other than the locally owned one does "
1619 "not make sense."));
1627 std::vector<IndexSet> dof_set =
1629 const ::types::subdomain_id n_subdomains = dof_set.size();
1635 subdomain_id < n_subdomains;
1648 std::vector<types::global_dof_index> local_dof_indices;
1649 std::set<types::global_dof_index> subdomain_halo_global_dof_indices;
1650 for (
typename std::vector<
1652 const_iterator it_cell = active_halo_layer.begin();
1653 it_cell != active_halo_layer.end();
1659 cell->subdomain_id() != subdomain_id,
1661 "The subdomain ID of the halo cell should not match that of the vector entry."));
1664 cell->get_dof_indices(local_dof_indices);
1668 subdomain_halo_global_dof_indices.insert(local_dof_index);
1671 dof_set[subdomain_id].add_indices(
1672 subdomain_halo_global_dof_indices.begin(),
1673 subdomain_halo_global_dof_indices.end());
1675 dof_set[subdomain_id].compress();
1681 template <
int dim,
int spacedim>
1685 std::vector<types::subdomain_id> &subdomain_association)
1693 "For parallel::distributed::Triangulation objects and "
1694 "associated DoF handler objects, asking for any subdomain other "
1695 "than the locally owned one does not make sense."));
1697 Assert(subdomain_association.size() == dof_handler.
n_dofs(),
1711 std::vector<types::subdomain_id> cell_owners(
1717 cell_owners = tr->get_true_subdomain_ids_of_cells();
1718 Assert(tr->get_true_subdomain_ids_of_cells().size() ==
1719 tr->n_active_cells(),
1726 cell_owners[cell->active_cell_index()] = cell->subdomain_id();
1730 std::fill_n(subdomain_association.begin(),
1734 std::vector<types::global_dof_index> local_dof_indices;
1735 local_dof_indices.reserve(
1743 cell_owners[cell->active_cell_index()];
1744 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1745 local_dof_indices.resize(dofs_per_cell);
1746 cell->get_dof_indices(local_dof_indices);
1752 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1753 if (subdomain_association[local_dof_indices[i]] ==
1755 subdomain_association[local_dof_indices[i]] = subdomain_id;
1756 else if (subdomain_association[local_dof_indices[i]] > subdomain_id)
1758 subdomain_association[local_dof_indices[i]] = subdomain_id;
1762 Assert(std::find(subdomain_association.begin(),
1763 subdomain_association.end(),
1765 subdomain_association.end(),
1771 template <
int dim,
int spacedim>
1777 std::vector<types::subdomain_id> subdomain_association(
1781 return std::count(subdomain_association.begin(),
1782 subdomain_association.end(),
1788 template <
int dim,
int spacedim>
1801 "For parallel::distributed::Triangulation objects and "
1802 "associated DoF handler objects, asking for any subdomain other "
1803 "than the locally owned one does not make sense."));
1807 std::vector<types::global_dof_index> local_dof_indices;
1808 local_dof_indices.reserve(
1816 std::vector<types::global_dof_index> subdomain_indices;
1819 if ((cell->is_artificial() ==
false) &&
1820 (cell->subdomain_id() == subdomain))
1822 const unsigned int dofs_per_cell = cell->get_fe().n_dofs_per_cell();
1823 local_dof_indices.resize(dofs_per_cell);
1824 cell->get_dof_indices(local_dof_indices);
1825 subdomain_indices.insert(subdomain_indices.end(),
1826 local_dof_indices.begin(),
1827 local_dof_indices.end());
1830 std::sort(subdomain_indices.begin(), subdomain_indices.end());
1831 index_set.
add_indices(subdomain_indices.begin(), subdomain_indices.end());
1839 template <
int dim,
int spacedim>
1844 std::vector<unsigned int> &n_dofs_on_subdomain)
1849 std::fill(n_dofs_on_subdomain.begin(), n_dofs_on_subdomain.end(), 0);
1858 return cell.subdomain_id() == subdomain;
1860 ExcMessage(
"There are no cells for the given subdomain!"));
1862 std::vector<types::subdomain_id> subdomain_association(
1866 std::vector<unsigned char> component_association(dof_handler.
n_dofs());
1869 component_association);
1871 for (
unsigned int c = 0; c < dof_handler.
get_fe(0).n_components(); ++c)
1874 if ((subdomain_association[i] == subdomain) &&
1875 (component_association[i] ==
static_cast<unsigned char>(c)))
1876 ++n_dofs_on_subdomain[c];
1888 template <
int dim,
int spacedim>
1891 const std::vector<unsigned char> &dofs_by_component,
1892 const std::vector<unsigned int> &target_component,
1893 const bool only_once,
1894 std::vector<types::global_dof_index> &dofs_per_component,
1895 unsigned int &component)
1914 for (
unsigned int dd = 0; dd < d; ++dd, ++component)
1915 dofs_per_component[target_component[component]] +=
1916 std::count(dofs_by_component.begin(),
1917 dofs_by_component.end(),
1924 for (
unsigned int dd = 1; dd < d; ++dd)
1925 dofs_per_component[target_component[component - d + dd]] =
1926 dofs_per_component[target_component[component - d]];
1933 template <
int dim,
int spacedim>
1936 const std::vector<unsigned char> &dofs_by_component,
1937 const std::vector<unsigned int> &target_component,
1938 const bool only_once,
1939 std::vector<types::global_dof_index> &dofs_per_component,
1940 unsigned int &component)
1945 for (
unsigned int fe = 1; fe < fe_collection.
size(); ++fe)
1947 Assert(fe_collection[fe].n_components() ==
1948 fe_collection[0].n_components(),
1950 Assert(fe_collection[fe].n_base_elements() ==
1951 fe_collection[0].n_base_elements(),
1953 for (
unsigned int b = 0; b < fe_collection[0].n_base_elements(); ++b)
1955 Assert(fe_collection[fe].base_element(b).n_components() ==
1956 fe_collection[0].base_element(b).n_components(),
1958 Assert(fe_collection[fe].base_element(b).n_base_elements() ==
1959 fe_collection[0].base_element(b).n_base_elements(),
1982 template <
int dim,
int spacedim>
1994 template <
int dim,
int spacedim>
1996 all_elements_are_primitive(
1997 const ::hp::FECollection<dim, spacedim> &fe_collection)
1999 for (
unsigned int i = 0; i < fe_collection.size(); ++i)
2000 if (fe_collection[i].is_primitive() ==
false)
2010 template <
int dim,
int spacedim>
2011 std::vector<types::global_dof_index>
2014 const bool only_once,
2015 const std::vector<unsigned int> &target_component_)
2021 std::vector<unsigned int> target_component = target_component_;
2022 if (target_component.empty())
2024 target_component.resize(n_components);
2025 for (
unsigned int i = 0; i < n_components; ++i)
2026 target_component[i] = i;
2029 Assert(target_component.size() == n_components,
2033 const unsigned int max_component =
2034 *std::max_element(target_component.begin(), target_component.end());
2035 const unsigned int n_target_components = max_component + 1;
2037 std::vector<types::global_dof_index> dofs_per_component(
2042 if (n_components == 1)
2045 return dofs_per_component;
2051 std::vector<unsigned char> dofs_by_component(
2058 unsigned int component = 0;
2069 Assert((internal::all_elements_are_primitive(
2071 (std::accumulate(dofs_per_component.begin(),
2072 dofs_per_component.end(),
2078#ifdef DEAL_II_WITH_MPI
2084 std::vector<types::global_dof_index> local_dof_count =
2087 const int ierr = MPI_Allreduce(local_dof_count.data(),
2088 dofs_per_component.data(),
2089 n_target_components,
2092 tria->get_communicator());
2097 return dofs_per_component;
2102 template <
int dim,
int spacedim>
2103 std::vector<types::global_dof_index>
2105 const std::vector<unsigned int> &target_block_)
2107 const ::hp::FECollection<dim, spacedim> &fe_collection =
2116 const unsigned int n_blocks = fe_collection[0].n_blocks();
2118 std::vector<unsigned int> target_block = target_block_;
2119 if (target_block.empty())
2121 target_block.resize(fe_collection[0].n_blocks());
2122 for (
unsigned int i = 0; i < n_blocks; ++i)
2123 target_block[i] = i;
2126 Assert(target_block.size() == n_blocks,
2128 for (
unsigned int f = 1; f < fe_collection.size(); ++f)
2129 Assert(fe_collection[0].n_blocks() == fe_collection[f].n_blocks(),
2130 ExcMessage(
"This function can only work if all elements in a "
2131 "collection have the same number of blocks."));
2137 std::vector<types::global_dof_index> dofs_per_block(1);
2138 dofs_per_block[0] = dof_handler.
n_dofs();
2139 return dofs_per_block;
2143 const unsigned int max_block =
2144 *std::max_element(target_block.begin(), target_block.end());
2145 const unsigned int n_target_blocks = max_block + 1;
2147 std::vector<types::global_dof_index> dofs_per_block(n_target_blocks);
2151 for (
unsigned int this_fe = fe_collection.size() - 1;
2152 this_fe < fe_collection.size();
2157 std::vector<unsigned char> dofs_by_block(
2162 for (
unsigned int block = 0; block < fe.
n_blocks(); ++block)
2163 dofs_per_block[target_block[block]] +=
2164 std::count(dofs_by_block.begin(), dofs_by_block.end(), block);
2166#ifdef DEAL_II_WITH_MPI
2173 std::vector<types::global_dof_index> local_dof_count =
2175 const int ierr = MPI_Allreduce(local_dof_count.data(),
2176 dofs_per_block.data(),
2180 tria->get_communicator());
2186 return dofs_per_block;
2191 template <
int dim,
int spacedim>
2194 std::vector<types::global_dof_index> &mapping)
2197 mapping.insert(mapping.end(),
2201 std::vector<types::global_dof_index> dofs_on_face;
2212 for (
const unsigned int f : cell->face_indices())
2213 if (cell->at_boundary(f))
2215 const unsigned int dofs_per_face =
2216 cell->get_fe().n_dofs_per_face(f);
2217 dofs_on_face.resize(dofs_per_face);
2218 cell->face(f)->get_dof_indices(dofs_on_face,
2219 cell->active_fe_index());
2220 for (
unsigned int i = 0; i < dofs_per_face; ++i)
2222 mapping[dofs_on_face[i]] = next_boundary_index++;
2230 template <
int dim,
int spacedim>
2233 const std::set<types::boundary_id> &boundary_ids,
2234 std::vector<types::global_dof_index> &mapping)
2241 mapping.insert(mapping.end(),
2246 if (boundary_ids.empty())
2249 std::vector<types::global_dof_index> dofs_on_face;
2254 for (
const unsigned int f : cell->face_indices())
2255 if (boundary_ids.find(cell->face(f)->boundary_id()) !=
2258 const unsigned int dofs_per_face =
2259 cell->get_fe().n_dofs_per_face(f);
2260 dofs_on_face.resize(dofs_per_face);
2261 cell->face(f)->get_dof_indices(dofs_on_face,
2262 cell->active_fe_index());
2263 for (
unsigned int i = 0; i < dofs_per_face; ++i)
2265 mapping[dofs_on_face[i]] = next_boundary_index++;
2276 template <
int dim,
int spacedim>
2277 std::map<types::global_dof_index, Point<spacedim>>
2283 std::map<types::global_dof_index, Point<spacedim>> support_points;
2289 for (
unsigned int fe_index = 0; fe_index < fe_collection.
size();
2293 Assert((fe_collection[fe_index].n_dofs_per_cell() == 0) ||
2294 (fe_collection[fe_index].has_support_points()),
2297 fe_collection[fe_index].get_unit_support_points()));
2302 (in_mask.
size() == 0 ?
2318 std::vector<types::global_dof_index> local_dof_indices;
2319 for (
const auto &cell : dof_handler.active_cell_iterators())
2321 if (cell->is_artificial() == false)
2323 hp_fe_values.reinit(cell);
2327 local_dof_indices.resize(cell->get_fe().n_dofs_per_cell());
2328 cell->get_dof_indices(local_dof_indices);
2330 const std::vector<Point<spacedim>> &points =
2332 for (
unsigned int i = 0; i < cell->get_fe().n_dofs_per_cell();
2335 const unsigned int dof_comp =
2336 cell->get_fe().system_to_component_index(i).first;
2340 support_points[local_dof_indices[i]] = points[i];
2344 return support_points;
2348 template <
int dim,
int spacedim>
2349 std::vector<Point<spacedim>>
2350 map_dofs_to_support_points_vector(
2355 std::vector<Point<spacedim>> support_points(dof_handler.
n_dofs());
2358 const std::map<types::global_dof_index, Point<spacedim>>
2366 Assert(x_support_points.find(i) != x_support_points.end(),
2369 support_points[i] = x_support_points.find(i)->second;
2372 return support_points;
2378 template <
int dim,
int spacedim>
2390 "This function can not be used with distributed triangulations. "
2391 "See the documentation for more information."));
2398 internal::map_dofs_to_support_points_vector(mapping_collection,
2404 template <
int dim,
int spacedim>
2417 "This function can not be used with distributed triangulations. "
2418 "See the documentation for more information."));
2423 internal::map_dofs_to_support_points_vector(mapping, dof_handler, mask);
2428 template <
int dim,
int spacedim>
2436 support_points.clear();
2442 support_points = internal::map_dofs_to_support_points(mapping_collection,
2449 template <
int dim,
int spacedim>
2457 support_points.clear();
2462 internal::map_dofs_to_support_points(mapping, dof_handler, mask);
2466 template <
int dim,
int spacedim>
2467 std::map<types::global_dof_index, Point<spacedim>>
2476 return internal::map_dofs_to_support_points(mapping_collection,
2482 template <
int dim,
int spacedim>
2483 std::map<types::global_dof_index, Point<spacedim>>
2489 return internal::map_dofs_to_support_points(mapping, dof_handler, mask);
2493 template <
int spacedim>
2501 std::map<Point<spacedim>,
2502 std::vector<types::global_dof_index>,
2507 for (
const auto &it : support_points)
2509 std::vector<types::global_dof_index> &v = point_map[it.second];
2510 v.push_back(it.first);
2514 for (
const auto &it : point_map)
2516 out << it.first <<
" \"";
2517 const std::vector<types::global_dof_index> &v = it.second;
2518 for (
unsigned int i = 0; i < v.size(); ++i)
2532 template <
int dim,
int spacedim>
2541 const unsigned int nb = fe.
n_blocks();
2543 tables_by_block.resize(1);
2544 tables_by_block[0].reinit(nb, nb);
2545 tables_by_block[0].fill(
none);
2553 tables_by_block[0](ib, jb) |= table(i, j);
2561 tables_by_block.resize(fe_collection.
size());
2563 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
2567 const unsigned int nb = fe.
n_blocks();
2568 tables_by_block[f].reinit(nb, nb);
2569 tables_by_block[f].fill(
none);
2576 tables_by_block[f](ib, jb) |= table(i, j);
2585 template <
int dim,
int spacedim>
2589 const unsigned int level,
2590 const std::vector<bool> &selected_dofs,
2593 std::vector<types::global_dof_index> indices;
2598 if (cell->is_locally_owned_on_level())
2605 if (cell->is_locally_owned_on_level())
2607 indices.resize(cell->get_fe().n_dofs_per_cell());
2608 cell->get_mg_dof_indices(indices);
2610 if (selected_dofs.size() != 0)
2615 if (selected_dofs.empty())
2616 block_list.
add(i, indices[j] - offset);
2619 if (selected_dofs[j])
2620 block_list.
add(i, indices[j] - offset);
2628 template <
int dim,
int spacedim>
2632 const unsigned int level,
2633 const bool interior_only)
2638 std::vector<types::global_dof_index> indices;
2639 std::vector<bool> exclude;
2643 indices.resize(cell->get_fe().n_dofs_per_cell());
2644 cell->get_mg_dof_indices(indices);
2650 std::fill(exclude.begin(), exclude.end(),
false);
2652 for (
const unsigned int face : cell->face_indices())
2653 if (cell->at_boundary(face) ||
2654 cell->neighbor(face)->level() != cell->level())
2659 block_list.
add(0, indices[j]);
2663 for (
const auto index : indices)
2664 block_list.
add(0, index);
2670 template <
int dim,
int spacedim>
2674 const unsigned int level,
2675 const bool interior_dofs_only,
2676 const bool boundary_dofs)
2681 std::vector<types::global_dof_index> indices;
2682 std::vector<bool> exclude;
2684 unsigned int block = 0;
2687 if (pcell->is_active())
2690 for (
unsigned int child = 0; child < pcell->n_children(); ++child)
2692 const auto cell = pcell->child(child);
2697 indices.resize(n_dofs);
2698 exclude.resize(n_dofs);
2699 std::fill(exclude.begin(), exclude.end(),
false);
2700 cell->get_mg_dof_indices(indices);
2702 if (interior_dofs_only)
2706 for (
unsigned int d = 0; d < dim; ++d)
2708 const unsigned int face =
2717 for (
const unsigned int face :
2719 if (cell->at_boundary(face))
2725 for (
unsigned int i = 0; i < n_dofs; ++i)
2727 block_list.
add(block, indices[i]);
2733 template <
int dim,
int spacedim>
2734 std::vector<unsigned int>
2737 const unsigned int level,
2738 const bool interior_only,
2739 const bool boundary_patches,
2740 const bool level_boundary_patches,
2741 const bool single_cell_patches,
2742 const bool invert_vertex_mapping)
2749 exclude_boundary_dofs,
2751 level_boundary_patches,
2752 single_cell_patches,
2753 invert_vertex_mapping);
2756 template <
int dim,
int spacedim>
2757 std::vector<unsigned int>
2760 const unsigned int level,
2762 const bool boundary_patches,
2763 const bool level_boundary_patches,
2764 const bool single_cell_patches,
2765 const bool invert_vertex_mapping)
2769 std::vector<unsigned int> vertex_cell_count(
2773 std::vector<bool> vertex_boundary(
2776 std::vector<unsigned int> vertex_mapping(
2781 std::vector<unsigned int> vertex_dof_count(
2787 for (
const unsigned int v : cell->vertex_indices())
2789 const unsigned int vg = cell->vertex_index(v);
2790 vertex_dof_count[vg] += cell->get_fe().n_dofs_per_cell();
2791 ++vertex_cell_count[vg];
2792 for (
unsigned int d = 0; d < dim; ++d)
2795 if (cell->at_boundary(face))
2796 vertex_boundary[vg] =
true;
2797 else if ((!level_boundary_patches) &&
2798 (cell->neighbor(face)->level() !=
2799 static_cast<int>(
level)))
2800 vertex_boundary[vg] =
true;
2806 for (
unsigned int vg = 0; vg < vertex_dof_count.size(); ++vg)
2807 if ((!single_cell_patches && vertex_cell_count[vg] < 2) ||
2808 (!boundary_patches && vertex_boundary[vg]))
2809 vertex_dof_count[vg] = 0;
2812 unsigned int n_vertex_count = 0;
2813 for (
unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2814 if (vertex_dof_count[vg] != 0)
2815 vertex_mapping[vg] = n_vertex_count++;
2818 for (
unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2819 if (vertex_dof_count[vg] != 0)
2820 vertex_dof_count[vertex_mapping[vg]] = vertex_dof_count[vg];
2824 vertex_dof_count.resize(n_vertex_count);
2828 block_list.
reinit(vertex_dof_count.size(),
2832 std::vector<types::global_dof_index> indices;
2833 std::vector<bool> exclude;
2839 cell->get_mg_dof_indices(indices);
2841 for (
const unsigned int v : cell->vertex_indices())
2843 const unsigned int vg = cell->vertex_index(v);
2844 const unsigned int block = vertex_mapping[vg];
2850 if (exclude_boundary_dofs.
size() == 0 ||
2856 std::fill(exclude.begin(), exclude.end(),
false);
2858 for (
unsigned int d = 0; d < dim; ++d)
2860 const unsigned int a_face =
2862 const unsigned int face =
2875 for (
unsigned int j = 0; j < indices.size(); ++j)
2877 block_list.
add(block, indices[j]);
2881 for (
const auto index : indices)
2882 block_list.
add(block, index);
2887 if (invert_vertex_mapping)
2890 unsigned int n_vertex_count = 0;
2891 for (
unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2893 vertex_mapping[n_vertex_count++] = vg;
2896 vertex_mapping.resize(n_vertex_count);
2899 return vertex_mapping;
2903 template <
int dim,
int spacedim>
2909 std::set<types::global_dof_index> dofs_on_patch;
2910 std::vector<types::global_dof_index> local_dof_indices;
2915 for (
unsigned int i = 0; i < patch.size(); ++i)
2919 Assert(cell->is_artificial() ==
false,
2920 ExcMessage(
"This function can not be called with cells that are "
2921 "not either locally owned or ghost cells."));
2923 cell->get_dof_indices(local_dof_indices);
2924 dofs_on_patch.insert(local_dof_indices.begin(),
2925 local_dof_indices.end());
2929 return dofs_on_patch.size();
2934 template <
int dim,
int spacedim>
2935 std::vector<types::global_dof_index>
2940 std::set<types::global_dof_index> dofs_on_patch;
2941 std::vector<types::global_dof_index> local_dof_indices;
2946 for (
unsigned int i = 0; i < patch.size(); ++i)
2950 Assert(cell->is_artificial() ==
false,
2951 ExcMessage(
"This function can not be called with cells that are "
2952 "not either locally owned or ghost cells."));
2954 cell->get_dof_indices(local_dof_indices);
2955 dofs_on_patch.insert(local_dof_indices.begin(),
2956 local_dof_indices.end());
2966 return std::vector<types::global_dof_index>(dofs_on_patch.begin(),
2967 dofs_on_patch.end());
2977#include "dof_tools.inst"
const std::vector< std::pair< size_type, number > > * get_constraint_entries(const size_type line_n) const
size_type n_constraints() const
unsigned int size() const
unsigned int n_selected_blocks(const unsigned int overall_number_of_blocks=numbers::invalid_unsigned_int) const
bool represents_n_components(const unsigned int n) const
bool represents_the_all_selected_mask() const
unsigned int size() const
unsigned int n_selected_components(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
unsigned int first_selected_component(const unsigned int overall_number_of_components=numbers::invalid_unsigned_int) const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const IndexSet & locally_owned_mg_dofs(const unsigned int level) const
types::global_dof_index n_boundary_dofs() const
const Triangulation< dim, spacedim > & get_triangulation() const
const IndexSet & locally_owned_dofs() const
active_cell_iterator begin_active(const unsigned int level=0) const
std::vector< types::fe_index > get_active_fe_indices() const
bool has_hp_capabilities() const
types::global_dof_index n_dofs() const
types::global_dof_index n_locally_owned_dofs() const
const std::vector< Point< spacedim > > & get_quadrature_points() const
const FEValues< dim, spacedim > & get_present_fe_values() const
unsigned int n_dofs_per_vertex() const
unsigned int n_dofs_per_cell() const
unsigned int n_dofs_per_line() const
unsigned int n_dofs_per_face(unsigned int face_no=0, unsigned int child=0) const
unsigned int n_blocks() const
unsigned int n_components() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
const ComponentMask & get_nonzero_components(const unsigned int i) const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const unsigned char combined_orientation=ReferenceCell::default_combined_face_orientation()) const
bool is_primitive() const
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
unsigned int component_to_block_index(const unsigned int component) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
unsigned int element_multiplicity(const unsigned int index) const
unsigned int n_base_elements() const
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index, const unsigned int face_no=0) const
size_type index_within_set(const size_type global_index) const
bool is_element(const size_type index) const
void add_index(const size_type index)
void fill_binary_vector(VectorType &vector) const
void subtract_set(const IndexSet &other)
size_type nth_index_in_set(const size_type local_index) const
void add_indices(const ForwardIterator &begin, const ForwardIterator &end)
Abstract base class for mapping classes.
void reinit(const size_type m, const size_type n, const ArrayView< const unsigned int > &row_lengths)
void add(const size_type i, const size_type j)
virtual types::subdomain_id locally_owned_subdomain() const
unsigned int n_active_cells() const
unsigned int n_levels() const
unsigned int n_vertices() const
virtual size_type size() const override
unsigned int size() const
unsigned int max_dofs_per_face() const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
unsigned int n_components() const
unsigned int max_dofs_per_cell() const
void push_back(const Quadrature< dim_in > &new_quadrature)
virtual MPI_Comm get_communicator() const override
#define DEAL_II_DEPRECATED
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
IteratorRange< FilteredIterator< BaseIterator > > filter_iterators(IteratorRange< BaseIterator > i, const Predicate &p)
IteratorRange< cell_iterator > cell_iterators_on_level(const unsigned int level) const
static ::ExceptionBase & ExcIO()
static ::ExceptionBase & ExcInvalidBoundaryIndicator()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
#define AssertThrowMPI(error_code)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcIndexRange(std::size_t arg1, std::size_t arg2, std::size_t arg3)
static ::ExceptionBase & ExcDimensionMismatch(std::size_t arg1, std::size_t arg2)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
typename ActiveSelector::line_iterator line_iterator
typename ActiveSelector::face_iterator face_iterator
typename ActiveSelector::active_cell_iterator active_cell_iterator
typename ActiveSelector::CellAccessor cell_accessor
@ update_quadrature_points
Transformed quadrature points.
unsigned int n_mpi_processes(const MPI_Comm mpi_communicator)
const types::boundary_id internal_face_boundary_id
const types::subdomain_id artificial_subdomain_id
const types::subdomain_id invalid_subdomain_id
static const unsigned int invalid_unsigned_int
const types::global_dof_index invalid_dof_index
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
unsigned short int fe_index
static std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices()
#define DEAL_II_DOF_INDEX_MPI_TYPE