17 #include <deal.II/base/timer.h> 18 #include <deal.II/base/utilities.h> 20 #include <deal.II/grid/grid_reordering.h> 21 #include <deal.II/grid/grid_tools.h> 29 DEAL_II_NAMESPACE_OPEN
44 CheapEdge(
const unsigned int v0,
const unsigned int v1)
54 operator<(
const CheapEdge &e)
const 56 return ((v0 <
e.v0) || ((v0 ==
e.v0) && (v1 <
e.v1)));
63 const unsigned int v0, v1;
79 std::set<CheapEdge> edges;
81 for (
typename std::vector<
CellData<dim>>::const_iterator c = cells.begin();
90 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
92 const CheapEdge reverse_edge(
95 if (edges.find(reverse_edge) != edges.end())
100 const CheapEdge correct_edge(
103 edges.insert(correct_edge);
126 static const unsigned int starter_edges[dim];
132 static const unsigned int n_other_parallel_edges = (1 << (dim - 1)) - 1;
134 [n_other_parallel_edges];
138 const unsigned int ParallelEdges<2>::starter_edges[2] = {0, 2};
141 const unsigned int ParallelEdges<2>::parallel_edges[4][1] = {{1},
147 const unsigned int ParallelEdges<3>::starter_edges[3] = {0, 2, 8};
150 const unsigned int ParallelEdges<3>::parallel_edges[12][3] = {
176 : cell_index(
numbers::invalid_unsigned_int)
177 , edge_within_cell(
numbers::invalid_unsigned_int)
183 AdjacentCell(
const unsigned int cell_index,
184 const unsigned int edge_within_cell)
185 : cell_index(cell_index)
186 , edge_within_cell(edge_within_cell)
190 unsigned int cell_index;
191 unsigned int edge_within_cell;
205 class AdjacentCells<2>
212 using const_iterator =
const AdjacentCell *;
223 push_back(
const AdjacentCell &adjacent_cell)
226 adjacent_cells[0] = adjacent_cell;
231 adjacent_cells[1] = adjacent_cell;
243 return adjacent_cells;
259 return adjacent_cells;
261 return adjacent_cells + 1;
263 return adjacent_cells + 2;
273 AdjacentCell adjacent_cells[2];
286 class AdjacentCells<3> :
public std::vector<AdjacentCell>
308 Edge(
const CellData<dim> &cell,
const unsigned int edge_number)
309 : orientation_status(not_oriented)
321 if (vertex_indices[0] > vertex_indices[1])
322 std::swap(vertex_indices[0], vertex_indices[1]);
330 operator<(const Edge<dim> &
e)
const 332 return ((vertex_indices[0] <
e.vertex_indices[0]) ||
333 ((vertex_indices[0] ==
e.vertex_indices[0]) &&
334 (vertex_indices[1] <
e.vertex_indices[1])));
341 operator==(
const Edge<dim> &e)
const 343 return ((vertex_indices[0] ==
e.vertex_indices[0]) &&
344 (vertex_indices[1] ==
e.vertex_indices[1]));
351 unsigned int vertex_indices[2];
357 enum OrientationStatus
364 OrientationStatus orientation_status;
370 AdjacentCells<dim> adjacent_cells;
387 Cell(
const CellData<dim> &c,
const std::vector<Edge<dim>> &edge_list)
389 for (
unsigned int i = 0; i < GeometryInfo<dim>::vertices_per_cell; ++i)
394 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
396 const Edge<dim>
e(c, l);
398 (std::lower_bound(edge_list.begin(), edge_list.end(),
e) -
432 class EdgeDeltaSet<2>
438 using const_iterator =
const unsigned int *;
464 insert(
const unsigned int edge_index)
467 edge_indices[0] = edge_index;
472 edge_indices[1] = edge_index;
499 return edge_indices + 1;
501 return edge_indices + 2;
508 unsigned int edge_indices[2];
525 class EdgeDeltaSet<3> :
public std::set<unsigned int>
535 std::vector<Edge<dim>>
542 std::vector<Edge<dim>> edge_list;
544 for (
unsigned int i = 0; i < cells.size(); ++i)
545 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
546 edge_list.emplace_back(cells[i], l);
549 std::sort(edge_list.begin(), edge_list.end());
550 edge_list.erase(std::unique(edge_list.begin(), edge_list.end()),
563 std::vector<Cell<dim>>
564 build_cells_and_connect_edges(
const std::vector<
CellData<dim>> &cells,
565 std::vector<Edge<dim>> & edges)
567 std::vector<Cell<dim>> cell_list;
568 cell_list.reserve(cells.size());
569 for (
unsigned int i = 0; i < cells.size(); ++i)
573 cell_list.emplace_back(cells[i], edges);
577 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
578 edges[cell_list.back().edge_indices[
l]].adjacent_cells.push_back(
594 get_next_unoriented_cell(
const std::vector<Cell<dim>> &cells,
595 const std::vector<Edge<dim>> &edges,
596 const unsigned int current_cell)
598 for (
unsigned int c = current_cell; c < cells.size(); ++c)
599 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
600 if (edges[cells[c].edge_indices[l]].orientation_status ==
601 Edge<dim>::not_oriented)
616 orient_one_set_of_parallel_edges(
const std::vector<Cell<dim>> &cells,
617 std::vector<Edge<dim>> & edges,
618 const unsigned int cell,
619 const unsigned int local_edge)
641 if (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0] ==
646 edges[cells[cell].edge_indices[local_edge]].orientation_status =
647 (dim == 2 ? Edge<dim>::backward : Edge<dim>::forward);
651 edges[cells[cell].edge_indices[local_edge]].vertex_indices[0] ==
656 edges[cells[cell].edge_indices[local_edge]].vertex_indices[1] ==
663 edges[cells[cell].edge_indices[local_edge]].orientation_status =
664 (dim == 2 ? Edge<dim>::forward : Edge<dim>::backward);
673 EdgeDeltaSet<dim> Delta_k;
674 EdgeDeltaSet<dim> Delta_k_minus_1;
675 Delta_k_minus_1.insert(cells[cell].edge_indices[local_edge]);
677 while (Delta_k_minus_1.begin() !=
678 Delta_k_minus_1.end())
682 for (
typename EdgeDeltaSet<dim>::const_iterator delta =
683 Delta_k_minus_1.begin();
684 delta != Delta_k_minus_1.end();
687 Assert(edges[*delta].orientation_status != Edge<dim>::not_oriented,
691 for (
typename AdjacentCells<dim>::const_iterator adjacent_cell =
692 edges[*delta].adjacent_cells.begin();
693 adjacent_cell != edges[*delta].adjacent_cells.end();
696 const unsigned int K = adjacent_cell->cell_index;
697 const unsigned int delta_is_edge_in_K =
698 adjacent_cell->edge_within_cell;
702 const unsigned int first_edge_vertex =
703 (edges[*delta].orientation_status == Edge<dim>::forward ?
704 edges[*delta].vertex_indices[0] :
705 edges[*delta].vertex_indices[1]);
706 const unsigned int first_edge_vertex_in_K =
709 delta_is_edge_in_K, 0)];
711 first_edge_vertex == first_edge_vertex_in_K ||
714 dim>::line_to_cell_vertices(delta_is_edge_in_K, 1)],
719 for (
unsigned int o_e = 0;
720 o_e < ParallelEdges<dim>::n_other_parallel_edges;
726 const unsigned int opposite_edge =
727 cells[K].edge_indices[ParallelEdges<
728 dim>::parallel_edges[delta_is_edge_in_K][o_e]];
729 const unsigned int first_opposite_edge_vertex =
730 cells[K].vertex_indices
732 ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K]
734 (first_edge_vertex == first_edge_vertex_in_K ? 0 :
741 const typename Edge<dim>::OrientationStatus
742 opposite_edge_orientation =
743 (edges[opposite_edge].vertex_indices[0] ==
744 first_opposite_edge_vertex ?
746 Edge<dim>::backward);
750 if (edges[opposite_edge].orientation_status ==
751 Edge<dim>::not_oriented)
755 edges[opposite_edge].orientation_status =
756 opposite_edge_orientation;
757 Delta_k.insert(opposite_edge);
772 Assert(edges[opposite_edge].orientation_status ==
773 opposite_edge_orientation,
778 if (edges[opposite_edge].orientation_status !=
779 opposite_edge_orientation)
792 Delta_k_minus_1 = Delta_k;
806 rotate_cell(
const std::vector<Cell<dim>> &cell_list,
807 const std::vector<Edge<dim>> &edge_list,
808 const unsigned int cell_index,
814 for (
unsigned int e = 0; e < GeometryInfo<dim>::lines_per_cell; ++
e)
816 Assert(edge_list[cell_list[cell_index].edge_indices[e]]
817 .orientation_status != Edge<dim>::not_oriented,
819 if (edge_list[cell_list[cell_index].edge_indices[e]]
820 .orientation_status == Edge<dim>::forward)
821 starting_vertex_of_edge[
e] =
822 edge_list[cell_list[cell_index].edge_indices[
e]].vertex_indices[0];
824 starting_vertex_of_edge[
e] =
825 edge_list[cell_list[cell_index].edge_indices[
e]].vertex_indices[1];
839 if ((starting_vertex_of_edge[0] == starting_vertex_of_edge[2]) ||
840 (starting_vertex_of_edge[0] == starting_vertex_of_edge[3]))
841 origin_vertex_of_cell = starting_vertex_of_edge[0];
842 else if ((starting_vertex_of_edge[1] ==
843 starting_vertex_of_edge[2]) ||
844 (starting_vertex_of_edge[1] == starting_vertex_of_edge[3]))
845 origin_vertex_of_cell = starting_vertex_of_edge[1];
857 for (origin_vertex_of_cell = 0;
858 origin_vertex_of_cell < GeometryInfo<dim>::vertices_per_cell;
859 ++origin_vertex_of_cell)
860 if (std::count(starting_vertex_of_edge,
861 starting_vertex_of_edge +
863 cell_list[cell_index]
864 .vertex_indices[origin_vertex_of_cell]) == dim)
888 while (raw_cells[cell_index].vertices[0] != origin_vertex_of_cell)
890 const unsigned int tmp = raw_cells[cell_index].vertices[0];
891 raw_cells[cell_index].vertices[0] =
892 raw_cells[cell_index].vertices[1];
893 raw_cells[cell_index].vertices[1] =
894 raw_cells[cell_index].vertices[3];
895 raw_cells[cell_index].vertices[3] =
896 raw_cells[cell_index].vertices[2];
897 raw_cells[cell_index].vertices[2] = tmp;
913 static const unsigned int cube_permutations[8][8] = {
914 {0, 1, 2, 3, 4, 5, 6, 7},
915 {1, 5, 3, 7, 0, 4, 2, 6},
916 {2, 6, 0, 4, 3, 7, 1, 5},
917 {3, 2, 1, 0, 7, 6, 5, 4},
918 {4, 0, 6, 2, 5, 1, 7, 3},
919 {5, 4, 7, 6, 1, 0, 3, 2},
920 {6, 7, 4, 5, 2, 3, 0, 1},
921 {7, 3, 5, 1, 6, 2, 4, 0}};
925 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell;
927 temp_vertex_indices[v] =
928 raw_cells[cell_index]
929 .vertices[cube_permutations[origin_vertex_of_cell][v]];
930 for (
unsigned int v = 0; v < GeometryInfo<dim>::vertices_per_cell;
932 raw_cells[cell_index].vertices[v] = temp_vertex_indices[v];
956 std::vector<Edge<dim>> edge_list = build_edges(cells);
957 std::vector<Cell<dim>> cell_list =
958 build_cells_and_connect_edges(cells, edge_list);
962 unsigned int next_cell_with_unoriented_edge = 0;
963 while ((next_cell_with_unoriented_edge = get_next_unoriented_cell(
964 cell_list, edge_list, next_cell_with_unoriented_edge)) !=
978 for (
unsigned int l = 0;
l < dim; ++
l)
979 if (edge_list[cell_list[next_cell_with_unoriented_edge]
980 .edge_indices[ParallelEdges<dim>::starter_edges[l]]]
981 .orientation_status == Edge<dim>::not_oriented)
982 orient_one_set_of_parallel_edges(
985 next_cell_with_unoriented_edge,
986 ParallelEdges<dim>::starter_edges[l]);
990 for (
unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++
l)
992 edge_list[cell_list[next_cell_with_unoriented_edge].edge_indices[l]]
993 .orientation_status != Edge<dim>::not_oriented,
1000 for (
unsigned int c = 0; c < cells.size(); ++c)
1001 rotate_cell(cell_list, edge_list, c, cells);
1024 void reorder_new_to_old_style(std::vector<
CellData<1>> &)
1028 void reorder_new_to_old_style(std::vector<
CellData<2>> &cells)
1030 for (
auto &cell : cells)
1031 std::swap(cell.vertices[2], cell.vertices[3]);
1035 void reorder_new_to_old_style(std::vector<
CellData<3>> &cells)
1038 for (
auto &cell : cells)
1040 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_cell; ++i)
1041 tmp[i] = cell.vertices[i];
1042 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_cell; ++i)
1051 void reorder_old_to_new_style(std::vector<
CellData<1>> &)
1055 void reorder_old_to_new_style(std::vector<
CellData<2>> &cells)
1058 reorder_new_to_old_style(cells);
1062 void reorder_old_to_new_style(std::vector<
CellData<3>> &cells)
1066 for (
auto &cell : cells)
1068 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_cell; ++i)
1069 tmp[i] = cell.vertices[i];
1070 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_cell; ++i)
1078 template <
int dim,
int spacedim>
1081 const bool use_new_style_ordering)
1083 Assert(cells.size() != 0,
1084 ExcMessage(
"List of elements to orient must have at least one cell"));
1091 if (use_new_style_ordering ==
false)
1092 reorder_old_to_new_style(cells);
1096 if (!is_consistent(cells))
1111 if (use_new_style_ordering ==
false)
1112 reorder_new_to_old_style(cells);
1152 const std::vector<
Point<2>> &all_vertices,
1156 unsigned int n_negative_cells = 0;
1157 for (
auto &cell : cells)
1162 for (
unsigned int i = 0; i < GeometryInfo<2>::vertices_per_cell; ++i)
1164 if (GridTools::cell_measure<2>(all_vertices, vertices_lex) < 0)
1167 std::swap(cell.vertices[1], cell.vertices[3]);
1175 for (
unsigned int i = 0; i < GeometryInfo<2>::vertices_per_cell; ++i)
1177 AssertThrow(GridTools::cell_measure<2>(all_vertices, vertices_lex) >
1189 AssertThrow(n_negative_cells == 0 || n_negative_cells == cells.size(),
1192 "This class assumes that either all cells have positive " 1193 "volume, or that all cells have been specified in an " 1194 "inverted vertex order so that their volume is negative. " 1195 "(In the latter case, this class automatically inverts " 1196 "every cell.) However, the mesh you have specified " 1197 "appears to have both cells with positive and cells with " 1198 "negative volume. You need to check your mesh which " 1199 "cells these are and how they got there.\n" 1200 "As a hint, of the total ") +
1203 " appear to have a negative volume."));
1222 const std::vector<
Point<3>> &all_vertices,
1226 unsigned int n_negative_cells = 0;
1227 for (
auto &cell : cells)
1232 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_cell; ++i)
1234 if (GridTools::cell_measure<3>(all_vertices, vertices_lex) < 0)
1238 for (
unsigned int i = 0; i < 4; ++i)
1239 std::swap(cell.vertices[i], cell.vertices[i + 4]);
1247 for (
unsigned int i = 0; i < GeometryInfo<3>::vertices_per_cell; ++i)
1249 AssertThrow(GridTools::cell_measure<3>(all_vertices, vertices_lex) >
1263 AssertThrow(n_negative_cells == 0 || n_negative_cells == cells.size(),
1264 ExcMessage(
"While sorting the cells that will be passed for " 1265 "creating a Triangulation object, deal.II found that " 1266 "some but not all cells have a negative volume. (If " 1267 "all cells had a negative volume, they would simply " 1268 "all have been inverted.) This usually happens in " 1269 "hand-generated meshes if one accidentally uses an " 1270 "incorrect convention for ordering the vertices in " 1271 "one or more cells; in that case, you may want to " 1272 "double check that you specified the vertex indices " 1273 "in their correct order. If you are reading a mesh " 1274 "that was created by a mesh generator, then this " 1275 "exception indicates that some of the cells created " 1276 "are so badly distorted that their volume becomes " 1277 "negative; this commonly occurs at complex geometric " 1278 "features, and you may see if the problem can be " 1279 "fixed by playing with the parameters that control " 1280 "mesh properties in your mesh generator, such as " 1281 "the number of cells, the mesh density, etc."));
1294 DEAL_II_NAMESPACE_CLOSE
static void reorder_cells(std::vector< CellData< dim >> &original_cells, const bool use_new_style_ordering=false)
static const unsigned int invalid_unsigned_int
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
static unsigned int line_to_cell_vertices(const unsigned int line, const unsigned int vertex)
#define AssertThrow(cond, exc)
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
void swap(Vector< Number > &u, Vector< Number > &v)
__global__ void set(Number *val, const Number s, const size_type N)
static ::ExceptionBase & ExcNotImplemented()
static void invert_all_cells_of_negative_grid(const std::vector< Point< spacedim >> &all_vertices, std::vector< CellData< dim >> &original_cells)
static ::ExceptionBase & ExcMeshNotOrientable()
Tensor< 2, dim, Number > l(const Tensor< 2, dim, Number > &F, const Tensor< 2, dim, Number > &dF_dt)
unsigned int vertices[GeometryInfo< structdim >::vertices_per_cell]
static ::ExceptionBase & ExcInternalError()