54 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
58 const
Point<spacedim> & p,
59 const
std::vector<
bool> &marked_vertices)
70 marked_vertices.size() == 0,
72 marked_vertices.size()));
82 marked_vertices.size() == 0 ||
83 std::equal(marked_vertices.begin(),
84 marked_vertices.end(),
86 [](
bool p,
bool q) { return !p || q; }),
88 "marked_vertices should be a subset of used vertices in the triangulation "
89 "but marked_vertices contains one or more vertices that are not used vertices!"));
97 const std::vector<bool> &used = (marked_vertices.size() == 0) ?
103 std::vector<bool>::const_iterator
first =
104 std::find(used.begin(), used.end(),
true);
110 unsigned int best_vertex = std::distance(used.begin(),
first);
111 double best_dist = (p -
vertices[best_vertex]).norm_square();
115 for (
unsigned int j = best_vertex + 1; j <
vertices.size(); ++j)
118 double dist = (p -
vertices[j]).norm_square();
119 if (dist < best_dist)
131 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
135 const MeshType<dim, spacedim> &mesh,
136 const
Point<spacedim> & p,
137 const
std::vector<
bool> &marked_vertices)
140 if (mapping.preserves_vertex_locations() ==
true)
152 marked_vertices.size() == 0,
154 marked_vertices.size()));
164 marked_vertices.size() == 0 ||
165 std::equal(marked_vertices.begin(),
166 marked_vertices.end(),
168 [](
bool p,
bool q) { return !p || q; }),
170 "marked_vertices should be a subset of used vertices in the triangulation "
171 "but marked_vertices contains one or more vertices that are not used vertices!"));
174 if (marked_vertices.size())
177 if (marked_vertices[it->first] ==
false)
192 template <
class MeshType>
195 std::vector<typename MeshType::active_cell_iterator>
198 typename ::internal::ActiveCellIterator<MeshType::dimension,
199 MeshType::space_dimension,
203 const unsigned int vertex)
205 const int dim = MeshType::dimension;
206 const int spacedim = MeshType::space_dimension;
212 Assert(mesh.get_triangulation().get_used_vertices()[vertex],
218 std::set<typename ::internal::
219 ActiveCellIterator<dim, spacedim, MeshType>::type>
222 typename ::internal::ActiveCellIterator<dim, spacedim, MeshType>::type
223 cell = mesh.begin_active(),
259 for (; cell != endc; ++cell)
261 for (
const unsigned int v : cell->vertex_indices())
262 if (cell->vertex_index(v) == vertex)
274 for (
const auto face :
275 cell->reference_cell().faces_for_given_vertex(v))
276 if (!cell->at_boundary(face) &&
277 cell->neighbor(face)->is_active())
299 for (
unsigned int e = 0; e < cell->n_lines(); ++e)
300 if (cell->line(e)->has_children())
304 if (cell->line(e)->child(0)->vertex_index(1) == vertex)
327 return std::vector<typename ::internal::
328 ActiveCellIterator<dim, spacedim, MeshType>::type>(
336 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
339 void find_active_cell_around_point_internal(
340 const MeshType<dim, spacedim> &mesh,
342 std::set<typename MeshType<dim, spacedim>::active_cell_iterator>
344 std::set<typename MeshType<dim, spacedim>::active_cell_iterator>
348 typename ::internal::
349 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
353 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
358 using cell_iterator =
359 typename MeshType<dim, spacedim>::active_cell_iterator;
361 using cell_iterator = typename ::internal::
362 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type;
370 std::set<cell_iterator> adjacent_cells_new;
374 std::vector<cell_iterator> active_neighbors;
375 get_active_neighbors<MeshType<dim, spacedim>>(cell, active_neighbors);
376 for (
unsigned int i = 0; i < active_neighbors.size(); ++i)
377 if (searched_cells.find(active_neighbors[i]) ==
378 searched_cells.end())
379 adjacent_cells_new.insert(active_neighbors[i]);
383 adjacent_cells_new.end());
392 cell_iterator it = mesh.begin_active();
393 for (; it != mesh.end(); ++it)
394 if (searched_cells.find(it) == searched_cells.end())
405 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
409 typename MeshType<dim, spacedim>::active_cell_iterator
411 typename ::internal::
412 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
416 const std::vector<bool> &marked_vertices,
417 const double tolerance)
419 return find_active_cell_around_point<dim, MeshType, spacedim>(
430 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
434 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
436 std::pair<typename ::internal::
437 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
441 const MeshType<dim, spacedim> &mesh,
443 const std::vector<bool> &marked_vertices,
444 const double tolerance)
446 using active_cell_iterator = typename ::internal::
447 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type;
453 double best_distance = tolerance;
455 std::pair<active_cell_iterator, Point<dim>> best_cell;
458 best_cell.first = mesh.end();
462 std::vector<active_cell_iterator> adjacent_cells_tmp =
471 std::set<active_cell_iterator>
adjacent_cells(adjacent_cells_tmp.begin(),
472 adjacent_cells_tmp.end());
473 std::set<active_cell_iterator> searched_cells;
481 const auto n_active_cells = mesh.get_triangulation().n_active_cells();
483 unsigned int cells_searched = 0;
484 while (!found && cells_searched < n_active_cells)
488 if (cell->is_artificial() ==
false)
491 if (marked_vertices.size() > 0)
493 bool any_vertex_marked =
false;
494 for (
const auto &v : cell->vertex_indices())
496 if (marked_vertices[cell->vertex_index(v)])
498 any_vertex_marked =
true;
502 if (!any_vertex_marked)
514 cell->reference_cell().closest_point(p_cell).distance(
521 if ((dist < best_distance) ||
522 ((dist == best_distance) &&
523 (cell->level() > best_level)))
526 best_distance = dist;
527 best_level = cell->level();
528 best_cell = std::make_pair(cell, p_cell);
559 if (!found && cells_searched < n_active_cells)
561 find_active_cell_around_point_internal<dim, MeshType, spacedim>(
571 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
575 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
578 std::vector<std::pair<
579 typename ::internal::
580 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
584 const MeshType<dim, spacedim> &mesh,
586 const double tolerance,
587 const std::vector<bool> &marked_vertices)
590 mapping, mesh, p, marked_vertices, tolerance);
592 if (cell_and_point.first == mesh.end())
596 mapping, mesh, p, tolerance, cell_and_point);
601 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
605 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
608 std::vector<std::pair<
609 typename ::internal::
610 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
615 const MeshType<dim, spacedim> &mesh,
617 const double tolerance,
618 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
622 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
627 cells_and_points.push_back(first_cell);
631 const Point<dim> unit_point = cells_and_points.front().second;
632 const auto my_cell = cells_and_points.front().first;
635 unsigned int n_dirs_at_threshold = 0;
637 for (
unsigned int d = 0; d < dim; ++d)
639 distance_to_center[d] =
std::abs(unit_point[d] - 0.5);
640 if (distance_to_center[d] > 0.5 - tolerance)
642 ++n_dirs_at_threshold;
643 last_point_at_threshold = d;
647 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
650 if (n_dirs_at_threshold == 1)
652 unsigned int neighbor_index =
653 2 * last_point_at_threshold +
654 (unit_point[last_point_at_threshold] > 0.5 ? 1 : 0);
655 if (!my_cell->at_boundary(neighbor_index))
657 const auto neighbor_cell = my_cell->neighbor(neighbor_index);
659 if (neighbor_cell->is_active())
660 cells_to_add.push_back(neighbor_cell);
662 for (
const auto &child_cell : neighbor_cell->child_iterators())
664 if (child_cell->is_active())
665 cells_to_add.push_back(child_cell);
670 else if (n_dirs_at_threshold == dim)
672 unsigned int local_vertex_index = 0;
673 for (
unsigned int d = 0; d < dim; ++d)
674 local_vertex_index += (unit_point[d] > 0.5 ? 1 : 0) << d;
675 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
677 mesh, my_cell->vertex_index(local_vertex_index));
678 for (
const auto &cell : cells)
681 cells_to_add.push_back(cell);
689 else if (n_dirs_at_threshold == 2)
692 unsigned int count_vertex_indices = 0;
694 for (
unsigned int d = 0; d < dim; ++d)
696 if (distance_to_center[d] > 0.5 - tolerance)
700 unit_point[d] > 0.5 ? 1 : 0;
701 count_vertex_indices++;
711 const unsigned int first_vertex =
714 for (
unsigned int d = 0; d < 2; ++d)
718 my_cell->vertex_index(first_vertex + (d << free_direction)));
719 for (
const auto &cell : tentative_cells)
721 bool cell_not_yet_present =
true;
722 for (
const auto &other_cell : cells_to_add)
723 if (cell == other_cell)
725 cell_not_yet_present =
false;
728 if (cell_not_yet_present)
729 cells_to_add.push_back(cell);
734 const double original_distance_to_unit_cell =
735 my_cell->reference_cell().closest_point(unit_point).distance(unit_point);
736 for (
const auto &cell : cells_to_add)
743 if (cell->reference_cell().closest_point(p_unit).distance(
744 p_unit) < original_distance_to_unit_cell + tolerance)
745 cells_and_points.emplace_back(cell, p_unit);
752 cells_and_points.begin(),
753 cells_and_points.end(),
754 [](
const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
756 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
757 Point<dim>> &b) { return a.first < b.first; });
759 return cells_and_points;
764 template <
class MeshType>
768 const MeshType &mesh,
769 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
772 std::vector<typename MeshType::active_cell_iterator> active_halo_layer;
773 std::vector<bool> locally_active_vertices_on_subdomain(
774 mesh.get_triangulation().n_vertices(),
false);
779 for (
const auto &cell : mesh.active_cell_iterators())
781 for (
const auto v : cell->vertex_indices())
782 locally_active_vertices_on_subdomain[cell->vertex_index(v)] =
true;
787 for (
const auto &cell : mesh.active_cell_iterators())
788 if (!predicate(cell))
789 for (
const auto v : cell->vertex_indices())
790 if (locally_active_vertices_on_subdomain[cell->vertex_index(v)] ==
793 active_halo_layer.push_back(cell);
797 return active_halo_layer;
802 template <
class MeshType>
806 const MeshType &mesh,
807 const std::function<
bool(
const typename MeshType::cell_iterator &)>
809 const unsigned int level)
811 std::vector<typename MeshType::cell_iterator> level_halo_layer;
812 std::vector<bool> locally_active_vertices_on_level_subdomain(
813 mesh.get_triangulation().n_vertices(),
false);
818 for (
typename MeshType::cell_iterator cell = mesh.begin(
level);
819 cell != mesh.end(
level);
822 for (
const unsigned int v : cell->vertex_indices())
823 locally_active_vertices_on_level_subdomain[cell->vertex_index(v)] =
829 for (
typename MeshType::cell_iterator cell = mesh.begin(
level);
830 cell != mesh.end(
level);
832 if (!predicate(cell))
833 for (
const unsigned int v : cell->vertex_indices())
834 if (locally_active_vertices_on_level_subdomain[cell->vertex_index(
837 level_halo_layer.push_back(cell);
841 return level_halo_layer;
847 template <
class MeshType>
849 bool contains_locally_owned_cells(
850 const std::vector<typename MeshType::active_cell_iterator> &cells)
852 for (
typename std::vector<
853 typename MeshType::active_cell_iterator>::const_iterator it =
858 if ((*it)->is_locally_owned())
864 template <
class MeshType>
866 bool contains_artificial_cells(
867 const std::vector<typename MeshType::active_cell_iterator> &cells)
869 for (
typename std::vector<
870 typename MeshType::active_cell_iterator>::const_iterator it =
875 if ((*it)->is_artificial())
884 template <
class MeshType>
890 std::function<
bool(
const typename MeshType::active_cell_iterator &)>
893 const std::vector<typename MeshType::active_cell_iterator>
898 Assert(contains_locally_owned_cells<MeshType>(active_halo_layer) ==
false,
899 ExcMessage(
"Halo layer contains locally owned cells"));
900 Assert(contains_artificial_cells<MeshType>(active_halo_layer) ==
false,
901 ExcMessage(
"Halo layer contains artificial cells"));
903 return active_halo_layer;
908 template <
class MeshType>
912 const MeshType &mesh,
913 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
915 const double layer_thickness)
917 std::vector<typename MeshType::active_cell_iterator>
918 subdomain_boundary_cells, active_cell_layer_within_distance;
919 std::vector<bool> vertices_outside_subdomain(
920 mesh.get_triangulation().n_vertices(),
false);
922 const unsigned int spacedim = MeshType::space_dimension;
924 unsigned int n_non_predicate_cells = 0;
932 for (
const auto &cell : mesh.active_cell_iterators())
933 if (!predicate(cell))
935 for (
const unsigned int v : cell->vertex_indices())
936 vertices_outside_subdomain[cell->vertex_index(v)] =
true;
937 n_non_predicate_cells++;
944 if (n_non_predicate_cells == 0 ||
945 n_non_predicate_cells == mesh.get_triangulation().n_active_cells())
946 return std::vector<typename MeshType::active_cell_iterator>();
950 for (
const auto &cell : mesh.active_cell_iterators())
953 for (
const unsigned int v : cell->vertex_indices())
954 if (vertices_outside_subdomain[cell->vertex_index(v)] ==
true)
956 subdomain_boundary_cells.push_back(cell);
966 const double DOUBLE_EPSILON = 100. * std::numeric_limits<double>::epsilon();
969 for (
unsigned int d = 0; d < spacedim; ++d)
971 bounding_box.first[d] -= (layer_thickness + DOUBLE_EPSILON);
972 bounding_box.second[d] += (layer_thickness + DOUBLE_EPSILON);
975 std::vector<Point<spacedim>>
976 subdomain_boundary_cells_centers;
979 subdomain_boundary_cells_radii;
981 subdomain_boundary_cells_centers.reserve(subdomain_boundary_cells.size());
982 subdomain_boundary_cells_radii.reserve(subdomain_boundary_cells.size());
984 for (
typename std::vector<typename MeshType::active_cell_iterator>::
985 const_iterator subdomain_boundary_cell_iterator =
986 subdomain_boundary_cells.begin();
987 subdomain_boundary_cell_iterator != subdomain_boundary_cells.end();
988 ++subdomain_boundary_cell_iterator)
990 const std::pair<Point<spacedim>,
double>
991 &subdomain_boundary_cell_enclosing_ball =
992 (*subdomain_boundary_cell_iterator)->enclosing_ball();
994 subdomain_boundary_cells_centers.push_back(
995 subdomain_boundary_cell_enclosing_ball.first);
996 subdomain_boundary_cells_radii.push_back(
997 subdomain_boundary_cell_enclosing_ball.second);
999 AssertThrow(subdomain_boundary_cells_radii.size() ==
1000 subdomain_boundary_cells_centers.size(),
1009 for (
const auto &cell : mesh.active_cell_iterators())
1012 if (predicate(cell))
1015 const std::pair<Point<spacedim>,
double> &cell_enclosing_ball =
1016 cell->enclosing_ball();
1019 cell_enclosing_ball.first;
1020 const double cell_enclosing_ball_radius = cell_enclosing_ball.second;
1022 bool cell_inside =
true;
1024 for (
unsigned int d = 0; d < spacedim; ++d)
1026 (cell_enclosing_ball_center[d] + cell_enclosing_ball_radius >
1027 bounding_box.first[d]) &&
1028 (cell_enclosing_ball_center[d] - cell_enclosing_ball_radius <
1029 bounding_box.second[d]);
1035 for (
unsigned int i = 0; i < subdomain_boundary_cells_radii.size();
1038 subdomain_boundary_cells_centers[i]) <
1039 Utilities::fixed_power<2>(cell_enclosing_ball_radius +
1040 subdomain_boundary_cells_radii[i] +
1041 layer_thickness + DOUBLE_EPSILON))
1043 active_cell_layer_within_distance.push_back(cell);
1048 return active_cell_layer_within_distance;
1053 template <
class MeshType>
1063 std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1064 predicate(locally_owned_cell_predicate);
1066 const std::vector<typename MeshType::active_cell_iterator>
1067 ghost_cell_layer_within_distance =
1075 contains_locally_owned_cells<MeshType>(
1076 ghost_cell_layer_within_distance) ==
false,
1078 "Ghost cells within layer_thickness contains locally owned cells."));
1080 contains_artificial_cells<MeshType>(ghost_cell_layer_within_distance) ==
1083 "Ghost cells within layer_thickness contains artificial cells. "
1084 "The function compute_ghost_cell_layer_within_distance "
1085 "is probably called while using parallel::distributed::Triangulation. "
1086 "In such case please refer to the description of this function."));
1088 return ghost_cell_layer_within_distance;
1093 template <
class MeshType>
1099 const std::function<
bool(
1100 const typename MeshType::
1101 active_cell_iterator &)>
1104 std::vector<bool> locally_active_vertices_on_subdomain(
1105 mesh.get_triangulation().n_vertices(),
false);
1107 const unsigned int spacedim = MeshType::space_dimension;
1114 for (
const auto &cell : mesh.active_cell_iterators())
1115 if (predicate(cell))
1117 minp = cell->center();
1118 maxp = cell->center();
1124 for (
const auto &cell : mesh.active_cell_iterators())
1125 if (predicate(cell))
1126 for (
const unsigned int v : cell->vertex_indices())
1127 if (locally_active_vertices_on_subdomain[cell->vertex_index(v)] ==
1130 locally_active_vertices_on_subdomain[cell->vertex_index(v)] =
1132 for (
unsigned int d = 0; d < spacedim; ++d)
1134 minp[d] =
std::min(minp[d], cell->vertex(v)[d]);
1135 maxp[d] =
std::max(maxp[d], cell->vertex(v)[d]);
1139 return std::make_pair(minp, maxp);
1144 template <
typename MeshType>
1146 std::list<std::pair<
1147 typename MeshType::cell_iterator,
1154 ExcMessage(
"The two meshes must be represent triangulations that "
1155 "have the same coarse meshes"));
1162 bool remove_ghost_cells =
false;
1163#ifdef DEAL_II_WITH_MPI
1165 constexpr int dim = MeshType::dimension;
1166 constexpr int spacedim = MeshType::space_dimension;
1168 *
>(&mesh_1.get_triangulation()) !=
nullptr ||
1170 *
>(&mesh_2.get_triangulation()) !=
nullptr)
1172 Assert(&mesh_1.get_triangulation() == &mesh_2.get_triangulation(),
1173 ExcMessage(
"This function can only be used with meshes "
1174 "corresponding to distributed Triangulations when "
1175 "both Triangulations are equal."));
1176 remove_ghost_cells =
true;
1188 using CellList = std::list<std::pair<
typename MeshType::cell_iterator,
1189 typename MeshType::cell_iterator>>;
1193 typename MeshType::cell_iterator cell_1 = mesh_1.begin(0),
1194 cell_2 = mesh_2.begin(0);
1195 for (; cell_1 != mesh_1.end(0); ++cell_1, ++cell_2)
1196 cell_list.emplace_back(cell_1, cell_2);
1199 typename CellList::iterator cell_pair = cell_list.begin();
1200 while (cell_pair != cell_list.end())
1204 if (cell_pair->first->has_children() &&
1205 cell_pair->second->has_children())
1207 Assert(cell_pair->first->refinement_case() ==
1208 cell_pair->second->refinement_case(),
1210 for (
unsigned int c = 0; c < cell_pair->first->n_children(); ++c)
1211 cell_list.emplace_back(cell_pair->first->child(c),
1212 cell_pair->second->child(c));
1217 const auto previous_cell_pair = cell_pair;
1219 cell_list.erase(previous_cell_pair);
1224 if (remove_ghost_cells &&
1225 ((cell_pair->first->is_active() &&
1226 !cell_pair->first->is_locally_owned()) ||
1227 (cell_pair->second->is_active() &&
1228 !cell_pair->second->is_locally_owned())))
1231 const auto previous_cell_pair = cell_pair;
1233 cell_list.erase(previous_cell_pair);
1242 for (cell_pair = cell_list.begin(); cell_pair != cell_list.end();
1244 Assert(cell_pair->first->is_active() || cell_pair->second->is_active() ||
1245 (cell_pair->first->refinement_case() !=
1246 cell_pair->second->refinement_case()),
1254 template <
int dim,
int spacedim>
1271 endc = mesh_1.
end(0);
1272 for (; cell_1 != endc; ++cell_1, ++cell_2)
1274 if (cell_1->n_vertices() != cell_2->n_vertices())
1276 for (
const unsigned int v : cell_1->vertex_indices())
1277 if (cell_1->vertex(v) != cell_2->vertex(v))
1290 template <
typename MeshType>
1295 mesh_2.get_triangulation());
1300 template <
int dim,
int spacedim>
1301 std::pair<typename DoFHandler<dim, spacedim>::active_cell_iterator,
1307 const double tolerance)
1311 ExcMessage(
"Mapping collection needs to have either size 1 "
1312 "or size equal to the number of elements in "
1313 "the FECollection."));
1315 using cell_iterator =
1318 std::pair<cell_iterator, Point<dim>> best_cell;
1322 if (mapping.
size() == 1)
1324 const std::vector<bool> marked_vertices = {};
1326 mapping[0], mesh, p, marked_vertices, tolerance);
1333 double best_distance = tolerance;
1334 int best_level = -1;
1341 std::vector<cell_iterator> adjacent_cells_tmp =
1349 std::set<cell_iterator>
adjacent_cells(adjacent_cells_tmp.begin(),
1350 adjacent_cells_tmp.end());
1351 std::set<cell_iterator> searched_cells;
1361 unsigned int cells_searched = 0;
1362 while (!found && cells_searched < n_cells)
1369 mapping[cell->active_fe_index()]
1370 .transform_real_to_unit_cell(cell, p);
1376 cell->reference_cell().closest_point(p_cell).distance(
1383 if (dist < best_distance ||
1384 (dist == best_distance && cell->level() > best_level))
1387 best_distance = dist;
1388 best_level = cell->level();
1389 best_cell = std::make_pair(cell, p_cell);
1415 if (!found && cells_searched < n_cells)
1417 find_active_cell_around_point_internal<dim,
1429 template <
class MeshType>
1432 const typename MeshType::active_cell_iterator &cell)
1434 Assert(cell->is_locally_owned(),
1435 ExcMessage(
"This function only makes sense if the cell for "
1436 "which you are asking for a patch, is locally "
1439 std::vector<typename MeshType::active_cell_iterator> patch;
1440 patch.push_back(cell);
1441 for (
const unsigned int face_number : cell->face_indices())
1442 if (cell->face(face_number)->at_boundary() ==
false)
1444 if (cell->neighbor(face_number)->has_children() ==
false)
1445 patch.push_back(cell->neighbor(face_number));
1450 if (MeshType::dimension > 1)
1452 for (
unsigned int subface = 0;
1453 subface < cell->face(face_number)->n_children();
1456 cell->neighbor_child_on_subface(face_number, subface));
1462 typename MeshType::cell_iterator neighbor =
1463 cell->neighbor(face_number);
1464 while (neighbor->has_children())
1465 neighbor = neighbor->child(1 - face_number);
1467 Assert(neighbor->neighbor(1 - face_number) == cell,
1469 patch.push_back(neighbor);
1477 template <
class Container>
1478 std::vector<typename Container::cell_iterator>
1480 const std::vector<typename Container::active_cell_iterator> &patch)
1484 "Vector containing patch cells should not be an empty vector!"));
1488 int min_level = patch[0]->level();
1490 for (
unsigned int i = 0; i < patch.size(); ++i)
1492 std::set<typename Container::cell_iterator> uniform_cells;
1493 typename std::vector<
1494 typename Container::active_cell_iterator>::const_iterator patch_cell;
1496 for (patch_cell = patch.begin(); patch_cell != patch.end(); ++patch_cell)
1501 if ((*patch_cell)->level() == min_level)
1502 uniform_cells.insert(*patch_cell);
1509 typename Container::cell_iterator parent = *patch_cell;
1511 while (parent->level() > min_level)
1512 parent = parent->parent();
1513 uniform_cells.insert(parent);
1517 return std::vector<typename Container::cell_iterator>(uniform_cells.begin(),
1518 uniform_cells.end());
1523 template <
class Container>
1526 const std::vector<typename Container::active_cell_iterator> &patch,
1528 &local_triangulation,
1531 Container::space_dimension>::active_cell_iterator,
1532 typename Container::active_cell_iterator> &patch_to_global_tria_map)
1535 const std::vector<typename Container::cell_iterator> uniform_cells =
1536 get_cells_at_coarsest_common_level<Container>(patch);
1538 local_triangulation.
clear();
1539 std::vector<Point<Container::space_dimension>>
vertices;
1540 const unsigned int n_uniform_cells = uniform_cells.size();
1541 std::vector<CellData<Container::dimension>> cells(n_uniform_cells);
1544 typename std::vector<typename Container::cell_iterator>::const_iterator
1546 for (uniform_cell = uniform_cells.begin();
1547 uniform_cell != uniform_cells.end();
1550 for (
const unsigned int v : (*uniform_cell)->vertex_indices())
1553 (*uniform_cell)->vertex(v);
1554 bool repeat_vertex =
false;
1556 for (
unsigned int m = 0; m < i; ++m)
1560 repeat_vertex =
true;
1561 cells[k].vertices[v] = m;
1565 if (repeat_vertex ==
false)
1568 cells[k].vertices[v] = i;
1579 unsigned int index = 0;
1582 Container::space_dimension>::cell_iterator,
1583 typename Container::cell_iterator>
1584 patch_to_global_tria_map_tmp;
1586 Container::space_dimension>::cell_iterator
1587 coarse_cell = local_triangulation.
begin();
1588 coarse_cell != local_triangulation.
end();
1589 ++coarse_cell, ++index)
1591 patch_to_global_tria_map_tmp.insert(
1592 std::make_pair(coarse_cell, uniform_cells[index]));
1596 Assert(coarse_cell->center().distance(uniform_cells[index]->center()) <=
1597 1e-15 * coarse_cell->diameter(),
1600 bool refinement_necessary;
1605 refinement_necessary =
false;
1606 for (
const auto &active_tria_cell :
1609 if (patch_to_global_tria_map_tmp[active_tria_cell]->has_children())
1611 active_tria_cell->set_refine_flag();
1612 refinement_necessary =
true;
1615 for (
unsigned int i = 0; i < patch.size(); ++i)
1620 if (patch_to_global_tria_map_tmp[active_tria_cell] ==
1625 for (
const unsigned int v :
1626 active_tria_cell->vertex_indices())
1627 active_tria_cell->vertex(v) = patch[i]->vertex(v);
1629 Assert(active_tria_cell->center().distance(
1630 patch_to_global_tria_map_tmp[active_tria_cell]
1632 1e-15 * active_tria_cell->diameter(),
1635 active_tria_cell->set_user_flag();
1641 if (refinement_necessary)
1646 Container::dimension,
1647 Container::space_dimension>::cell_iterator cell =
1648 local_triangulation.
begin();
1649 cell != local_triangulation.
end();
1652 if (patch_to_global_tria_map_tmp.find(cell) !=
1653 patch_to_global_tria_map_tmp.end())
1655 if (cell->has_children())
1662 for (
unsigned int c = 0; c < cell->n_children(); ++c)
1664 if (patch_to_global_tria_map_tmp.find(cell->child(
1665 c)) == patch_to_global_tria_map_tmp.end())
1667 patch_to_global_tria_map_tmp.insert(
1670 patch_to_global_tria_map_tmp[cell]->child(
1692 patch_to_global_tria_map_tmp.erase(cell);
1698 while (refinement_necessary);
1704 Container::space_dimension>::cell_iterator
1705 cell = local_triangulation.
begin();
1706 cell != local_triangulation.
end();
1709 if (cell->user_flag_set())
1711 Assert(patch_to_global_tria_map_tmp.find(cell) !=
1712 patch_to_global_tria_map_tmp.end(),
1715 Assert(cell->center().distance(
1716 patch_to_global_tria_map_tmp[cell]->center()) <=
1717 1e-15 * cell->diameter(),
1725 Container::space_dimension>::cell_iterator,
1726 typename Container::cell_iterator>::iterator
1727 map_tmp_it = patch_to_global_tria_map_tmp.
begin(),
1728 map_tmp_end = patch_to_global_tria_map_tmp.end();
1732 for (; map_tmp_it != map_tmp_end; ++map_tmp_it)
1733 patch_to_global_tria_map[map_tmp_it->first] = map_tmp_it->second;
1738 template <
int dim,
int spacedim>
1741 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
1755 std::set<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
1756 dof_to_set_of_cells_map;
1758 std::vector<types::global_dof_index> local_dof_indices;
1759 std::vector<types::global_dof_index> local_face_dof_indices;
1760 std::vector<types::global_dof_index> local_line_dof_indices;
1764 std::vector<bool> user_flags;
1769 std::map<typename DoFHandler<dim, spacedim>::active_line_iterator,
1771 lines_to_parent_lines_map;
1778 .clear_user_flags();
1783 endc = dof_handler.
end();
1784 for (; cell != endc; ++cell)
1790 if (cell->is_artificial() ==
false)
1792 for (
unsigned int l = 0; l < cell->n_lines(); ++l)
1793 if (cell->line(l)->has_children())
1794 for (
unsigned int c = 0; c < cell->line(l)->n_children();
1797 lines_to_parent_lines_map[cell->line(l)->child(c)] =
1801 cell->line(l)->child(c)->set_user_flag();
1816 endc = dof_handler.
end();
1817 for (; cell != endc; ++cell)
1822 if (cell->is_artificial() ==
false)
1824 const unsigned int n_dofs_per_cell =
1826 local_dof_indices.resize(n_dofs_per_cell);
1830 cell->get_dof_indices(local_dof_indices);
1831 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
1832 dof_to_set_of_cells_map[local_dof_indices[i]].insert(cell);
1842 for (
const unsigned int f : cell->face_indices())
1844 if (cell->face(f)->has_children())
1846 for (
unsigned int c = 0; c < cell->face(f)->n_children();
1861 Assert(cell->face(f)->child(c)->has_children() ==
false,
1864 const unsigned int n_dofs_per_face =
1866 local_face_dof_indices.resize(n_dofs_per_face);
1868 cell->face(f)->child(c)->get_dof_indices(
1869 local_face_dof_indices);
1870 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
1871 dof_to_set_of_cells_map[local_face_dof_indices[i]]
1875 else if ((cell->face(f)->at_boundary() ==
false) &&
1876 (cell->neighbor_is_coarser(f)))
1893 std::pair<unsigned int, unsigned int>
1894 neighbor_face_no_subface_no =
1895 cell->neighbor_of_coarser_neighbor(f);
1896 unsigned int face_no = neighbor_face_no_subface_no.first;
1897 unsigned int subface = neighbor_face_no_subface_no.second;
1899 const unsigned int n_dofs_per_face =
1901 local_face_dof_indices.resize(n_dofs_per_face);
1903 cell->neighbor(f)->face(face_no)->get_dof_indices(
1904 local_face_dof_indices);
1905 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
1906 dof_to_set_of_cells_map[local_face_dof_indices[i]].insert(
1911 for (
unsigned int c = 0;
1912 c < cell->neighbor(f)->face(face_no)->n_children();
1918 const unsigned int n_dofs_per_face =
1920 local_face_dof_indices.resize(n_dofs_per_face);
1925 ->has_children() ==
false,
1930 ->get_dof_indices(local_face_dof_indices);
1931 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
1932 dof_to_set_of_cells_map[local_face_dof_indices[i]]
1949 for (
unsigned int l = 0; l < cell->n_lines(); ++l)
1951 if (cell->line(l)->has_children())
1953 for (
unsigned int c = 0;
1954 c < cell->line(l)->n_children();
1957 Assert(cell->line(l)->child(c)->has_children() ==
1963 const unsigned int n_dofs_per_line =
1966 local_line_dof_indices.resize(n_dofs_per_line);
1968 cell->line(l)->child(c)->get_dof_indices(
1969 local_line_dof_indices);
1970 for (
unsigned int i = 0; i < n_dofs_per_line; ++i)
1971 dof_to_set_of_cells_map[local_line_dof_indices[i]]
1979 else if (cell->line(l)->user_flag_set() ==
true)
1983 lines_to_parent_lines_map[cell->line(l)];
1988 const unsigned int n_dofs_per_line =
1991 local_line_dof_indices.resize(n_dofs_per_line);
1993 parent_line->get_dof_indices(local_line_dof_indices);
1994 for (
unsigned int i = 0; i < n_dofs_per_line; ++i)
1995 dof_to_set_of_cells_map[local_line_dof_indices[i]]
1998 for (
unsigned int c = 0; c < parent_line->n_children();
2001 Assert(parent_line->child(c)->has_children() ==
2005 const unsigned int n_dofs_per_line =
2008 local_line_dof_indices.resize(n_dofs_per_line);
2010 parent_line->child(c)->get_dof_indices(
2011 local_line_dof_indices);
2012 for (
unsigned int i = 0; i < n_dofs_per_line; ++i)
2013 dof_to_set_of_cells_map[local_line_dof_indices[i]]
2030 .load_user_flags(user_flags);
2037 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
2038 dof_to_cell_patches;
2042 std::set<typename DoFHandler<dim, spacedim>::active_cell_iterator>>::
2043 iterator it = dof_to_set_of_cells_map.begin(),
2044 it_end = dof_to_set_of_cells_map.end();
2045 for (; it != it_end; ++it)
2046 dof_to_cell_patches[it->first].assign(it->second.begin(),
2049 return dof_to_cell_patches;
2055 template <
typename CellIterator>
2058 std::set<std::pair<CellIterator, unsigned int>> &pairs1,
2061 const unsigned int direction,
2063 const ::Tensor<1, CellIterator::AccessorType::space_dimension>
2067 static const int space_dim = CellIterator::AccessorType::space_dimension;
2073 constexpr int dim = CellIterator::AccessorType::dimension;
2074 constexpr int spacedim = CellIterator::AccessorType::space_dimension;
2079 if (!(((pairs1.size() > 0) &&
2080 (
dynamic_cast<const parallel::fullydistributed::
2081 Triangulation<dim, spacedim> *
>(
2082 &pairs1.begin()->first->get_triangulation()) !=
nullptr)) ||
2083 ((pairs2.size() > 0) &&
2086 *
>(&pairs2.begin()->first->get_triangulation()) !=
nullptr))))
2087 Assert(pairs1.size() == pairs2.size(),
2088 ExcMessage(
"Unmatched faces on periodic boundaries"));
2092 unsigned int n_matches = 0;
2095 std::bitset<3> orientation;
2096 using PairIterator =
2097 typename std::set<std::pair<CellIterator, unsigned int>>::const_iterator;
2098 for (PairIterator it1 = pairs1.begin(); it1 != pairs1.end(); ++it1)
2100 for (PairIterator it2 = pairs2.begin(); it2 != pairs2.end(); ++it2)
2102 const CellIterator cell1 = it1->first;
2103 const CellIterator cell2 = it2->first;
2104 const unsigned int face_idx1 = it1->second;
2105 const unsigned int face_idx2 = it2->second;
2107 cell1->face(face_idx1),
2108 cell2->face(face_idx2),
2117 {cell1, cell2}, {face_idx1, face_idx2}, orientation, matrix};
2118 matched_pairs.push_back(matched_face);
2132 constexpr int dim = CellIterator::AccessorType::dimension;
2133 constexpr int spacedim = CellIterator::AccessorType::space_dimension;
2134 if (!(((pairs1.size() > 0) &&
2135 (
dynamic_cast<const parallel::fullydistributed::
2136 Triangulation<dim, spacedim> *
>(
2137 &pairs1.begin()->first->get_triangulation()) !=
nullptr)) ||
2138 ((pairs2.size() > 0) &&
2141 *
>(&pairs2.begin()->first->get_triangulation()) !=
nullptr))))
2142 AssertThrow(n_matches == pairs1.size() && pairs2.size() == 0,
2143 ExcMessage(
"Unmatched faces on periodic boundaries"));
2149 template <
typename MeshType>
2152 const MeshType & mesh,
2153 const
types::boundary_id b_id,
2154 const
unsigned int direction,
2157 const
Tensor<1, MeshType::space_dimension> &offset,
2160 static const int dim = MeshType::dimension;
2161 static const int space_dim = MeshType::space_dimension;
2171 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs1;
2172 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs2;
2174 for (
typename MeshType::cell_iterator cell = mesh.begin(0);
2175 cell != mesh.end(0);
2178 const typename MeshType::face_iterator face_1 =
2179 cell->face(2 * direction);
2180 const typename MeshType::face_iterator face_2 =
2181 cell->face(2 * direction + 1);
2183 if (face_1->at_boundary() && face_1->boundary_id() == b_id)
2185 const std::pair<typename MeshType::cell_iterator, unsigned int>
2186 pair1 = std::make_pair(cell, 2 * direction);
2187 pairs1.insert(pair1);
2190 if (face_2->at_boundary() && face_2->boundary_id() == b_id)
2192 const std::pair<typename MeshType::cell_iterator, unsigned int>
2193 pair2 = std::make_pair(cell, 2 * direction + 1);
2194 pairs2.insert(pair2);
2198 Assert(pairs1.size() == pairs2.size(),
2199 ExcMessage(
"Unmatched faces on periodic boundaries"));
2201 Assert(pairs1.size() > 0,
2202 ExcMessage(
"No new periodic face pairs have been found. "
2203 "Are you sure that you've selected the correct boundary "
2204 "id's and that the coarsest level mesh is colorized?"));
2207 const unsigned int size_old = matched_pairs.size();
2212 pairs1, pairs2, direction, matched_pairs, offset, matrix);
2216 const unsigned int size_new = matched_pairs.size();
2217 for (
unsigned int i = size_old; i < size_new; ++i)
2219 Assert(matched_pairs[i].orientation == 1,
2221 "Found a face match with non standard orientation. "
2222 "This function is only suitable for meshes with cells "
2223 "in default orientation"));
2230 template <
typename MeshType>
2233 const MeshType & mesh,
2234 const
types::boundary_id b_id1,
2235 const
types::boundary_id b_id2,
2236 const
unsigned int direction,
2239 const
Tensor<1, MeshType::space_dimension> &offset,
2242 static const int dim = MeshType::dimension;
2243 static const int space_dim = MeshType::space_dimension;
2251 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs1;
2252 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs2;
2254 for (
typename MeshType::cell_iterator cell = mesh.begin(0);
2255 cell != mesh.end(0);
2258 for (
const unsigned int i : cell->face_indices())
2260 const typename MeshType::face_iterator face = cell->face(i);
2261 if (face->at_boundary() && face->boundary_id() == b_id1)
2263 const std::pair<typename MeshType::cell_iterator, unsigned int>
2264 pair1 = std::make_pair(cell, i);
2265 pairs1.insert(pair1);
2268 if (face->at_boundary() && face->boundary_id() == b_id2)
2270 const std::pair<typename MeshType::cell_iterator, unsigned int>
2271 pair2 = std::make_pair(cell, i);
2272 pairs2.insert(pair2);
2283 if (!(((pairs1.size() > 0) &&
2286 *
>(&pairs1.begin()->first->get_triangulation()) !=
nullptr)) ||
2287 ((pairs2.size() > 0) &&
2290 *
>(&pairs2.begin()->first->get_triangulation()) !=
nullptr))))
2291 Assert(pairs1.size() == pairs2.size(),
2292 ExcMessage(
"Unmatched faces on periodic boundaries"));
2295 (pairs1.size() > 0 ||
2298 &mesh.begin()->get_triangulation()) !=
nullptr)),
2299 ExcMessage(
"No new periodic face pairs have been found. "
2300 "Are you sure that you've selected the correct boundary "
2301 "id's and that the coarsest level mesh is colorized?"));
2305 pairs1, pairs2, direction, matched_pairs, offset, matrix);
2319 template <
int spacedim>
2323 const unsigned int direction,
2333 if (matrix.m() == spacedim)
2334 for (
unsigned int i = 0; i < spacedim; ++i)
2335 for (
unsigned int j = 0; j < spacedim; ++j)
2336 distance(i) += matrix(i, j) * point1(j);
2340 distance += offset - point2;
2342 for (
unsigned int i = 0; i < spacedim; ++i)
2348 if (
std::abs(distance(i)) > 1.e-10)
2374 std::array<unsigned int, GeometryInfo<1>::vertices_per_face>;
2375 static inline std::bitset<3>
2387 std::array<unsigned int, GeometryInfo<2>::vertices_per_face>;
2388 static inline std::bitset<3>
2396 static const MATCH_T m_tff = {{0, 1}};
2397 if (matching == m_tff)
2399 static const MATCH_T m_ttf = {{1, 0}};
2400 if (matching == m_ttf)
2413 std::array<unsigned int, GeometryInfo<3>::vertices_per_face>;
2415 static inline std::bitset<3>
2423 static const MATCH_T m_tff = {{0, 1, 2, 3}};
2424 if (matching == m_tff)
2426 static const MATCH_T m_tft = {{1, 3, 0, 2}};
2427 if (matching == m_tft)
2429 static const MATCH_T m_ttf = {{3, 2, 1, 0}};
2430 if (matching == m_ttf)
2432 static const MATCH_T m_ttt = {{2, 0, 3, 1}};
2433 if (matching == m_ttt)
2435 static const MATCH_T m_fff = {{0, 2, 1, 3}};
2436 if (matching == m_fff)
2438 static const MATCH_T m_fft = {{2, 3, 0, 1}};
2439 if (matching == m_fft)
2441 static const MATCH_T m_ftf = {{3, 1, 2, 0}};
2442 if (matching == m_ftf)
2444 static const MATCH_T m_ftt = {{1, 0, 3, 2}};
2445 if (matching == m_ftt)
2456 template <
typename FaceIterator>
2459 std::bitset<3> & orientation,
2460 const FaceIterator & face1,
2461 const FaceIterator & face2,
2462 const unsigned int direction,
2466 Assert(matrix.m() == matrix.n(),
2467 ExcMessage(
"The supplied matrix must be a square matrix"));
2469 static const int dim = FaceIterator::AccessorType::dimension;
2473 std::array<unsigned int, GeometryInfo<dim>::vertices_per_face> matching;
2477 std::set<unsigned int> face2_vertices;
2478 for (
unsigned int i = 0; i < face1->n_vertices(); ++i)
2479 face2_vertices.insert(i);
2481 for (
unsigned int i = 0; i < face1->n_vertices(); ++i)
2483 for (std::set<unsigned int>::iterator it = face2_vertices.begin();
2484 it != face2_vertices.end();
2494 face2_vertices.erase(it);
2501 if (face2_vertices.empty())
2504 return face2_vertices.empty();
2509 template <
typename FaceIterator>
2512 const FaceIterator & face1,
2513 const FaceIterator & face2,
2514 const unsigned int direction,
2519 std::bitset<3> dummy;
2528#include "grid_tools_dof_handlers.inst"
cell_iterator end() const
const hp::FECollection< dim, spacedim > & get_fe_collection() const
const FiniteElement< dim, spacedim > & get_fe(const types::fe_index index=0) const
const Triangulation< dim, spacedim > & get_triangulation() const
active_cell_iterator begin_active(const unsigned int level=0) 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
Abstract base class for mapping classes.
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
numbers::NumberTraits< Number >::real_type distance_square(const Point< dim, Number > &p) const
cell_iterator begin(const unsigned int level=0) const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
unsigned int n_active_cells() const
void save_user_flags(std::ostream &out) const
const std::vector< Point< spacedim > > & get_vertices() const
cell_iterator end() const
virtual void execute_coarsening_and_refinement()
unsigned int n_cells() const
const std::vector< bool > & get_used_vertices() const
Triangulation< dim, spacedim > & get_triangulation()
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_CXX20_REQUIRES(condition)
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcVertexNotUsed(unsigned int arg1)
#define AssertDimension(dim1, dim2)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
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::active_cell_iterator active_cell_iterator
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
static const unsigned int invalid_unsigned_int
typename type_identity< T >::type type_identity_t
::VectorizedArray< Number, width > min(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
unsigned int global_dof_index
const ::Triangulation< dim, spacedim > & tria