76 template <
int dim,
typename Number =
double>
88 double downstream_size = 0;
90 for (
unsigned int d = 0;
d < dim; ++
d)
92 downstream_size += (rhs[
d] - lhs[
d]) * weight;
95 if (downstream_size < 0)
97 else if (downstream_size > 0)
101 for (
unsigned int d = 0;
d < dim; ++
d)
103 if (lhs[
d] == rhs[
d])
105 return lhs[
d] < rhs[
d];
125 template <
int dim,
int spacedim>
126 std::vector<unsigned char>
130 std::vector<unsigned char> local_component_association(
139 local_component_association[i] =
150 const unsigned int first_comp =
155 local_component_association[i] = first_comp;
162 for (
unsigned int c = first_comp; c < fe.
n_components(); ++c)
163 if (component_mask[c] ==
true)
165 local_component_association[i] = c;
170 Assert(std::find(local_component_association.begin(),
171 local_component_association.end(),
172 static_cast<unsigned char>(-1)) ==
173 local_component_association.end(),
176 return local_component_association;
196 template <
typename DoFHandlerType>
200 std::vector<unsigned char> &dofs_by_component)
202 const ::hp::FECollection<DoFHandlerType::dimension,
203 DoFHandlerType::space_dimension>
204 &fe_collection = dof.get_fe_collection();
206 Assert(dofs_by_component.size() == dof.n_locally_owned_dofs(),
208 dof.n_locally_owned_dofs()));
217 std::vector<std::vector<unsigned char>> local_component_association(
218 fe_collection.size());
219 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
222 DoFHandlerType::space_dimension> &fe =
224 local_component_association[f] =
229 std::vector<types::global_dof_index> indices;
230 for (
const auto &c : dof.active_cell_iterators())
231 if (c->is_locally_owned())
233 const unsigned int fe_index = c->active_fe_index();
234 const unsigned int dofs_per_cell = c->get_fe().dofs_per_cell;
235 indices.resize(dofs_per_cell);
236 c->get_dof_indices(indices);
237 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
238 if (dof.locally_owned_dofs().is_element(indices[i]))
239 dofs_by_component[dof.locally_owned_dofs().index_within_set(
240 indices[i])] = local_component_association[fe_index][i];
252 template <
typename DoFHandlerType>
255 std::vector<unsigned char> &dofs_by_block)
257 const ::hp::FECollection<DoFHandlerType::dimension,
258 DoFHandlerType::space_dimension>
259 &fe_collection = dof.get_fe_collection();
261 Assert(dofs_by_block.size() == dof.n_locally_owned_dofs(),
263 dof.n_locally_owned_dofs()));
272 std::vector<std::vector<unsigned char>> local_block_association(
273 fe_collection.size());
274 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
277 DoFHandlerType::space_dimension> &fe =
279 local_block_association[f].resize(fe.dofs_per_cell,
280 static_cast<unsigned char>(-1));
281 for (
unsigned int i = 0; i < fe.dofs_per_cell; ++i)
282 local_block_association[f][i] = fe.system_to_block_index(i).first;
284 Assert(std::find(local_block_association[f].
begin(),
285 local_block_association[f].
end(),
286 static_cast<unsigned char>(-1)) ==
287 local_block_association[f].
end(),
292 std::vector<types::global_dof_index> indices;
293 for (
const auto &cell : dof.active_cell_iterators())
294 if (cell->is_locally_owned())
296 const unsigned int fe_index = cell->active_fe_index();
297 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
298 indices.resize(dofs_per_cell);
299 cell->get_dof_indices(indices);
300 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
301 if (dof.locally_owned_dofs().is_element(indices[i]))
302 dofs_by_block[dof.locally_owned_dofs().index_within_set(
303 indices[i])] = local_block_association[fe_index][i];
310 template <
typename DoFHandlerType,
typename Number>
313 const Vector<Number> &cell_data,
315 const unsigned int component)
317 const unsigned int dim = DoFHandlerType::dimension;
318 const unsigned int spacedim = DoFHandlerType::space_dimension;
331 const bool consider_components = (
n_components(dof_handler) != 1);
334 if (consider_components ==
false)
338 std::vector<unsigned char> component_dofs(
339 dof_handler.n_locally_owned_dofs());
342 dof_handler.get_fe_collection().component_mask(
346 for (
unsigned int i = 0; i < dof_data.
size(); ++i)
347 if (component_dofs[i] ==
static_cast<unsigned char>(component))
352 std::vector<unsigned char> touch_count(dof_handler.n_dofs(), 0);
354 typename DoFHandlerType::active_cell_iterator cell =
355 dof_handler.begin_active(),
356 endc = dof_handler.end();
357 std::vector<types::global_dof_index> dof_indices;
360 for (
unsigned int present_cell = 0; cell != endc; ++cell, ++present_cell)
362 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
363 dof_indices.resize(dofs_per_cell);
364 cell->get_dof_indices(dof_indices);
366 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
369 if (!consider_components ||
370 (cell->get_fe().system_to_component_index(i).first == component))
373 dof_data(dof_indices[i]) += cell_data(present_cell);
376 ++touch_count[dof_indices[i]];
386 Assert(consider_components || (touch_count[i] != 0),
388 if (touch_count[i] != 0)
389 dof_data(i) /= touch_count[i];
395 template <
typename DoFHandlerType>
399 std::vector<bool> & selected_dofs)
402 dof.get_fe_collection().n_components()),
404 "The given component mask is not sized correctly to represent the "
405 "components of the given finite element."));
406 Assert(selected_dofs.size() == dof.n_locally_owned_dofs(),
408 dof.n_locally_owned_dofs()));
413 dof.get_fe_collection().n_components()) == 0)
415 std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(),
false);
419 dof.get_fe_collection().n_components()) ==
420 dof.get_fe_collection().n_components())
422 std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(),
true);
428 std::fill_n(selected_dofs.begin(), dof.n_locally_owned_dofs(),
false);
432 std::vector<unsigned char> dofs_by_component(dof.n_locally_owned_dofs());
436 if (component_mask[dofs_by_component[i]] ==
true)
437 selected_dofs[i] =
true;
442 template <
typename DoFHandlerType>
447 dof.get_fe_collection().n_components()),
449 "The given component mask is not sized correctly to represent the "
450 "components of the given finite element."));
455 dof.get_fe_collection().n_components()) == 0)
458 dof.get_fe_collection().n_components()) ==
459 dof.get_fe_collection().n_components())
460 return dof.locally_owned_dofs();
466 const auto & fe_collection = dof.get_fe_collection();
467 std::vector<std::vector<bool>> local_selected_dofs;
469 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
472 const std::vector<unsigned char> local_component_association =
477 std::vector<bool> this_selected_dofs(fe_collection[f].dofs_per_cell);
478 for (
unsigned int i = 0; i < fe_collection[f].dofs_per_cell; ++i)
479 this_selected_dofs[i] =
480 component_mask[local_component_association[i]];
482 local_selected_dofs.emplace_back(std::move(this_selected_dofs));
486 IndexSet selected_dofs(dof.n_dofs());
487 std::vector<types::global_dof_index> local_dof_indices;
488 for (
const auto &c : dof.active_cell_iterators())
489 if (c->is_locally_owned())
491 const unsigned int fe_index = c->active_fe_index();
492 const unsigned int dofs_per_cell = c->get_fe().dofs_per_cell;
493 local_dof_indices.resize(dofs_per_cell);
494 c->get_dof_indices(local_dof_indices);
495 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
496 if (dof.locally_owned_dofs().is_element(local_dof_indices[i]) &&
497 local_selected_dofs[fe_index][i])
498 selected_dofs.
add_index(local_dof_indices[i]);
501 return selected_dofs;
506 template <
typename DoFHandlerType>
510 std::vector<bool> & selected_dofs)
513 extract_dofs<DoFHandlerType>(
514 dof, dof.get_fe_collection().component_mask(block_mask), selected_dofs);
519 template <
typename DoFHandlerType>
524 return extract_dofs<DoFHandlerType>(
525 dof, dof.get_fe_collection().component_mask(block_mask));
530 template <
typename DoFHandlerType>
531 std::vector<IndexSet>
535 const auto n_comps = dof.get_fe_collection().n_components();
538 "The given component mask is not sized correctly to represent the "
539 "components of the given finite element."));
541 const auto &locally_owned_dofs = dof.locally_owned_dofs();
545 std::vector<unsigned char> dofs_by_component(dof.n_locally_owned_dofs());
548 std::vector<IndexSet> index_per_comp(n_comps,
IndexSet(dof.n_dofs()));
552 const auto &comp_i = dofs_by_component[i];
553 if (component_mask[comp_i])
554 index_per_comp[comp_i].add_index(
555 locally_owned_dofs.nth_index_in_set(i));
557 for (
auto &c : index_per_comp)
559 return index_per_comp;
564 template <
typename DoFHandlerType>
567 const DoFHandlerType &dof,
569 std::vector<bool> & selected_dofs)
572 DoFHandlerType::space_dimension> &fe = dof.get_fe();
576 "The given component mask is not sized correctly to represent the "
577 "components of the given finite element."));
585 std::fill_n(selected_dofs.begin(), dof.n_dofs(
level),
false);
591 std::fill_n(selected_dofs.begin(), dof.n_dofs(
level),
true);
596 std::fill_n(selected_dofs.begin(), dof.n_dofs(
level),
false);
600 std::vector<unsigned char> local_component_asssociation =
602 std::vector<bool> local_selected_dofs(fe.dofs_per_cell);
603 for (
unsigned int i = 0; i < fe.dofs_per_cell; ++i)
604 local_selected_dofs[i] = component_mask[local_component_asssociation[i]];
607 std::vector<types::global_dof_index> indices(fe.dofs_per_cell);
608 typename DoFHandlerType::level_cell_iterator c;
609 for (c = dof.begin(
level); c != dof.end(
level); ++c)
611 c->get_mg_dof_indices(indices);
612 for (
unsigned int i = 0; i < fe.dofs_per_cell; ++i)
613 selected_dofs[indices[i]] = local_selected_dofs[i];
619 template <
typename DoFHandlerType>
622 const DoFHandlerType &dof,
624 std::vector<bool> & selected_dofs)
629 dof.get_fe().component_mask(block_mask),
635 template <
typename DoFHandlerType>
639 std::vector<bool> & selected_dofs,
640 const std::set<types::boundary_id> &boundary_ids)
643 DoFHandlerType::dimension,
644 DoFHandlerType::space_dimension
> *>(
645 &dof_handler.get_triangulation()) ==
nullptr),
647 "This function can not be used with distributed triangulations. "
648 "See the documentation for more information."));
654 selected_dofs.clear();
655 selected_dofs.resize(dof_handler.n_dofs(),
false);
662 template <
typename DoFHandlerType>
667 const std::set<types::boundary_id> &boundary_ids)
670 ExcMessage(
"Component mask has invalid size."));
674 const unsigned int dim = DoFHandlerType::dimension;
677 selected_dofs.
clear();
678 selected_dofs.
set_size(dof_handler.n_dofs());
682 const bool check_boundary_id = (boundary_ids.size() != 0);
686 const bool check_vector_component =
691 std::vector<types::global_dof_index> face_dof_indices;
700 for (
const auto &cell : dof_handler.active_cell_iterators())
704 if (cell->is_artificial() ==
false)
705 for (
const unsigned int face :
707 if (cell->at_boundary(face))
708 if (!check_boundary_id ||
709 (boundary_ids.find(cell->face(face)->boundary_id()) !=
713 DoFHandlerType::space_dimension> &fe =
717 face_dof_indices.resize(dofs_per_face);
718 cell->face(face)->get_dof_indices(face_dof_indices,
719 cell->active_fe_index());
721 for (
unsigned int i = 0; i < fe.dofs_per_face; ++i)
722 if (!check_vector_component)
723 selected_dofs.
add_index(face_dof_indices[i]);
731 const unsigned int cell_index =
735 (i < 2 * fe.dofs_per_vertex ?
737 i + 2 * fe.dofs_per_vertex) :
738 (dim == 3 ? (i < 4 * fe.dofs_per_vertex ?
740 (i < 4 * fe.dofs_per_vertex +
741 4 * fe.dofs_per_line ?
742 i + 4 * fe.dofs_per_vertex :
743 i + 4 * fe.dofs_per_vertex +
744 8 * fe.dofs_per_line)) :
746 if (fe.is_primitive(cell_index))
749 [fe.face_system_to_component_index(i).first] ==
751 selected_dofs.
add_index(face_dof_indices[i]);
755 const unsigned int first_nonzero_comp =
756 fe.get_nonzero_components(cell_index)
757 .first_selected_component();
758 Assert(first_nonzero_comp < fe.n_components(),
761 if (component_mask[first_nonzero_comp] ==
true)
762 selected_dofs.
add_index(face_dof_indices[i]);
770 template <
typename DoFHandlerType>
773 const DoFHandlerType & dof_handler,
775 std::vector<bool> & selected_dofs,
776 const std::set<types::boundary_id> &boundary_ids)
779 ExcMessage(
"This component mask has the wrong size."));
786 const bool check_boundary_id = (boundary_ids.size() != 0);
790 const bool check_vector_component =
794 selected_dofs.clear();
795 selected_dofs.resize(dof_handler.n_dofs(),
false);
796 std::vector<types::global_dof_index> cell_dof_indices;
805 for (
const auto &cell : dof_handler.active_cell_iterators())
806 for (
const unsigned int face :
808 if (cell->at_boundary(face))
809 if (!check_boundary_id ||
810 (boundary_ids.find(cell->face(face)->boundary_id()) !=
814 DoFHandlerType::space_dimension> &fe =
818 cell_dof_indices.resize(dofs_per_cell);
819 cell->get_dof_indices(cell_dof_indices);
821 for (
unsigned int i = 0; i < fe.dofs_per_cell; ++i)
822 if (fe.has_support_on_face(i, face))
824 if (!check_vector_component)
825 selected_dofs[cell_dof_indices[i]] =
true;
831 if (fe.is_primitive(i))
832 selected_dofs[cell_dof_indices[i]] =
833 (component_mask[fe.system_to_component_index(i)
837 const unsigned int first_nonzero_comp =
838 fe.get_nonzero_components(i)
839 .first_selected_component();
840 Assert(first_nonzero_comp < fe.n_components(),
843 selected_dofs[cell_dof_indices[i]] =
844 (component_mask[first_nonzero_comp] ==
true);
853 template <
typename DoFHandlerType,
typename number>
856 const DoFHandlerType &dof_handler,
858 bool(
const typename DoFHandlerType::active_cell_iterator &)> &predicate,
861 const std::function<
bool(
862 const typename DoFHandlerType::active_cell_iterator &)>
864 [=](
const typename DoFHandlerType::active_cell_iterator &cell) ->
bool {
865 return cell->is_locally_owned() && predicate(cell);
868 std::vector<types::global_dof_index> local_dof_indices;
872 std::set<types::global_dof_index> predicate_dofs;
874 typename DoFHandlerType::active_cell_iterator cell =
875 dof_handler.begin_active(),
876 endc = dof_handler.end();
877 for (; cell != endc; ++cell)
878 if (!cell->is_artificial() && predicate(cell))
880 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
881 cell->get_dof_indices(local_dof_indices);
882 predicate_dofs.insert(local_dof_indices.begin(),
883 local_dof_indices.end());
887 std::set<types::global_dof_index> dofs_with_support_on_halo_cells;
889 const std::vector<typename DoFHandlerType::active_cell_iterator>
892 for (
typename std::vector<
893 typename DoFHandlerType::active_cell_iterator>::const_iterator it =
895 it != halo_cells.end();
907 const unsigned int dofs_per_cell = (*it)->get_fe().dofs_per_cell;
908 local_dof_indices.resize(dofs_per_cell);
909 (*it)->get_dof_indices(local_dof_indices);
910 dofs_with_support_on_halo_cells.insert(local_dof_indices.begin(),
911 local_dof_indices.end());
923 std::set<types::global_dof_index> extra_halo;
924 for (
const auto dof : dofs_with_support_on_halo_cells)
930 const unsigned int line_size = line_ptr->size();
931 for (
unsigned int j = 0; j < line_size; ++j)
932 extra_halo.insert((*line_ptr)[j].first);
935 dofs_with_support_on_halo_cells.insert(extra_halo.begin(),
941 IndexSet support_set(dof_handler.n_dofs());
942 support_set.
add_indices(predicate_dofs.begin(), predicate_dofs.end());
945 IndexSet halo_set(dof_handler.n_dofs());
946 halo_set.
add_indices(dofs_with_support_on_halo_cells.begin(),
947 dofs_with_support_on_halo_cells.end());
962 template <
int spacedim>
965 const ::DoFHandler<1, spacedim> &dof_handler)
968 return IndexSet(dof_handler.n_dofs());
972 template <
int spacedim>
975 const ::DoFHandler<2, spacedim> &dof_handler)
977 const unsigned int dim = 2;
979 IndexSet selected_dofs(dof_handler.n_dofs());
985 for (
const auto &cell : dof_handler.active_cell_iterators())
986 if (!cell->is_artificial())
989 if (cell->face(face)->has_children())
991 const typename ::DoFHandler<dim,
992 spacedim>::line_iterator
993 line = cell->face(face);
996 selected_dofs.add_index(
997 line->child(0)->vertex_dof_index(1, dof));
999 for (
unsigned int child = 0; child < 2; ++child)
1001 if (cell->neighbor_child_on_subface(face, child)
1006 selected_dofs.add_index(
1007 line->child(child)->dof_index(dof));
1012 selected_dofs.compress();
1013 return selected_dofs;
1017 template <
int spacedim>
1020 const ::DoFHandler<3, spacedim> &dof_handler)
1022 const unsigned int dim = 3;
1024 IndexSet selected_dofs(dof_handler.n_dofs());
1025 IndexSet unconstrained_dofs(dof_handler.n_dofs());
1029 for (
const auto &cell : dof_handler.active_cell_iterators())
1030 if (!cell->is_artificial())
1033 const typename ::DoFHandler<dim, spacedim>::face_iterator
1034 face = cell->face(f);
1035 if (cell->face(f)->has_children())
1037 for (
unsigned int child = 0; child < 4; ++child)
1038 if (!cell->neighbor_child_on_subface(f, child)
1042 std::vector<types::global_dof_index> ldi(
1044 face->child(child)->get_dof_indices(ldi);
1045 selected_dofs.add_indices(ldi.begin(), ldi.end());
1050 for (
unsigned int vertex = 0; vertex < 4; ++vertex)
1053 unconstrained_dofs.add_index(
1054 face->vertex_dof_index(vertex, dof));
1057 selected_dofs.subtract_set(unconstrained_dofs);
1058 return selected_dofs;
1065 template <
int dim,
int spacedim>
1068 std::vector<bool> & selected_dofs)
1070 const IndexSet selected_dofs_as_index_set =
1075 std::fill(selected_dofs.begin(), selected_dofs.end(),
false);
1076 for (
const auto index : selected_dofs_as_index_set)
1077 selected_dofs[index] =
true;
1082 template <
int dim,
int spacedim>
1091 template <
typename DoFHandlerType>
1095 std::vector<bool> & selected_dofs)
1097 Assert(selected_dofs.size() == dof_handler.n_dofs(),
1101 std::fill_n(selected_dofs.begin(), dof_handler.n_dofs(),
false);
1103 std::vector<types::global_dof_index> local_dof_indices;
1108 typename DoFHandlerType::active_cell_iterator cell =
1109 dof_handler.begin_active(),
1110 endc = dof_handler.end();
1111 for (; cell != endc; ++cell)
1114 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1115 local_dof_indices.resize(dofs_per_cell);
1116 cell->get_dof_indices(local_dof_indices);
1117 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1118 selected_dofs[local_dof_indices[i]] =
true;
1124 template <
typename DoFHandlerType>
1130 dof_set = dof_handler.locally_owned_dofs();
1136 template <
typename DoFHandlerType>
1142 dof_set = dof_handler.locally_owned_dofs();
1147 std::vector<types::global_dof_index> dof_indices;
1148 std::set<types::global_dof_index> global_dof_indices;
1150 typename DoFHandlerType::active_cell_iterator cell =
1151 dof_handler.begin_active(),
1152 endc = dof_handler.
end();
1153 for (; cell != endc; ++cell)
1154 if (cell->is_locally_owned())
1156 dof_indices.resize(cell->get_fe().dofs_per_cell);
1157 cell->get_dof_indices(dof_indices);
1161 global_dof_indices.insert(dof_index);
1164 dof_set.
add_indices(global_dof_indices.begin(), global_dof_indices.end());
1171 template <
typename DoFHandlerType>
1177 dof_set = dof_handler.locally_owned_dofs();
1187 std::vector<types::global_dof_index> dof_indices;
1188 std::vector<types::global_dof_index> dofs_on_ghosts;
1190 typename DoFHandlerType::active_cell_iterator cell =
1191 dof_handler.begin_active(),
1192 endc = dof_handler.
end();
1193 for (; cell != endc; ++cell)
1194 if (cell->is_ghost())
1196 dof_indices.resize(cell->get_fe().dofs_per_cell);
1197 cell->get_dof_indices(dof_indices);
1198 for (
const auto dof_index : dof_indices)
1200 dofs_on_ghosts.push_back(dof_index);
1204 std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1206 std::unique(dofs_on_ghosts.begin(),
1207 dofs_on_ghosts.end()));
1213 template <
typename DoFHandlerType>
1216 const unsigned int level,
1220 dof_set = dof_handler.locally_owned_mg_dofs(
level);
1230 std::vector<types::global_dof_index> dof_indices;
1231 std::vector<types::global_dof_index> dofs_on_ghosts;
1233 typename DoFHandlerType::cell_iterator cell = dof_handler.
begin(
level),
1234 endc = dof_handler.end(
level);
1235 for (; cell != endc; ++cell)
1240 if (
id == dof_handler.get_triangulation().locally_owned_subdomain() ||
1244 dof_indices.resize(cell->get_fe().dofs_per_cell);
1245 cell->get_mg_dof_indices(dof_indices);
1246 for (
const auto dof_index : dof_indices)
1248 dofs_on_ghosts.push_back(dof_index);
1252 std::sort(dofs_on_ghosts.begin(), dofs_on_ghosts.end());
1254 std::unique(dofs_on_ghosts.begin(),
1255 dofs_on_ghosts.end()));
1262 template <
typename DoFHandlerType>
1266 std::vector<std::vector<bool>> &constant_modes)
1270 if (dof_handler.n_locally_owned_dofs() == 0)
1272 constant_modes = std::vector<std::vector<bool>>(0);
1276 const unsigned int n_components = dof_handler.get_fe(0).n_components();
1280 std::vector<unsigned char> dofs_by_component(
1281 dof_handler.n_locally_owned_dofs());
1285 unsigned int n_selected_dofs = 0;
1287 if (component_mask[i] ==
true)
1289 std::count(dofs_by_component.begin(), dofs_by_component.end(), i);
1292 const IndexSet &locally_owned_dofs = dof_handler.locally_owned_dofs();
1293 std::vector<unsigned int> component_numbering(
1295 for (
unsigned int i = 0, count = 0; i < locally_owned_dofs.
n_elements();
1297 if (component_mask[dofs_by_component[i]])
1298 component_numbering[i] = count++;
1305 const ::hp::FECollection<DoFHandlerType::dimension,
1306 DoFHandlerType::space_dimension>
1307 & fe_collection = dof_handler.get_fe_collection();
1308 std::vector<Table<2, bool>> element_constant_modes;
1309 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
1311 unsigned int n_constant_modes = 0;
1312 for (
unsigned int f = 0; f < fe_collection.size(); ++f)
1314 std::pair<Table<2, bool>, std::vector<unsigned int>> data =
1315 fe_collection[f].get_constant_modes();
1316 element_constant_modes.push_back(data.first);
1318 for (
unsigned int i = 0; i < data.second.size(); ++i)
1319 if (component_mask[data.second[i]])
1320 constant_mode_to_component_translation[data.second[i]]
1321 .emplace_back(n_constant_modes++, i);
1323 element_constant_modes[0].n_rows());
1327 constant_modes.clear();
1328 constant_modes.resize(n_constant_modes,
1329 std::vector<bool>(n_selected_dofs,
false));
1332 std::vector<types::global_dof_index> dof_indices;
1333 for (
const auto &cell : dof_handler.active_cell_iterators())
1334 if (cell->is_locally_owned())
1336 dof_indices.resize(cell->get_fe().dofs_per_cell);
1337 cell->get_dof_indices(dof_indices);
1339 for (
unsigned int i = 0; i < dof_indices.size(); ++i)
1340 if (locally_owned_dofs.
is_element(dof_indices[i]))
1342 const unsigned int loc_index =
1344 const unsigned int comp = dofs_by_component[loc_index];
1345 if (component_mask[comp])
1346 for (
auto &indices :
1347 constant_mode_to_component_translation[comp])
1348 constant_modes[indices
1349 .first][component_numbering[loc_index]] =
1350 element_constant_modes[cell->active_fe_index()](
1358 template <
typename DoFHandlerType>
1361 std::vector<unsigned int> &active_fe_indices)
1364 dof_handler.get_triangulation().n_active_cells());
1366 typename DoFHandlerType::active_cell_iterator cell =
1367 dof_handler.begin_active(),
1368 endc = dof_handler.end();
1369 for (; cell != endc; ++cell)
1370 active_fe_indices[cell->active_cell_index()] = cell->active_fe_index();
1373 template <
typename DoFHandlerType>
1374 std::vector<IndexSet>
1377 Assert(dof_handler.n_dofs() > 0,
1378 ExcMessage(
"The given DoFHandler has no DoFs."));
1383 DoFHandlerType::dimension,
1384 DoFHandlerType::space_dimension
> *>(
1385 &dof_handler.get_triangulation()) ==
nullptr),
1387 "For parallel::distributed::Triangulation objects and "
1388 "associated DoF handler objects, asking for any information "
1389 "related to a subdomain other than the locally owned one does "
1390 "not make sense."));
1394 std::vector<::types::subdomain_id> subdomain_association(
1395 dof_handler.n_dofs());
1397 subdomain_association);
1411 const unsigned int n_subdomains =
1413 DoFHandlerType::dimension,
1414 DoFHandlerType::space_dimension
> *>(
1415 &dof_handler.get_triangulation()) ==
nullptr ?
1417 unsigned int max_subdomain_id = 0;
1418 for (
const auto &cell : dof_handler.active_cell_iterators())
1420 std::max(max_subdomain_id, cell->subdomain_id());
1421 return max_subdomain_id + 1;
1425 DoFHandlerType::dimension,
1426 DoFHandlerType::space_dimension
> *>(
1427 &dof_handler.get_triangulation())
1428 ->get_communicator()));
1429 Assert(n_subdomains > *std::max_element(subdomain_association.begin(),
1430 subdomain_association.end()),
1433 std::vector<::IndexSet> index_sets(
1434 n_subdomains, ::
IndexSet(dof_handler.n_dofs()));
1442 index < subdomain_association.size();
1446 if (subdomain_association[index] != this_subdomain)
1448 index_sets[this_subdomain].add_range(i_min, index);
1450 this_subdomain = subdomain_association[index];
1455 if (i_min == subdomain_association.size() - 1)
1457 index_sets[this_subdomain].add_index(i_min);
1463 index_sets[this_subdomain].add_range(i_min,
1464 subdomain_association.size());
1467 for (
unsigned int i = 0; i < n_subdomains; i++)
1473 template <
typename DoFHandlerType>
1474 std::vector<IndexSet>
1480 DoFHandlerType::dimension,
1481 DoFHandlerType::space_dimension
> *>(
1482 &dof_handler.get_triangulation()) ==
nullptr),
1484 "For parallel::distributed::Triangulation objects and "
1485 "associated DoF handler objects, asking for any information "
1486 "related to a subdomain other than the locally owned one does "
1487 "not make sense."));
1495 std::vector<IndexSet> dof_set =
1508 const typename DoFHandlerType::active_cell_iterator &)>
1510 const std::vector<typename DoFHandlerType::active_cell_iterator>
1515 std::vector<types::global_dof_index> local_dof_indices;
1516 std::set<types::global_dof_index> subdomain_halo_global_dof_indices;
1517 for (
typename std::vector<
1518 typename DoFHandlerType::active_cell_iterator>::const_iterator
1519 it_cell = active_halo_layer.begin();
1520 it_cell != active_halo_layer.end();
1523 const typename DoFHandlerType::active_cell_iterator &cell =
1528 "The subdomain ID of the halo cell should not match that of the vector entry."));
1530 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
1531 cell->get_dof_indices(local_dof_indices);
1535 subdomain_halo_global_dof_indices.insert(local_dof_index);
1539 subdomain_halo_global_dof_indices.begin(),
1540 subdomain_halo_global_dof_indices.end());
1548 template <
typename DoFHandlerType>
1551 const DoFHandlerType & dof_handler,
1552 std::vector<types::subdomain_id> &subdomain_association)
1557 DoFHandlerType::dimension,
1558 DoFHandlerType::space_dimension
> *>(
1559 &dof_handler.get_triangulation()) ==
nullptr),
1561 "For parallel::distributed::Triangulation objects and "
1562 "associated DoF handler objects, asking for any subdomain other "
1563 "than the locally owned one does not make sense."));
1565 Assert(subdomain_association.size() == dof_handler.n_dofs(),
1567 dof_handler.n_dofs()));
1579 std::vector<types::subdomain_id> cell_owners(
1580 dof_handler.get_triangulation().n_active_cells());
1582 DoFHandlerType::space_dimension>
1584 DoFHandlerType::dimension,
1585 DoFHandlerType::space_dimension
> *>(
1586 &dof_handler.get_triangulation())))
1588 cell_owners = tr->get_true_subdomain_ids_of_cells();
1589 Assert(tr->get_true_subdomain_ids_of_cells().size() ==
1590 tr->n_active_cells(),
1595 for (
typename DoFHandlerType::active_cell_iterator cell =
1596 dof_handler.begin_active();
1597 cell != dof_handler.end();
1599 if (cell->is_locally_owned())
1600 cell_owners[cell->active_cell_index()] = cell->subdomain_id();
1604 std::fill_n(subdomain_association.begin(),
1605 dof_handler.n_dofs(),
1608 std::vector<types::global_dof_index> local_dof_indices;
1613 typename DoFHandlerType::active_cell_iterator cell =
1614 dof_handler.begin_active(),
1615 endc = dof_handler.end();
1616 for (; cell != endc; ++cell)
1619 cell_owners[cell->active_cell_index()];
1620 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1621 local_dof_indices.resize(dofs_per_cell);
1622 cell->get_dof_indices(local_dof_indices);
1628 for (
unsigned int i = 0; i < dofs_per_cell; ++i)
1629 if (subdomain_association[local_dof_indices[i]] ==
1631 subdomain_association[local_dof_indices[i]] =
subdomain_id;
1632 else if (subdomain_association[local_dof_indices[i]] >
subdomain_id)
1634 subdomain_association[local_dof_indices[i]] =
subdomain_id;
1638 Assert(std::find(subdomain_association.begin(),
1639 subdomain_association.end(),
1641 subdomain_association.end(),
1647 template <
typename DoFHandlerType>
1652 std::vector<types::subdomain_id> subdomain_association(
1653 dof_handler.n_dofs());
1656 return std::count(subdomain_association.begin(),
1657 subdomain_association.end(),
1663 template <
typename DoFHandlerType>
1670 Assert((dof_handler.get_triangulation().locally_owned_subdomain() ==
1673 dof_handler.get_triangulation().locally_owned_subdomain()),
1675 "For parallel::distributed::Triangulation objects and "
1676 "associated DoF handler objects, asking for any subdomain other "
1677 "than the locally owned one does not make sense."));
1679 IndexSet index_set(dof_handler.n_dofs());
1681 std::vector<types::global_dof_index> local_dof_indices;
1689 std::vector<types::global_dof_index> subdomain_indices;
1691 typename DoFHandlerType::active_cell_iterator cell =
1692 dof_handler.begin_active(),
1693 endc = dof_handler.end();
1694 for (; cell != endc; ++cell)
1695 if ((cell->is_artificial() ==
false) &&
1696 (cell->subdomain_id() == subdomain))
1698 const unsigned int dofs_per_cell = cell->get_fe().dofs_per_cell;
1699 local_dof_indices.resize(dofs_per_cell);
1700 cell->get_dof_indices(local_dof_indices);
1701 subdomain_indices.insert(subdomain_indices.end(),
1702 local_dof_indices.begin(),
1703 local_dof_indices.end());
1706 std::sort(subdomain_indices.begin(), subdomain_indices.end());
1707 subdomain_indices.erase(std::unique(subdomain_indices.begin(),
1708 subdomain_indices.end()),
1709 subdomain_indices.end());
1712 index_set.
add_indices(subdomain_indices.begin(), subdomain_indices.end());
1720 template <
typename DoFHandlerType>
1723 const DoFHandlerType & dof_handler,
1725 std::vector<unsigned int> &n_dofs_on_subdomain)
1727 Assert(n_dofs_on_subdomain.size() == dof_handler.get_fe(0).n_components(),
1729 dof_handler.get_fe(0).n_components()));
1730 std::fill(n_dofs_on_subdomain.begin(), n_dofs_on_subdomain.end(), 0);
1734 dof_handler.begin_active(),
1735 typename DoFHandlerType::active_cell_iterator{dof_handler.end()},
1736 [subdomain](
const typename DoFHandlerType::cell_accessor &cell) {
1737 return cell.subdomain_id() == subdomain;
1739 ExcMessage(
"There are no cells for the given subdomain!"));
1741 std::vector<types::subdomain_id> subdomain_association(
1742 dof_handler.n_dofs());
1745 std::vector<unsigned char> component_association(dof_handler.n_dofs());
1747 std::vector<bool>(),
1748 component_association);
1750 for (
unsigned int c = 0; c < dof_handler.get_fe(0).
n_components(); ++c)
1753 if ((subdomain_association[i] == subdomain) &&
1754 (component_association[i] ==
static_cast<unsigned char>(c)))
1755 ++n_dofs_on_subdomain[c];
1767 template <
int dim,
int spacedim>
1770 const std::vector<unsigned char> & dofs_by_component,
1771 const std::vector<unsigned int> & target_component,
1772 const bool only_once,
1773 std::vector<types::global_dof_index> &dofs_per_component,
1774 unsigned int & component)
1793 for (
unsigned int dd = 0; dd <
d; ++dd, ++component)
1794 dofs_per_component[target_component[component]] +=
1795 std::count(dofs_by_component.begin(),
1796 dofs_by_component.end(),
1803 for (
unsigned int dd = 1; dd <
d; ++dd)
1804 dofs_per_component[target_component[component -
d + dd]] =
1805 dofs_per_component[target_component[component -
d]];
1812 template <
int dim,
int spacedim>
1815 const std::vector<unsigned char> & dofs_by_component,
1816 const std::vector<unsigned int> & target_component,
1817 const bool only_once,
1818 std::vector<types::global_dof_index> &dofs_per_component,
1819 unsigned int & component)
1824 for (
unsigned int fe = 1; fe < fe_collection.
size(); ++fe)
1829 Assert(fe_collection[fe].n_base_elements() ==
1830 fe_collection[0].n_base_elements(),
1832 for (
unsigned int b = 0;
b < fe_collection[0].n_base_elements(); ++
b)
1837 Assert(fe_collection[fe].base_element(
b).n_base_elements() ==
1838 fe_collection[0].base_element(
b).n_base_elements(),
1861 template <
int dim,
int spacedim>
1873 template <
int dim,
int spacedim>
1875 all_elements_are_primitive(
1876 const ::hp::FECollection<dim, spacedim> &fe_collection)
1878 for (
unsigned int i = 0; i < fe_collection.size(); ++i)
1879 if (fe_collection[i].is_primitive() ==
false)
1890 template <
typename DoFHandlerType>
1893 const DoFHandlerType & dof_handler,
1894 std::vector<types::global_dof_index> &dofs_per_component,
1895 const bool only_once,
1896 const std::vector<unsigned int> & target_component)
1898 dofs_per_component =
1904 template <
typename DoFHandlerType>
1905 std::vector<types::global_dof_index>
1907 const DoFHandlerType & dof_handler,
1908 const bool only_once,
1909 const std::vector<unsigned int> &target_component_)
1911 const unsigned int n_components = dof_handler.get_fe(0).n_components();
1915 std::vector<unsigned int> target_component = target_component_;
1916 if (target_component.size() == 0)
1920 target_component[i] = i;
1927 const unsigned int max_component =
1928 *std::max_element(target_component.begin(), target_component.end());
1929 const unsigned int n_target_components = max_component + 1;
1931 std::vector<types::global_dof_index> dofs_per_component(
1938 dofs_per_component[0] = dof_handler.n_locally_owned_dofs();
1939 return dofs_per_component;
1945 std::vector<unsigned char> dofs_by_component(
1946 dof_handler.n_locally_owned_dofs());
1952 unsigned int component = 0;
1963 Assert((internal::all_elements_are_primitive(
1964 dof_handler.get_fe_collection()) ==
false) ||
1965 (std::accumulate(dofs_per_component.begin(),
1966 dofs_per_component.end(),
1968 dof_handler.n_locally_owned_dofs()),
1972 #ifdef DEAL_II_WITH_MPI
1973 const unsigned int dim = DoFHandlerType::dimension;
1974 const unsigned int spacedim = DoFHandlerType::space_dimension;
1978 &dof_handler.get_triangulation())))
1980 std::vector<types::global_dof_index> local_dof_count =
1983 const int ierr = MPI_Allreduce(local_dof_count.data(),
1984 dofs_per_component.data(),
1985 n_target_components,
1988 tria->get_communicator());
1993 return dofs_per_component;
1999 template <
typename DoFHandlerType>
2002 std::vector<types::global_dof_index> &dofs_per_block,
2003 const std::vector<unsigned int> & target_block)
2010 template <
typename DoFHandlerType>
2011 std::vector<types::global_dof_index>
2013 const std::vector<unsigned int> &target_block_)
2015 const ::hp::FECollection<DoFHandlerType::dimension,
2016 DoFHandlerType::space_dimension>
2017 &fe_collection = dof_handler.get_fe_collection();
2025 const unsigned int n_blocks = fe_collection[0].n_blocks();
2027 std::vector<unsigned int> target_block = target_block_;
2028 if (target_block.size() == 0)
2030 target_block.resize(fe_collection[0].
n_blocks());
2031 for (
unsigned int i = 0; i <
n_blocks; ++i)
2032 target_block[i] = i;
2037 for (
unsigned int f = 1; f < fe_collection.size(); ++f)
2039 ExcMessage(
"This function can only work if all elements in a "
2040 "collection have the same number of blocks."));
2046 std::vector<types::global_dof_index> dofs_per_block(1);
2047 dofs_per_block[0] = dof_handler.n_dofs();
2048 return dofs_per_block;
2052 const unsigned int max_block =
2053 *std::max_element(target_block.begin(), target_block.end());
2054 const unsigned int n_target_blocks = max_block + 1;
2056 std::vector<types::global_dof_index> dofs_per_block(n_target_blocks);
2060 for (
unsigned int this_fe = fe_collection.size() - 1;
2061 this_fe < fe_collection.size();
2065 DoFHandlerType::space_dimension> &fe =
2066 fe_collection[this_fe];
2068 std::vector<unsigned char> dofs_by_block(
2069 dof_handler.n_locally_owned_dofs());
2073 for (
unsigned int block = 0; block < fe.n_blocks(); ++block)
2074 dofs_per_block[target_block[block]] +=
2075 std::count(dofs_by_block.begin(), dofs_by_block.end(), block);
2077 #ifdef DEAL_II_WITH_MPI
2081 DoFHandlerType::space_dimension>
2083 DoFHandlerType::dimension,
2084 DoFHandlerType::space_dimension
> *>(
2085 &dof_handler.get_triangulation())))
2087 std::vector<types::global_dof_index> local_dof_count =
2089 const int ierr = MPI_Allreduce(local_dof_count.data(),
2090 dofs_per_block.data(),
2094 tria->get_communicator());
2100 return dofs_per_block;
2105 template <
typename DoFHandlerType>
2108 std::vector<types::global_dof_index> &mapping)
2111 mapping.insert(mapping.end(),
2112 dof_handler.n_dofs(),
2115 std::vector<types::global_dof_index> dofs_on_face;
2125 typename DoFHandlerType::active_cell_iterator cell =
2126 dof_handler.begin_active(),
2127 endc = dof_handler.end();
2128 for (; cell != endc; ++cell)
2129 for (
const unsigned int f :
2131 if (cell->at_boundary(f))
2133 const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
2134 dofs_on_face.resize(dofs_per_face);
2135 cell->face(f)->get_dof_indices(dofs_on_face,
2136 cell->active_fe_index());
2137 for (
unsigned int i = 0; i < dofs_per_face; ++i)
2139 mapping[dofs_on_face[i]] = next_boundary_index++;
2147 template <
typename DoFHandlerType>
2150 const std::set<types::boundary_id> &boundary_ids,
2151 std::vector<types::global_dof_index> &mapping)
2158 mapping.insert(mapping.end(),
2159 dof_handler.n_dofs(),
2163 if (boundary_ids.size() == 0)
2166 std::vector<types::global_dof_index> dofs_on_face;
2170 typename DoFHandlerType::active_cell_iterator cell =
2171 dof_handler.begin_active(),
2172 endc = dof_handler.end();
2173 for (; cell != endc; ++cell)
2174 for (
const unsigned int f :
2176 if (boundary_ids.find(cell->face(f)->boundary_id()) !=
2179 const unsigned int dofs_per_face = cell->get_fe().dofs_per_face;
2180 dofs_on_face.resize(dofs_per_face);
2181 cell->face(f)->get_dof_indices(dofs_on_face,
2182 cell->active_fe_index());
2183 for (
unsigned int i = 0; i < dofs_per_face; ++i)
2185 mapping[dofs_on_face[i]] = next_boundary_index++;
2189 dof_handler.n_boundary_dofs(boundary_ids));
2196 template <
typename DoFHandlerType>
2200 DoFHandlerType::space_dimension> &mapping,
2201 const DoFHandlerType & dof_handler,
2206 const unsigned int dim = DoFHandlerType::dimension;
2207 const unsigned int spacedim = DoFHandlerType::space_dimension;
2210 dof_handler.get_fe_collection();
2213 for (
unsigned int fe_index = 0; fe_index < fe_collection.
size();
2217 Assert(fe_collection[fe_index].has_support_points(),
2220 fe_collection[fe_index].get_unit_support_points()));
2225 (in_mask.
size() == 0 ?
2240 typename DoFHandlerType::active_cell_iterator cell = dof_handler
2242 endc = dof_handler.end();
2244 std::vector<types::global_dof_index> local_dof_indices;
2245 for (; cell != endc; ++cell)
2247 if (cell->is_artificial() ==
false)
2249 hp_fe_values.reinit(cell);
2253 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
2254 cell->get_dof_indices(local_dof_indices);
2256 const std::vector<Point<spacedim>> &points =
2257 fe_values.get_quadrature_points();
2258 for (
unsigned int i = 0; i < cell->get_fe().dofs_per_cell; ++i)
2260 const unsigned int dof_comp =
2261 cell->get_fe().system_to_component_index(i).first;
2265 support_points[local_dof_indices[i]] = points[i];
2271 template <
typename DoFHandlerType>
2275 DoFHandlerType::space_dimension> &mapping,
2276 const DoFHandlerType & dof_handler,
2293 Assert(x_support_points.find(i) != x_support_points.end(),
2296 support_points[i] = x_support_points[i];
2302 template <
int dim,
int spacedim>
2314 "This function can not be used with distributed triangulations. "
2315 "See the documentation for more information."));
2328 template <
int dim,
int spacedim>
2341 "This function can not be used with distributed triangulations. "
2342 "See the documentation for more information."));
2353 template <
int dim,
int spacedim>
2361 support_points.clear();
2374 template <
int dim,
int spacedim>
2382 support_points.clear();
2392 template <
int spacedim>
2400 using dof_map_t = std::map<types::global_dof_index, Point<spacedim>>;
2402 using point_map_t = std::map<Point<spacedim>,
2403 std::vector<types::global_dof_index>,
2406 point_map_t point_map;
2409 for (
typename dof_map_t::const_iterator it = support_points.begin();
2410 it != support_points.end();
2413 std::vector<types::global_dof_index> &v = point_map[it->second];
2414 v.push_back(it->first);
2418 for (
typename point_map_t::iterator it = point_map.begin();
2419 it != point_map.end();
2422 out << it->first <<
" \"";
2423 const std::vector<types::global_dof_index> &v = it->second;
2424 for (
unsigned int i = 0; i < v.size(); ++i)
2438 template <
int dim,
int spacedim>
2445 const unsigned int nb = fe.
n_blocks();
2447 tables_by_block.resize(1);
2448 tables_by_block[0].reinit(nb, nb);
2449 tables_by_block[0].fill(
none);
2457 tables_by_block[0](ib, jb) |= table(i, j);
2463 template <
int dim,
int spacedim>
2471 tables_by_block.resize(fe_collection.
size());
2473 for (
unsigned int f = 0; f < fe_collection.
size(); ++f)
2477 const unsigned int nb = fe.
n_blocks();
2478 tables_by_block[f].reinit(nb, nb);
2479 tables_by_block[f].fill(
none);
2486 tables_by_block[f](ib, jb) |= table(i, j);
2494 template <
int dim,
int spacedim>
2498 const unsigned int level,
2499 const std::vector<bool> & selected_dofs,
2505 std::vector<types::global_dof_index> indices;
2509 for (cell = dof_handler.
begin(
level); cell != endc; ++cell)
2510 if (cell->is_locally_owned_on_level())
2514 dof_handler.
get_fe().dofs_per_cell);
2516 for (cell = dof_handler.
begin(
level); cell != endc; ++cell)
2517 if (cell->is_locally_owned_on_level())
2519 indices.resize(cell->
get_fe().dofs_per_cell);
2520 cell->get_mg_dof_indices(indices);
2522 if (selected_dofs.size() != 0)
2527 if (selected_dofs.size() == 0)
2528 block_list.
add(i, indices[j] - offset);
2531 if (selected_dofs[j])
2532 block_list.
add(i, indices[j] - offset);
2540 template <
typename DoFHandlerType>
2543 const DoFHandlerType &dof_handler,
2544 const unsigned int level,
2545 const bool interior_only)
2548 block_list.
reinit(1, dof_handler.n_dofs(
level), dof_handler.n_dofs(
level));
2549 typename DoFHandlerType::level_cell_iterator cell;
2550 typename DoFHandlerType::level_cell_iterator endc = dof_handler.end(
level);
2552 std::vector<types::global_dof_index> indices;
2553 std::vector<bool> exclude;
2555 for (cell = dof_handler.begin(
level); cell != endc; ++cell)
2557 indices.resize(cell->get_fe().dofs_per_cell);
2558 cell->get_mg_dof_indices(indices);
2564 std::fill(exclude.begin(), exclude.end(),
false);
2567 for (
const unsigned int face :
2569 if (cell->at_boundary(face) ||
2570 cell->neighbor(face)->level() != cell->level())
2571 for (
unsigned int i = 0; i < dpf; ++i)
2575 block_list.
add(0, indices[j]);
2579 for (
const auto index : indices)
2580 block_list.
add(0, index);
2586 template <
typename DoFHandlerType>
2589 const DoFHandlerType &dof_handler,
2590 const unsigned int level,
2591 const bool interior_dofs_only,
2592 const bool boundary_dofs)
2597 typename DoFHandlerType::level_cell_iterator pcell =
2598 dof_handler.begin(
level - 1);
2599 typename DoFHandlerType::level_cell_iterator endc =
2600 dof_handler.end(
level - 1);
2602 std::vector<types::global_dof_index> indices;
2603 std::vector<bool> exclude;
2605 for (
unsigned int block = 0; pcell != endc; ++pcell)
2607 if (pcell->is_active())
2610 for (
unsigned int child = 0; child < pcell->n_children(); ++child)
2612 const typename DoFHandlerType::level_cell_iterator cell =
2613 pcell->child(child);
2617 dof_handler.get_fe();
2619 indices.resize(n_dofs);
2620 exclude.resize(n_dofs);
2621 std::fill(exclude.begin(), exclude.end(),
false);
2622 cell->get_mg_dof_indices(indices);
2624 if (interior_dofs_only)
2630 for (
unsigned int d = 0;
d < DoFHandlerType::dimension; ++
d)
2633 DoFHandlerType::dimension>::vertex_to_face[child][
d];
2634 for (
unsigned int i = 0; i < dpf; ++i)
2641 for (
const unsigned int face :
2643 if (cell->at_boundary(face))
2644 for (
unsigned int i = 0; i < dpf; ++i)
2648 for (
unsigned int i = 0; i < n_dofs; ++i)
2650 block_list.
add(block, indices[i]);
2656 template <
typename DoFHandlerType>
2657 std::vector<unsigned int>
2659 const DoFHandlerType &dof_handler,
2660 const unsigned int level,
2661 const bool interior_only,
2662 const bool boundary_patches,
2663 const bool level_boundary_patches,
2664 const bool single_cell_patches,
2665 const bool invert_vertex_mapping)
2667 const unsigned int n_blocks = dof_handler.get_fe().n_blocks();
2672 exclude_boundary_dofs,
2674 level_boundary_patches,
2675 single_cell_patches,
2676 invert_vertex_mapping);
2679 template <
typename DoFHandlerType>
2680 std::vector<unsigned int>
2682 const DoFHandlerType &dof_handler,
2683 const unsigned int level,
2684 const BlockMask & exclude_boundary_dofs,
2685 const bool boundary_patches,
2686 const bool level_boundary_patches,
2687 const bool single_cell_patches,
2688 const bool invert_vertex_mapping)
2690 typename DoFHandlerType::level_cell_iterator cell;
2691 typename DoFHandlerType::level_cell_iterator endc = dof_handler.end(
level);
2695 std::vector<unsigned int> vertex_cell_count(
2696 dof_handler.get_triangulation().n_vertices(), 0);
2699 std::vector<bool> vertex_boundary(
2700 dof_handler.get_triangulation().n_vertices(),
false);
2702 std::vector<unsigned int> vertex_mapping(
2703 dof_handler.get_triangulation().n_vertices(),
2707 std::vector<unsigned int> vertex_dof_count(
2708 dof_handler.get_triangulation().n_vertices(), 0);
2712 for (cell = dof_handler.begin(
level); cell != endc; ++cell)
2713 for (
const unsigned int v :
2716 const unsigned int vg = cell->vertex_index(v);
2717 vertex_dof_count[vg] += cell->get_fe().dofs_per_cell;
2718 ++vertex_cell_count[vg];
2719 for (
unsigned int d = 0;
d < DoFHandlerType::dimension; ++
d)
2721 const unsigned int face =
2723 if (cell->at_boundary(face))
2724 vertex_boundary[vg] =
true;
2725 else if ((!level_boundary_patches) &&
2726 (cell->neighbor(face)->level() !=
2727 static_cast<int>(
level)))
2728 vertex_boundary[vg] =
true;
2734 for (
unsigned int vg = 0; vg < vertex_dof_count.size(); ++vg)
2735 if ((!single_cell_patches && vertex_cell_count[vg] < 2) ||
2736 (!boundary_patches && vertex_boundary[vg]))
2737 vertex_dof_count[vg] = 0;
2740 unsigned int n_vertex_count = 0;
2741 for (
unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2742 if (vertex_dof_count[vg] != 0)
2743 vertex_mapping[vg] = n_vertex_count++;
2746 for (
unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2747 if (vertex_dof_count[vg] != 0)
2748 vertex_dof_count[vertex_mapping[vg]] = vertex_dof_count[vg];
2752 vertex_dof_count.resize(n_vertex_count);
2756 block_list.
reinit(vertex_dof_count.size(),
2757 dof_handler.n_dofs(
level),
2760 std::vector<types::global_dof_index> indices;
2761 std::vector<bool> exclude;
2763 for (cell = dof_handler.begin(
level); cell != endc; ++cell)
2767 cell->get_mg_dof_indices(indices);
2769 for (
const unsigned int v :
2772 const unsigned int vg = cell->vertex_index(v);
2773 const unsigned int block = vertex_mapping[vg];
2779 if (exclude_boundary_dofs.
size() == 0 ||
2785 std::fill(exclude.begin(), exclude.end(),
false);
2788 for (
unsigned int d = 0;
d < DoFHandlerType::dimension; ++
d)
2791 DoFHandlerType::dimension>::vertex_to_face[v][
d];
2793 DoFHandlerType::dimension>::opposite_face[a_face];
2794 for (
unsigned int i = 0; i < dpf; ++i)
2805 for (
unsigned int j = 0; j < indices.size(); ++j)
2807 block_list.
add(block, indices[j]);
2811 for (
const auto index : indices)
2812 block_list.
add(block, index);
2817 if (invert_vertex_mapping)
2820 unsigned int n_vertex_count = 0;
2821 for (
unsigned int vg = 0; vg < vertex_mapping.size(); ++vg)
2823 vertex_mapping[n_vertex_count++] = vg;
2826 vertex_mapping.resize(n_vertex_count);
2829 return vertex_mapping;
2833 template <
typename DoFHandlerType>
2836 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch)
2838 std::set<types::global_dof_index> dofs_on_patch;
2839 std::vector<types::global_dof_index> local_dof_indices;
2844 for (
unsigned int i = 0; i < patch.size(); ++i)
2846 const typename DoFHandlerType::active_cell_iterator cell = patch[i];
2847 Assert(cell->is_artificial() ==
false,
2848 ExcMessage(
"This function can not be called with cells that are "
2849 "not either locally owned or ghost cells."));
2850 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
2851 cell->get_dof_indices(local_dof_indices);
2852 dofs_on_patch.insert(local_dof_indices.begin(),
2853 local_dof_indices.end());
2857 return dofs_on_patch.size();
2862 template <
typename DoFHandlerType>
2863 std::vector<types::global_dof_index>
2865 const std::vector<typename DoFHandlerType::active_cell_iterator> &patch)
2867 std::set<types::global_dof_index> dofs_on_patch;
2868 std::vector<types::global_dof_index> local_dof_indices;
2873 for (
unsigned int i = 0; i < patch.size(); ++i)
2875 const typename DoFHandlerType::active_cell_iterator cell = patch[i];
2876 Assert(cell->is_artificial() ==
false,
2877 ExcMessage(
"This function can not be called with cells that are "
2878 "not either locally owned or ghost cells."));
2879 local_dof_indices.resize(cell->get_fe().dofs_per_cell);
2880 cell->get_dof_indices(local_dof_indices);
2881 dofs_on_patch.insert(local_dof_indices.begin(),
2882 local_dof_indices.end());
2885 Assert(dofs_on_patch.size() == count_dofs_on_patch<DoFHandlerType>(patch),
2892 return std::vector<types::global_dof_index>(dofs_on_patch.begin(),
2893 dofs_on_patch.end());
2903 #include "dof_tools.inst"