55 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
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>
134 const MeshType<dim, spacedim> &mesh,
136 const std::vector<bool> & marked_vertices)
151 marked_vertices.size() == 0,
153 marked_vertices.size()));
163 marked_vertices.size() == 0 ||
164 std::equal(marked_vertices.begin(),
165 marked_vertices.end(),
167 [](
bool p,
bool q) { return !p || q; }),
169 "marked_vertices should be a subset of used vertices in the triangulation "
170 "but marked_vertices contains one or more vertices that are not used vertices!"));
173 if (marked_vertices.size())
176 if (marked_vertices[it->first] ==
false)
191 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
193 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
196 typename ::internal::
197 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
200 const unsigned int vertex)
206 Assert(mesh.get_triangulation().get_used_vertices()[vertex],
213 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
216 typename ::internal::
217 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
218 cell = mesh.begin_active(),
254 for (; cell != endc; ++cell)
256 for (
const unsigned int v : cell->vertex_indices())
257 if (cell->vertex_index(v) == vertex)
269 for (
const auto face :
270 cell->reference_cell().faces_for_given_vertex(v))
271 if (!cell->at_boundary(face) &&
272 cell->neighbor(face)->is_active())
294 for (
unsigned int e = 0;
e < cell->n_lines(); ++
e)
295 if (cell->line(
e)->has_children())
299 if (cell->line(
e)->child(0)->vertex_index(1) == vertex)
323 typename ::internal::
324 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>(
332 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
334 find_active_cell_around_point_internal(
335 const MeshType<dim, spacedim> &mesh,
337 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>
339 std::set<
typename MeshType<dim, spacedim>::active_cell_iterator>
343 typename ::internal::
344 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
347 typename ::internal::
348 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type>
353 using cell_iterator =
354 typename MeshType<dim, spacedim>::active_cell_iterator;
356 using cell_iterator = typename ::internal::
357 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type;
365 std::set<cell_iterator> adjacent_cells_new;
367 typename std::set<cell_iterator>::const_iterator cell =
371 for (; cell != endc; ++cell)
373 std::vector<cell_iterator> active_neighbors;
374 get_active_neighbors<MeshType<dim, spacedim>>(*cell,
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>
407 typename MeshType<dim, spacedim>::active_cell_iterator
409 typename ::internal::
410 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type
414 const std::vector<bool> & marked_vertices,
415 const double tolerance)
417 return find_active_cell_around_point<dim, MeshType, spacedim>(
428 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
430 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
Point<dim>>
432 std::pair<typename ::internal::
433 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
437 const MeshType<dim, spacedim> &mesh,
439 const std::vector<bool> & marked_vertices,
440 const double tolerance)
442 using active_cell_iterator = typename ::internal::
443 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type;
449 double best_distance = tolerance;
451 std::pair<active_cell_iterator, Point<dim>> best_cell;
454 best_cell.first = mesh.end();
458 std::vector<active_cell_iterator> adjacent_cells_tmp =
467 std::set<active_cell_iterator>
adjacent_cells(adjacent_cells_tmp.begin(),
468 adjacent_cells_tmp.end());
469 std::set<active_cell_iterator> searched_cells;
478 mesh.get_triangulation().n_active_cells();
480 unsigned int cells_searched = 0;
483 typename std::set<active_cell_iterator>::const_iterator
486 for (; cell != endc; ++cell)
488 if ((*cell)->is_artificial() ==
false)
504 if ((dist < best_distance) ||
505 ((dist == best_distance) &&
506 ((*cell)->level() > best_level)))
509 best_distance = dist;
510 best_level = (*cell)->level();
511 best_cell = std::make_pair(*cell, p_cell);
516 spacedim>::ExcTransformationFailed &)
539 if (marked_vertices.size() > 0)
551 find_active_cell_around_point_internal<dim, MeshType, spacedim>(
561 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
563 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
566 std::vector<std::pair<
567 typename ::internal::
568 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
572 const MeshType<dim, spacedim> &mesh,
574 const double tolerance,
575 const std::vector<bool> &marked_vertices)
578 mapping, mesh, p, marked_vertices, tolerance);
580 if (cell_and_point.first == mesh.end())
584 mapping, mesh, p, tolerance, cell_and_point);
589 template <
int dim,
template <
int,
int>
class MeshType,
int spacedim>
591 std::vector<std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
594 std::vector<std::pair<
595 typename ::internal::
596 ActiveCellIterator<dim, spacedim, MeshType<dim, spacedim>>::type,
601 const MeshType<dim, spacedim> &mesh,
603 const double tolerance,
604 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
608 std::pair<typename MeshType<dim, spacedim>::active_cell_iterator,
613 cells_and_points.push_back(first_cell);
617 const Point<dim> unit_point = cells_and_points.front().second;
618 const auto my_cell = cells_and_points.front().first;
620 unsigned int n_dirs_at_threshold = 0;
622 for (
unsigned int d = 0;
d < dim; ++
d)
624 distance_to_center[
d] =
std::abs(unit_point[
d] - 0.5);
625 if (distance_to_center[
d] > 0.5 - tolerance)
627 ++n_dirs_at_threshold;
628 last_point_at_threshold =
d;
632 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
635 if (n_dirs_at_threshold == 1)
637 unsigned int neighbor_index =
638 2 * last_point_at_threshold +
639 (unit_point[last_point_at_threshold] > 0.5 ? 1 : 0);
640 if (!my_cell->at_boundary(neighbor_index))
641 cells_to_add.push_back(my_cell->neighbor(neighbor_index));
644 else if (n_dirs_at_threshold == dim)
646 unsigned int local_vertex_index = 0;
647 for (
unsigned int d = 0;
d < dim; ++
d)
648 local_vertex_index += (unit_point[
d] > 0.5 ? 1 : 0) <<
d;
649 std::vector<typename MeshType<dim, spacedim>::active_cell_iterator>
651 mesh, my_cell->vertex_index(local_vertex_index));
652 for (
const auto &cell : cells)
654 cells_to_add.push_back(cell);
661 else if (n_dirs_at_threshold == 2)
664 unsigned int count_vertex_indices = 0;
666 for (
unsigned int d = 0;
d < dim; ++
d)
668 if (distance_to_center[
d] > 0.5 - tolerance)
672 unit_point[
d] > 0.5 ? 1 : 0;
673 count_vertex_indices++;
683 const unsigned int first_vertex =
686 for (
unsigned int d = 0;
d < 2; ++
d)
690 my_cell->vertex_index(first_vertex + (
d << free_direction)));
691 for (
const auto &cell : tentative_cells)
693 bool cell_not_yet_present =
true;
694 for (
const auto &other_cell : cells_to_add)
695 if (cell == other_cell)
697 cell_not_yet_present =
false;
700 if (cell_not_yet_present)
701 cells_to_add.push_back(cell);
706 const double original_distance_to_unit_cell =
708 for (
const auto &cell : cells_to_add)
716 original_distance_to_unit_cell + tolerance)
717 cells_and_points.emplace_back(cell, p_unit);
724 cells_and_points.begin(),
725 cells_and_points.end(),
726 [](
const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
728 const std::pair<
typename MeshType<dim, spacedim>::active_cell_iterator,
731 return cells_and_points;
736 template <
class MeshType>
737 std::vector<typename MeshType::active_cell_iterator>
739 const MeshType &mesh,
740 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
743 std::vector<typename MeshType::active_cell_iterator> active_halo_layer;
744 std::vector<bool> locally_active_vertices_on_subdomain(
745 mesh.get_triangulation().n_vertices(),
false);
750 for (
const auto &cell : mesh.active_cell_iterators())
752 for (
const auto v : cell->vertex_indices())
753 locally_active_vertices_on_subdomain[cell->vertex_index(v)] =
true;
758 for (
const auto &cell : mesh.active_cell_iterators())
759 if (!predicate(cell))
760 for (
const auto v : cell->vertex_indices())
761 if (locally_active_vertices_on_subdomain[cell->vertex_index(v)] ==
764 active_halo_layer.push_back(cell);
768 return active_halo_layer;
773 template <
class MeshType>
774 std::vector<typename MeshType::cell_iterator>
776 const MeshType &mesh,
777 const std::function<
bool(
const typename MeshType::cell_iterator &)>
779 const unsigned int level)
781 std::vector<typename MeshType::cell_iterator> level_halo_layer;
782 std::vector<bool> locally_active_vertices_on_level_subdomain(
783 mesh.get_triangulation().n_vertices(),
false);
788 for (
typename MeshType::cell_iterator cell = mesh.begin(
level);
789 cell != mesh.end(
level);
792 for (
const unsigned int v : cell->vertex_indices())
793 locally_active_vertices_on_level_subdomain[cell->vertex_index(v)] =
799 for (
typename MeshType::cell_iterator cell = mesh.begin(
level);
800 cell != mesh.end(
level);
802 if (!predicate(cell))
803 for (
const unsigned int v : cell->vertex_indices())
804 if (locally_active_vertices_on_level_subdomain[cell->vertex_index(
807 level_halo_layer.push_back(cell);
811 return level_halo_layer;
817 template <
class MeshType>
819 contains_locally_owned_cells(
820 const std::vector<typename MeshType::active_cell_iterator> &cells)
822 for (
typename std::vector<
823 typename MeshType::active_cell_iterator>::const_iterator it =
828 if ((*it)->is_locally_owned())
834 template <
class MeshType>
836 contains_artificial_cells(
837 const std::vector<typename MeshType::active_cell_iterator> &cells)
839 for (
typename std::vector<
840 typename MeshType::active_cell_iterator>::const_iterator it =
845 if ((*it)->is_artificial())
854 template <
class MeshType>
855 std::vector<typename MeshType::active_cell_iterator>
858 std::function<
bool(
const typename MeshType::active_cell_iterator &)>
861 const std::vector<typename MeshType::active_cell_iterator>
866 Assert(contains_locally_owned_cells<MeshType>(active_halo_layer) ==
false,
867 ExcMessage(
"Halo layer contains locally owned cells"));
868 Assert(contains_artificial_cells<MeshType>(active_halo_layer) ==
false,
869 ExcMessage(
"Halo layer contains artificial cells"));
871 return active_halo_layer;
876 template <
class MeshType>
877 std::vector<typename MeshType::active_cell_iterator>
879 const MeshType &mesh,
880 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
882 const double layer_thickness)
884 std::vector<typename MeshType::active_cell_iterator>
885 subdomain_boundary_cells, active_cell_layer_within_distance;
886 std::vector<bool> vertices_outside_subdomain(
887 mesh.get_triangulation().n_vertices(),
false);
889 const unsigned int spacedim = MeshType::space_dimension;
891 unsigned int n_non_predicate_cells = 0;
899 for (
const auto &cell : mesh.active_cell_iterators())
900 if (!predicate(cell))
902 for (
const unsigned int v : cell->vertex_indices())
903 vertices_outside_subdomain[cell->vertex_index(v)] =
true;
904 n_non_predicate_cells++;
911 if (n_non_predicate_cells == 0 ||
912 n_non_predicate_cells == mesh.get_triangulation().n_active_cells())
913 return std::vector<typename MeshType::active_cell_iterator>();
917 for (
const auto &cell : mesh.active_cell_iterators())
920 for (
const unsigned int v : cell->vertex_indices())
921 if (vertices_outside_subdomain[cell->vertex_index(v)] ==
true)
923 subdomain_boundary_cells.push_back(cell);
936 for (
unsigned int d = 0;
d < spacedim; ++
d)
938 bounding_box.first[
d] -= (layer_thickness + DOUBLE_EPSILON);
939 bounding_box.second[
d] += (layer_thickness + DOUBLE_EPSILON);
942 std::vector<Point<spacedim>>
943 subdomain_boundary_cells_centers;
946 subdomain_boundary_cells_radii;
948 subdomain_boundary_cells_centers.reserve(subdomain_boundary_cells.size());
949 subdomain_boundary_cells_radii.reserve(subdomain_boundary_cells.size());
951 for (
typename std::vector<typename MeshType::active_cell_iterator>::
952 const_iterator subdomain_boundary_cell_iterator =
953 subdomain_boundary_cells.begin();
954 subdomain_boundary_cell_iterator != subdomain_boundary_cells.end();
955 ++subdomain_boundary_cell_iterator)
957 const std::pair<Point<spacedim>,
double>
958 &subdomain_boundary_cell_enclosing_ball =
959 (*subdomain_boundary_cell_iterator)->enclosing_ball();
961 subdomain_boundary_cells_centers.push_back(
962 subdomain_boundary_cell_enclosing_ball.first);
963 subdomain_boundary_cells_radii.push_back(
964 subdomain_boundary_cell_enclosing_ball.second);
966 AssertThrow(subdomain_boundary_cells_radii.size() ==
967 subdomain_boundary_cells_centers.size(),
976 for (
const auto &cell : mesh.active_cell_iterators())
982 const std::pair<Point<spacedim>,
double> &cell_enclosing_ball =
983 cell->enclosing_ball();
986 cell_enclosing_ball.first;
987 const double cell_enclosing_ball_radius = cell_enclosing_ball.second;
989 bool cell_inside =
true;
991 for (
unsigned int d = 0;
d < spacedim; ++
d)
993 (cell_enclosing_ball_center[
d] + cell_enclosing_ball_radius >
994 bounding_box.first[
d]) &&
995 (cell_enclosing_ball_center[
d] - cell_enclosing_ball_radius <
996 bounding_box.second[
d]);
1002 for (
unsigned int i = 0; i < subdomain_boundary_cells_radii.size();
1005 subdomain_boundary_cells_centers[i]) <
1006 Utilities::fixed_power<2>(cell_enclosing_ball_radius +
1007 subdomain_boundary_cells_radii[i] +
1008 layer_thickness + DOUBLE_EPSILON))
1010 active_cell_layer_within_distance.push_back(cell);
1015 return active_cell_layer_within_distance;
1020 template <
class MeshType>
1021 std::vector<typename MeshType::active_cell_iterator>
1023 const double layer_thickness)
1026 std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1027 predicate(locally_owned_cell_predicate);
1029 const std::vector<typename MeshType::active_cell_iterator>
1030 ghost_cell_layer_within_distance =
1038 contains_locally_owned_cells<MeshType>(
1039 ghost_cell_layer_within_distance) ==
false,
1041 "Ghost cells within layer_thickness contains locally owned cells."));
1043 contains_artificial_cells<MeshType>(ghost_cell_layer_within_distance) ==
1046 "Ghost cells within layer_thickness contains artificial cells. "
1047 "The function compute_ghost_cell_layer_within_distance "
1048 "is probably called while using parallel::distributed::Triangulation. "
1049 "In such case please refer to the description of this function."));
1051 return ghost_cell_layer_within_distance;
1056 template <
class MeshType>
1059 const MeshType &mesh,
1060 const std::function<
bool(
const typename MeshType::active_cell_iterator &)>
1063 std::vector<bool> locally_active_vertices_on_subdomain(
1064 mesh.get_triangulation().n_vertices(),
false);
1066 const unsigned int spacedim = MeshType::space_dimension;
1073 for (
const auto &cell : mesh.active_cell_iterators())
1074 if (predicate(cell))
1076 minp = cell->center();
1077 maxp = cell->center();
1083 for (
const auto &cell : mesh.active_cell_iterators())
1084 if (predicate(cell))
1085 for (
const unsigned int v : cell->vertex_indices())
1086 if (locally_active_vertices_on_subdomain[cell->vertex_index(v)] ==
1089 locally_active_vertices_on_subdomain[cell->vertex_index(v)] =
1091 for (
unsigned int d = 0;
d < spacedim; ++
d)
1098 return std::make_pair(minp, maxp);
1103 template <
typename MeshType>
1104 std::list<std::pair<
typename MeshType::cell_iterator,
1105 typename MeshType::cell_iterator>>
1109 ExcMessage(
"The two meshes must be represent triangulations that "
1110 "have the same coarse meshes"));
1117 bool remove_ghost_cells =
false;
1118#ifdef DEAL_II_WITH_MPI
1120 constexpr int dim = MeshType::dimension;
1121 constexpr int spacedim = MeshType::space_dimension;
1123 *
>(&mesh_1.get_triangulation()) !=
nullptr ||
1125 *
>(&mesh_2.get_triangulation()) !=
nullptr)
1127 Assert(&mesh_1.get_triangulation() == &mesh_2.get_triangulation(),
1128 ExcMessage(
"This function can only be used with meshes "
1129 "corresponding to distributed Triangulations when "
1130 "both Triangulations are equal."));
1131 remove_ghost_cells =
true;
1143 using CellList = std::list<std::pair<
typename MeshType::cell_iterator,
1144 typename MeshType::cell_iterator>>;
1148 typename MeshType::cell_iterator cell_1 = mesh_1.begin(0),
1149 cell_2 = mesh_2.begin(0);
1150 for (; cell_1 != mesh_1.end(0); ++cell_1, ++cell_2)
1151 cell_list.emplace_back(cell_1, cell_2);
1154 typename CellList::iterator cell_pair = cell_list.begin();
1155 while (cell_pair != cell_list.end())
1159 if (cell_pair->first->has_children() &&
1160 cell_pair->second->has_children())
1162 Assert(cell_pair->first->refinement_case() ==
1163 cell_pair->second->refinement_case(),
1165 for (
unsigned int c = 0; c < cell_pair->first->n_children(); ++c)
1166 cell_list.emplace_back(cell_pair->first->child(c),
1167 cell_pair->second->child(c));
1172 const auto previous_cell_pair = cell_pair;
1174 cell_list.erase(previous_cell_pair);
1179 if (remove_ghost_cells &&
1180 ((cell_pair->first->is_active() &&
1181 !cell_pair->first->is_locally_owned()) ||
1182 (cell_pair->second->is_active() &&
1183 !cell_pair->second->is_locally_owned())))
1186 const auto previous_cell_pair = cell_pair;
1188 cell_list.erase(previous_cell_pair);
1197 for (cell_pair = cell_list.begin(); cell_pair != cell_list.end();
1199 Assert(cell_pair->first->is_active() || cell_pair->second->is_active() ||
1200 (cell_pair->first->refinement_case() !=
1201 cell_pair->second->refinement_case()),
1209 template <
int dim,
int spacedim>
1226 endc = mesh_1.
end(0);
1227 for (; cell_1 != endc; ++cell_1, ++cell_2)
1229 if (cell_1->
n_vertices() != cell_2->n_vertices())
1231 for (
const unsigned int v : cell_1->vertex_indices())
1232 if (cell_1->vertex(v) != cell_2->vertex(v))
1245 template <
typename MeshType>
1250 mesh_2.get_triangulation());
1255 template <
int dim,
int spacedim>
1256 std::pair<typename ::DoFHandler<dim, spacedim>::active_cell_iterator,
1262 const double tolerance)
1266 ExcMessage(
"Mapping collection needs to have either size 1 "
1267 "or size equal to the number of elements in "
1268 "the FECollection."));
1270 using cell_iterator =
1271 typename ::DoFHandler<dim, spacedim>::active_cell_iterator;
1273 std::pair<cell_iterator, Point<dim>> best_cell;
1277 if (mapping.
size() == 1)
1279 const std::vector<bool> marked_vertices = {};
1281 mapping[0], mesh, p, marked_vertices, tolerance);
1288 double best_distance = tolerance;
1289 int best_level = -1;
1296 std::vector<cell_iterator> adjacent_cells_tmp =
1304 std::set<cell_iterator>
adjacent_cells(adjacent_cells_tmp.begin(),
1305 adjacent_cells_tmp.end());
1306 std::set<cell_iterator> searched_cells;
1316 unsigned int cells_searched = 0;
1317 while (!found && cells_searched <
n_cells)
1319 typename std::set<cell_iterator>::const_iterator
1322 for (; cell != endc; ++cell)
1327 mapping[(*cell)->active_fe_index()]
1328 .transform_real_to_unit_cell(*cell, p);
1340 if (dist < best_distance || (dist == best_distance &&
1341 (*cell)->level() > best_level))
1344 best_distance = dist;
1345 best_level = (*cell)->level();
1346 best_cell = std::make_pair(*cell, p_cell);
1351 spacedim>::ExcTransformationFailed &)
1373 if (!found && cells_searched <
n_cells)
1375 find_active_cell_around_point_internal<dim,
1387 template <
class MeshType>
1388 std::vector<typename MeshType::active_cell_iterator>
1391 Assert(cell->is_locally_owned(),
1392 ExcMessage(
"This function only makes sense if the cell for "
1393 "which you are asking for a patch, is locally "
1396 std::vector<typename MeshType::active_cell_iterator> patch;
1397 patch.push_back(cell);
1398 for (
const unsigned int face_number : cell->face_indices())
1399 if (cell->face(face_number)->at_boundary() ==
false)
1401 if (cell->neighbor(face_number)->has_children() ==
false)
1402 patch.push_back(cell->neighbor(face_number));
1407 if (MeshType::dimension > 1)
1409 for (
unsigned int subface = 0;
1410 subface < cell->face(face_number)->n_children();
1413 cell->neighbor_child_on_subface(face_number, subface));
1419 typename MeshType::cell_iterator neighbor =
1420 cell->neighbor(face_number);
1421 while (neighbor->has_children())
1422 neighbor = neighbor->child(1 - face_number);
1424 Assert(neighbor->neighbor(1 - face_number) == cell,
1426 patch.push_back(neighbor);
1434 template <
class Container>
1435 std::vector<typename Container::cell_iterator>
1437 const std::vector<typename Container::active_cell_iterator> &patch)
1441 "Vector containing patch cells should not be an empty vector!"));
1445 int min_level = patch[0]->level();
1447 for (
unsigned int i = 0; i < patch.size(); ++i)
1449 std::set<typename Container::cell_iterator> uniform_cells;
1450 typename std::vector<
1451 typename Container::active_cell_iterator>::const_iterator patch_cell;
1453 for (patch_cell = patch.begin(); patch_cell != patch.end(); ++patch_cell)
1458 if ((*patch_cell)->level() == min_level)
1459 uniform_cells.insert(*patch_cell);
1466 typename Container::cell_iterator parent = *patch_cell;
1468 while (parent->level() > min_level)
1469 parent = parent->parent();
1470 uniform_cells.insert(parent);
1474 return std::vector<typename Container::cell_iterator>(uniform_cells.begin(),
1475 uniform_cells.end());
1480 template <
class Container>
1483 const std::vector<typename Container::active_cell_iterator> &patch,
1485 &local_triangulation,
1488 Container::space_dimension>::active_cell_iterator,
1489 typename Container::active_cell_iterator> &patch_to_global_tria_map)
1492 const std::vector<typename Container::cell_iterator> uniform_cells =
1493 get_cells_at_coarsest_common_level<Container>(patch);
1495 local_triangulation.
clear();
1496 std::vector<Point<Container::space_dimension>>
vertices;
1497 const unsigned int n_uniform_cells = uniform_cells.size();
1498 std::vector<CellData<Container::dimension>> cells(n_uniform_cells);
1501 typename std::vector<typename Container::cell_iterator>::const_iterator
1503 for (uniform_cell = uniform_cells.begin();
1504 uniform_cell != uniform_cells.end();
1507 for (
const unsigned int v : (*uniform_cell)->vertex_indices())
1510 (*uniform_cell)->vertex(v);
1511 bool repeat_vertex =
false;
1513 for (
unsigned int m = 0; m < i; ++m)
1517 repeat_vertex =
true;
1518 cells[k].vertices[v] = m;
1522 if (repeat_vertex ==
false)
1525 cells[k].vertices[v] = i;
1536 unsigned int index = 0;
1539 Container::space_dimension>::cell_iterator,
1540 typename Container::cell_iterator>
1541 patch_to_global_tria_map_tmp;
1543 Container::space_dimension>::cell_iterator
1544 coarse_cell = local_triangulation.
begin();
1545 coarse_cell != local_triangulation.
end();
1546 ++coarse_cell, ++index)
1548 patch_to_global_tria_map_tmp.insert(
1549 std::make_pair(coarse_cell, uniform_cells[index]));
1553 Assert(coarse_cell->center().distance(uniform_cells[index]->center()) <=
1554 1
e-15 * coarse_cell->diameter(),
1557 bool refinement_necessary;
1562 refinement_necessary =
false;
1563 for (
const auto &active_tria_cell :
1566 if (patch_to_global_tria_map_tmp[active_tria_cell]->has_children())
1568 active_tria_cell->set_refine_flag();
1569 refinement_necessary =
true;
1572 for (
unsigned int i = 0; i < patch.size(); ++i)
1577 if (patch_to_global_tria_map_tmp[active_tria_cell] ==
1582 for (
const unsigned int v :
1583 active_tria_cell->vertex_indices())
1584 active_tria_cell->vertex(v) = patch[i]->vertex(v);
1586 Assert(active_tria_cell->center().distance(
1587 patch_to_global_tria_map_tmp[active_tria_cell]
1589 1
e-15 * active_tria_cell->diameter(),
1592 active_tria_cell->set_user_flag();
1598 if (refinement_necessary)
1603 Container::dimension,
1604 Container::space_dimension>::cell_iterator cell =
1605 local_triangulation.
begin();
1606 cell != local_triangulation.
end();
1609 if (patch_to_global_tria_map_tmp.find(cell) !=
1610 patch_to_global_tria_map_tmp.end())
1612 if (cell->has_children())
1619 for (
unsigned int c = 0; c < cell->n_children(); ++c)
1621 if (patch_to_global_tria_map_tmp.find(cell->child(
1622 c)) == patch_to_global_tria_map_tmp.end())
1624 patch_to_global_tria_map_tmp.insert(
1627 patch_to_global_tria_map_tmp[cell]->child(
1649 patch_to_global_tria_map_tmp.erase(cell);
1655 while (refinement_necessary);
1661 Container::space_dimension>::cell_iterator
1662 cell = local_triangulation.
begin();
1663 cell != local_triangulation.
end();
1666 if (cell->user_flag_set())
1668 Assert(patch_to_global_tria_map_tmp.find(cell) !=
1669 patch_to_global_tria_map_tmp.end(),
1672 Assert(cell->center().distance(
1673 patch_to_global_tria_map_tmp[cell]->center()) <=
1674 1
e-15 * cell->diameter(),
1682 Container::space_dimension>::cell_iterator,
1683 typename Container::cell_iterator>::iterator
1684 map_tmp_it = patch_to_global_tria_map_tmp.
begin(),
1685 map_tmp_end = patch_to_global_tria_map_tmp.end();
1689 for (; map_tmp_it != map_tmp_end; ++map_tmp_it)
1690 patch_to_global_tria_map[map_tmp_it->first] = map_tmp_it->second;
1695 template <
int dim,
int spacedim>
1698 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
1712 std::set<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
1713 dof_to_set_of_cells_map;
1715 std::vector<types::global_dof_index> local_dof_indices;
1716 std::vector<types::global_dof_index> local_face_dof_indices;
1717 std::vector<types::global_dof_index> local_line_dof_indices;
1721 std::vector<bool> user_flags;
1726 std::map<typename DoFHandler<dim, spacedim>::active_line_iterator,
1728 lines_to_parent_lines_map;
1735 .clear_user_flags();
1740 endc = dof_handler.
end();
1741 for (; cell != endc; ++cell)
1747 if (cell->is_artificial() ==
false)
1749 for (
unsigned int l = 0;
l < cell->n_lines(); ++
l)
1750 if (cell->line(
l)->has_children())
1751 for (
unsigned int c = 0; c < cell->line(
l)->n_children();
1754 lines_to_parent_lines_map[cell->line(
l)->child(c)] =
1758 cell->line(
l)->child(c)->set_user_flag();
1773 endc = dof_handler.
end();
1774 for (; cell != endc; ++cell)
1779 if (cell->is_artificial() ==
false)
1781 const unsigned int n_dofs_per_cell =
1783 local_dof_indices.resize(n_dofs_per_cell);
1787 cell->get_dof_indices(local_dof_indices);
1788 for (
unsigned int i = 0; i < n_dofs_per_cell; ++i)
1789 dof_to_set_of_cells_map[local_dof_indices[i]].
insert(cell);
1799 for (
const unsigned int f : cell->face_indices())
1801 if (cell->face(f)->has_children())
1803 for (
unsigned int c = 0; c < cell->face(f)->n_children();
1818 Assert(cell->face(f)->child(c)->has_children() ==
false,
1821 const unsigned int n_dofs_per_face =
1823 local_face_dof_indices.resize(n_dofs_per_face);
1825 cell->face(f)->child(c)->get_dof_indices(
1826 local_face_dof_indices);
1827 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
1828 dof_to_set_of_cells_map[local_face_dof_indices[i]]
1832 else if ((cell->face(f)->at_boundary() ==
false) &&
1833 (cell->neighbor_is_coarser(f)))
1850 std::pair<unsigned int, unsigned int>
1851 neighbor_face_no_subface_no =
1852 cell->neighbor_of_coarser_neighbor(f);
1853 unsigned int face_no = neighbor_face_no_subface_no.first;
1854 unsigned int subface = neighbor_face_no_subface_no.second;
1856 const unsigned int n_dofs_per_face =
1858 local_face_dof_indices.resize(n_dofs_per_face);
1860 cell->neighbor(f)->face(face_no)->get_dof_indices(
1861 local_face_dof_indices);
1862 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
1863 dof_to_set_of_cells_map[local_face_dof_indices[i]].
insert(
1868 for (
unsigned int c = 0;
1869 c < cell->neighbor(f)->face(face_no)->n_children();
1875 const unsigned int n_dofs_per_face =
1877 local_face_dof_indices.resize(n_dofs_per_face);
1882 ->has_children() ==
false,
1887 ->get_dof_indices(local_face_dof_indices);
1888 for (
unsigned int i = 0; i < n_dofs_per_face; ++i)
1889 dof_to_set_of_cells_map[local_face_dof_indices[i]]
1906 for (
unsigned int l = 0;
l < cell->n_lines(); ++
l)
1908 if (cell->line(
l)->has_children())
1910 for (
unsigned int c = 0;
1911 c < cell->line(
l)->n_children();
1914 Assert(cell->line(
l)->child(c)->has_children() ==
1920 const unsigned int n_dofs_per_line =
1923 local_line_dof_indices.resize(n_dofs_per_line);
1925 cell->line(
l)->child(c)->get_dof_indices(
1926 local_line_dof_indices);
1927 for (
unsigned int i = 0; i < n_dofs_per_line; ++i)
1928 dof_to_set_of_cells_map[local_line_dof_indices[i]]
1936 else if (cell->line(
l)->user_flag_set() ==
true)
1940 lines_to_parent_lines_map[cell->line(
l)];
1945 const unsigned int n_dofs_per_line =
1948 local_line_dof_indices.resize(n_dofs_per_line);
1950 parent_line->get_dof_indices(local_line_dof_indices);
1951 for (
unsigned int i = 0; i < n_dofs_per_line; ++i)
1952 dof_to_set_of_cells_map[local_line_dof_indices[i]]
1955 for (
unsigned int c = 0; c < parent_line->n_children();
1958 Assert(parent_line->child(c)->has_children() ==
1962 const unsigned int n_dofs_per_line =
1965 local_line_dof_indices.resize(n_dofs_per_line);
1967 parent_line->child(c)->get_dof_indices(
1968 local_line_dof_indices);
1969 for (
unsigned int i = 0; i < n_dofs_per_line; ++i)
1970 dof_to_set_of_cells_map[local_line_dof_indices[i]]
1987 .load_user_flags(user_flags);
1994 std::vector<typename DoFHandler<dim, spacedim>::active_cell_iterator>>
1995 dof_to_cell_patches;
1999 std::set<typename DoFHandler<dim, spacedim>::active_cell_iterator>>::
2000 iterator it = dof_to_set_of_cells_map.begin(),
2001 it_end = dof_to_set_of_cells_map.end();
2002 for (; it != it_end; ++it)
2003 dof_to_cell_patches[it->first].assign(it->second.begin(),
2006 return dof_to_cell_patches;
2012 template <
typename CellIterator>
2015 std::set<std::pair<CellIterator, unsigned int>> &pairs1,
2018 const int direction,
2020 const ::Tensor<1, CellIterator::AccessorType::space_dimension>
2024 static const int space_dim = CellIterator::AccessorType::space_dimension;
2030 constexpr int dim = CellIterator::AccessorType::dimension;
2031 constexpr int spacedim = CellIterator::AccessorType::space_dimension;
2036 if (!(((pairs1.size() > 0) &&
2037 (
dynamic_cast<const parallel::fullydistributed::
2038 Triangulation<dim, spacedim> *
>(
2039 &pairs1.begin()->first->get_triangulation()) !=
nullptr)) ||
2040 ((pairs2.size() > 0) &&
2043 *
>(&pairs2.begin()->first->get_triangulation()) !=
nullptr))))
2044 Assert(pairs1.size() == pairs2.size(),
2045 ExcMessage(
"Unmatched faces on periodic boundaries"));
2049 unsigned int n_matches = 0;
2052 std::bitset<3> orientation;
2053 using PairIterator =
2054 typename std::set<std::pair<CellIterator, unsigned int>>::const_iterator;
2055 for (PairIterator it1 = pairs1.begin(); it1 != pairs1.end(); ++it1)
2057 for (PairIterator it2 = pairs2.begin(); it2 != pairs2.end(); ++it2)
2059 const CellIterator cell1 = it1->first;
2060 const CellIterator cell2 = it2->first;
2061 const unsigned int face_idx1 = it1->second;
2062 const unsigned int face_idx2 = it2->second;
2064 cell1->face(face_idx1),
2065 cell2->face(face_idx2),
2074 {cell1, cell2}, {face_idx1, face_idx2}, orientation,
matrix};
2075 matched_pairs.push_back(matched_face);
2089 constexpr int dim = CellIterator::AccessorType::dimension;
2090 constexpr int spacedim = CellIterator::AccessorType::space_dimension;
2091 if (!(((pairs1.size() > 0) &&
2092 (
dynamic_cast<const parallel::fullydistributed::
2093 Triangulation<dim, spacedim> *
>(
2094 &pairs1.begin()->first->get_triangulation()) !=
nullptr)) ||
2095 ((pairs2.size() > 0) &&
2098 *
>(&pairs2.begin()->first->get_triangulation()) !=
nullptr))))
2099 AssertThrow(n_matches == pairs1.size() && pairs2.size() == 0,
2100 ExcMessage(
"Unmatched faces on periodic boundaries"));
2106 template <
typename MeshType>
2109 const MeshType & mesh,
2111 const int direction,
2117 static const int dim = MeshType::dimension;
2118 static const int space_dim = MeshType::space_dimension;
2128 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs1;
2129 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs2;
2131 for (
typename MeshType::cell_iterator cell = mesh.begin(0);
2132 cell != mesh.end(0);
2135 const typename MeshType::face_iterator face_1 =
2136 cell->face(2 * direction);
2137 const typename MeshType::face_iterator face_2 =
2138 cell->face(2 * direction + 1);
2140 if (face_1->at_boundary() && face_1->boundary_id() == b_id)
2142 const std::pair<typename MeshType::cell_iterator, unsigned int>
2143 pair1 = std::make_pair(cell, 2 * direction);
2144 pairs1.insert(pair1);
2147 if (face_2->at_boundary() && face_2->boundary_id() == b_id)
2149 const std::pair<typename MeshType::cell_iterator, unsigned int>
2150 pair2 = std::make_pair(cell, 2 * direction + 1);
2151 pairs2.insert(pair2);
2155 Assert(pairs1.size() == pairs2.size(),
2156 ExcMessage(
"Unmatched faces on periodic boundaries"));
2158 Assert(pairs1.size() > 0,
2159 ExcMessage(
"No new periodic face pairs have been found. "
2160 "Are you sure that you've selected the correct boundary "
2161 "id's and that the coarsest level mesh is colorized?"));
2164 const unsigned int size_old = matched_pairs.size();
2169 pairs1, pairs2, direction, matched_pairs, offset,
matrix);
2173 const unsigned int size_new = matched_pairs.size();
2174 for (
unsigned int i = size_old; i < size_new; ++i)
2176 Assert(matched_pairs[i].orientation == 1,
2178 "Found a face match with non standard orientation. "
2179 "This function is only suitable for meshes with cells "
2180 "in default orientation"));
2187 template <
typename MeshType>
2190 const MeshType & mesh,
2193 const int direction,
2199 static const int dim = MeshType::dimension;
2200 static const int space_dim = MeshType::space_dimension;
2208 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs1;
2209 std::set<std::pair<typename MeshType::cell_iterator, unsigned int>> pairs2;
2211 for (
typename MeshType::cell_iterator cell = mesh.begin(0);
2212 cell != mesh.end(0);
2215 for (
unsigned int i : cell->face_indices())
2217 const typename MeshType::face_iterator face = cell->face(i);
2218 if (face->at_boundary() && face->boundary_id() == b_id1)
2220 const std::pair<typename MeshType::cell_iterator, unsigned int>
2221 pair1 = std::make_pair(cell, i);
2222 pairs1.insert(pair1);
2225 if (face->at_boundary() && face->boundary_id() == b_id2)
2227 const std::pair<typename MeshType::cell_iterator, unsigned int>
2228 pair2 = std::make_pair(cell, i);
2229 pairs2.insert(pair2);
2240 if (!(((pairs1.size() > 0) &&
2243 *
>(&pairs1.begin()->first->get_triangulation()) !=
nullptr)) ||
2244 ((pairs2.size() > 0) &&
2247 *
>(&pairs2.begin()->first->get_triangulation()) !=
nullptr))))
2248 Assert(pairs1.size() == pairs2.size(),
2249 ExcMessage(
"Unmatched faces on periodic boundaries"));
2252 (pairs1.size() > 0 ||
2255 &mesh.begin()->get_triangulation()) !=
nullptr)),
2256 ExcMessage(
"No new periodic face pairs have been found. "
2257 "Are you sure that you've selected the correct boundary "
2258 "id's and that the coarsest level mesh is colorized?"));
2262 pairs1, pairs2, direction, matched_pairs, offset,
matrix);
2276 template <
int spacedim>
2280 const int direction,
2290 if (
matrix.m() == spacedim)
2291 for (
int i = 0; i < spacedim; ++i)
2292 for (
int j = 0; j < spacedim; ++j)
2293 distance(i) +=
matrix(i, j) * point1(j);
2297 distance += offset - point2;
2299 for (
int i = 0; i < spacedim; ++i)
2305 if (
std::abs(distance(i)) > 1.e-10)
2331 std::array<unsigned int, GeometryInfo<1>::vertices_per_face>;
2332 static inline std::bitset<3>
2344 std::array<unsigned int, GeometryInfo<2>::vertices_per_face>;
2345 static inline std::bitset<3>
2353 static const MATCH_T m_tff = {{0, 1}};
2354 if (matching == m_tff)
2356 static const MATCH_T m_ttf = {{1, 0}};
2357 if (matching == m_ttf)
2370 std::array<unsigned int, GeometryInfo<3>::vertices_per_face>;
2372 static inline std::bitset<3>
2380 static const MATCH_T m_tff = {{0, 1, 2, 3}};
2381 if (matching == m_tff)
2383 static const MATCH_T m_tft = {{1, 3, 0, 2}};
2384 if (matching == m_tft)
2386 static const MATCH_T m_ttf = {{3, 2, 1, 0}};
2387 if (matching == m_ttf)
2389 static const MATCH_T m_ttt = {{2, 0, 3, 1}};
2390 if (matching == m_ttt)
2392 static const MATCH_T m_fff = {{0, 2, 1, 3}};
2393 if (matching == m_fff)
2395 static const MATCH_T m_fft = {{2, 3, 0, 1}};
2396 if (matching == m_fft)
2398 static const MATCH_T m_ftf = {{3, 1, 2, 0}};
2399 if (matching == m_ftf)
2401 static const MATCH_T m_ftt = {{1, 0, 3, 2}};
2402 if (matching == m_ftt)
2413 template <
typename FaceIterator>
2415 std::bitset<3> & orientation,
2416 const FaceIterator & face1,
2417 const FaceIterator & face2,
2418 const int direction,
2423 ExcMessage(
"The supplied matrix must be a square matrix"));
2425 static const int dim = FaceIterator::AccessorType::dimension;
2429 std::array<unsigned int, GeometryInfo<dim>::vertices_per_face> matching;
2433 std::set<unsigned int> face2_vertices;
2434 for (
unsigned int i = 0; i < face1->n_vertices(); ++i)
2435 face2_vertices.insert(i);
2437 for (
unsigned int i = 0; i < face1->n_vertices(); ++i)
2439 for (std::set<unsigned int>::iterator it = face2_vertices.begin();
2440 it != face2_vertices.end();
2450 face2_vertices.erase(it);
2457 if (face2_vertices.empty())
2460 return face2_vertices.empty();
2465 template <
typename FaceIterator>
2468 const FaceIterator & face1,
2469 const FaceIterator & face2,
2470 const int direction,
2475 std::bitset<3> dummy;
2484#include "grid_tools_dof_handlers.inst"
cell_iterator end() const
const FiniteElement< dim, spacedim > & get_fe(const unsigned int index=0) const
const hp::FECollection< dim, spacedim > & get_fe_collection() 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
virtual Point< dim > transform_real_to_unit_cell(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const Point< spacedim > &p) const =0
virtual bool preserves_vertex_locations() 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 n_vertices() const
unsigned int size() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
IteratorRange< active_cell_iterator > active_cell_iterators() const
__global__ void set(Number *val, const Number s, const size_type N)
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)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
@ matrix
Contents is actually a matrix.
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
SymmetricTensor< 2, dim, Number > d(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
SymmetricTensor< 2, dim, Number > epsilon(const Tensor< 2, dim, Number > &Grad_u)
unsigned int n_active_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
unsigned int n_cells(const internal::TriangulationImplementation::NumberCache< 1 > &c)
static const unsigned int invalid_unsigned_int
::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
static double distance_to_unit_cell(const Point< dim > &p)