17 #include <deal.II/grid/grid_reordering.h> 18 #include <deal.II/grid/grid_tools.h> 19 #include <deal.II/base/utilities.h> 20 #include <deal.II/base/timer.h> 28 DEAL_II_NAMESPACE_OPEN
43 CheapEdge (
const unsigned int v0,
44 const unsigned int v1)
53 bool operator < (
const CheapEdge &e)
const 55 return ((v0 <
e.v0) || ((v0 ==
e.v0) && (v1 <
e.v1)));
62 const unsigned int v0, v1;
78 std::set<CheapEdge> edges;
80 for (
typename std::vector<
CellData<dim> >::const_iterator c = cells.begin();
81 c != cells.end(); ++c)
88 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
92 if (edges.find (reverse_edge) != edges.end())
99 edges.insert (correct_edge);
122 static const unsigned int starter_edges[dim];
128 static const unsigned int n_other_parallel_edges = (1<<(dim-1)) - 1;
133 const unsigned int ParallelEdges<2>::starter_edges[2] = { 0, 2 };
136 const unsigned int ParallelEdges<2>::parallel_edges[4][1] = { {1}, {0}, {3}, {2} };
139 const unsigned int ParallelEdges<3>::starter_edges[3] = { 0, 2, 8 };
142 const unsigned int ParallelEdges<3>::parallel_edges[12][3] = { {1, 4, 5},
168 cell_index (
numbers::invalid_unsigned_int),
169 edge_within_cell (
numbers::invalid_unsigned_int)
175 AdjacentCell (
const unsigned int cell_index,
176 const unsigned int edge_within_cell)
178 cell_index (cell_index),
179 edge_within_cell (edge_within_cell)
183 unsigned int cell_index;
184 unsigned int edge_within_cell;
189 template <
int dim>
class AdjacentCells;
197 class AdjacentCells<2>
204 typedef const AdjacentCell *const_iterator;
214 void push_back (
const AdjacentCell &adjacent_cell)
217 adjacent_cells[0] = adjacent_cell;
222 adjacent_cells[1] = adjacent_cell;
231 const_iterator begin ()
const 233 return adjacent_cells;
242 const_iterator end ()
const 248 return adjacent_cells;
250 return adjacent_cells + 1;
252 return adjacent_cells + 2;
262 AdjacentCell adjacent_cells[2];
275 class AdjacentCells<3> :
public std::vector<AdjacentCell>
298 const unsigned int edge_number)
300 orientation_status (not_oriented)
309 if (vertex_indices[0] > vertex_indices[1])
310 std::swap (vertex_indices[0], vertex_indices[1]);
317 bool operator< (const Edge<dim> &
e)
const 319 return ((vertex_indices[0] <
e.vertex_indices[0])
321 ((vertex_indices[0] ==
e.vertex_indices[0]) && (vertex_indices[1] <
e.vertex_indices[1])));
327 bool operator== (
const Edge<dim> &e)
const 329 return ((vertex_indices[0] ==
e.vertex_indices[0])
331 (vertex_indices[1] ==
e.vertex_indices[1]));
338 unsigned int vertex_indices[2];
344 enum OrientationStatus
351 OrientationStatus orientation_status;
357 AdjacentCells<dim> adjacent_cells;
375 const std::vector<Edge<dim> > &edge_list)
377 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
382 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
384 const Edge<dim>
e (c, l);
385 edge_indices[
l] = (std::lower_bound (edge_list.begin(), edge_list.end(),
e)
407 template <
int dim>
class EdgeDeltaSet;
419 class EdgeDeltaSet<2>
425 typedef const unsigned int *const_iterator;
449 void insert (
const unsigned int edge_index)
452 edge_indices[0] = edge_index;
457 edge_indices[1] = edge_index;
465 const_iterator begin ()
const 474 const_iterator end ()
const 482 return edge_indices + 1;
484 return edge_indices + 2;
491 unsigned int edge_indices[2];
508 class EdgeDeltaSet<3> :
public std::set<unsigned int>
519 std::vector<Edge<dim> >
526 std::vector<Edge<dim> > edge_list;
528 for (
unsigned int i=0; i<cells.size(); ++i)
529 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
530 edge_list.emplace_back (cells[i], l);
533 std::sort (edge_list.begin(), edge_list.end());
534 edge_list.erase(std::unique(edge_list.begin(),edge_list.end()),
547 std::vector<Cell<dim> >
548 build_cells_and_connect_edges (
const std::vector<
CellData<dim> > &cells,
549 std::vector<Edge<dim> > &edges)
551 std::vector<Cell<dim> > cell_list;
552 cell_list.reserve(cells.size());
553 for (
unsigned int i=0; i<cells.size(); ++i)
557 cell_list.emplace_back (cells[i], edges);
561 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
562 edges[cell_list.back().edge_indices[
l]].adjacent_cells.push_back (AdjacentCell (i, l));
577 get_next_unoriented_cell(
const std::vector<Cell<dim> > &cells,
578 const std::vector<Edge<dim> > &edges,
579 const unsigned int current_cell)
581 for (
unsigned int c=current_cell; c<cells.size(); ++c)
582 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
583 if (edges[cells[c].edge_indices[l]].orientation_status == Edge<dim>::not_oriented)
598 orient_one_set_of_parallel_edges (
const std::vector<Cell<dim> > &cells,
599 std::vector<Edge<dim> > &edges,
600 const unsigned int cell,
601 const unsigned int local_edge)
623 if (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
628 edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
629 Edge<dim>::backward :
633 Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
637 Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[1]
644 edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
646 Edge<dim>::backward);
655 EdgeDeltaSet<dim> Delta_k;
656 EdgeDeltaSet<dim> Delta_k_minus_1;
657 Delta_k_minus_1.insert (cells[cell].edge_indices[local_edge]);
659 while (Delta_k_minus_1.begin() != Delta_k_minus_1.end())
663 for (
typename EdgeDeltaSet<dim>::const_iterator delta = Delta_k_minus_1.begin();
664 delta != Delta_k_minus_1.end(); ++delta)
666 Assert (edges[*delta].orientation_status != Edge<dim>::not_oriented,
670 for (
typename AdjacentCells<dim>::const_iterator
671 adjacent_cell = edges[*delta].adjacent_cells.begin();
672 adjacent_cell != edges[*delta].adjacent_cells.end(); ++adjacent_cell)
674 const unsigned int K = adjacent_cell->cell_index;
675 const unsigned int delta_is_edge_in_K = adjacent_cell->edge_within_cell;
679 const unsigned int first_edge_vertex
680 = (edges[*delta].orientation_status == Edge<dim>::forward
682 edges[*delta].vertex_indices[0]
684 edges[*delta].vertex_indices[1]);
685 const unsigned int first_edge_vertex_in_K
687 Assert (first_edge_vertex == first_edge_vertex_in_K
694 for (
unsigned int o_e=0; o_e<ParallelEdges<dim>::n_other_parallel_edges; ++o_e)
699 const unsigned int opposite_edge
700 = cells[K].edge_indices[ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K][o_e]];
701 const unsigned int first_opposite_edge_vertex
703 ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K][o_e],
704 (first_edge_vertex == first_edge_vertex_in_K
714 const typename Edge<dim>::OrientationStatus opposite_edge_orientation
715 = (edges[opposite_edge].vertex_indices[0]
717 first_opposite_edge_vertex
721 Edge<dim>::backward);
725 if (edges[opposite_edge].orientation_status == Edge<dim>::not_oriented)
729 edges[opposite_edge].orientation_status = opposite_edge_orientation;
730 Delta_k.insert (opposite_edge);
744 Assert (edges[opposite_edge].orientation_status == opposite_edge_orientation,
749 if (edges[opposite_edge].orientation_status != opposite_edge_orientation)
762 Delta_k_minus_1 = Delta_k;
776 rotate_cell (
const std::vector<Cell<dim> > &cell_list,
777 const std::vector<Edge<dim> > &edge_list,
778 const unsigned int cell_index,
784 for (
unsigned int e=0; e<GeometryInfo<dim>::lines_per_cell; ++
e)
786 Assert (edge_list[cell_list[cell_index].edge_indices[e]].orientation_status
787 != Edge<dim>::not_oriented,
789 if (edge_list[cell_list[cell_index].edge_indices[e]].orientation_status == Edge<dim>::forward)
790 starting_vertex_of_edge[
e] = edge_list[cell_list[cell_index].edge_indices[
e]].vertex_indices[0];
792 starting_vertex_of_edge[
e] = edge_list[cell_list[cell_index].edge_indices[
e]].vertex_indices[1];
806 if ((starting_vertex_of_edge[0] == starting_vertex_of_edge[2])
808 (starting_vertex_of_edge[0] == starting_vertex_of_edge[3]))
809 origin_vertex_of_cell = starting_vertex_of_edge[0];
810 else if ((starting_vertex_of_edge[1] == starting_vertex_of_edge[2])
812 (starting_vertex_of_edge[1] == starting_vertex_of_edge[3]))
813 origin_vertex_of_cell = starting_vertex_of_edge[1];
825 for (origin_vertex_of_cell=0;
826 origin_vertex_of_cell<GeometryInfo<dim>::vertices_per_cell;
827 ++origin_vertex_of_cell)
828 if (std::count (starting_vertex_of_edge,
830 cell_list[cell_index].vertex_indices[origin_vertex_of_cell])
855 while (raw_cells[cell_index].vertices[0] != origin_vertex_of_cell)
857 const unsigned int tmp = raw_cells[cell_index].vertices[0];
858 raw_cells[cell_index].vertices[0] = raw_cells[cell_index].vertices[1];
859 raw_cells[cell_index].vertices[1] = raw_cells[cell_index].vertices[3];
860 raw_cells[cell_index].vertices[3] = raw_cells[cell_index].vertices[2];
861 raw_cells[cell_index].vertices[2] = tmp;
877 static const unsigned int cube_permutations[8][8] =
890 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
891 temp_vertex_indices[v]
892 = raw_cells[cell_index].vertices[cube_permutations[origin_vertex_of_cell][v]];
893 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
894 raw_cells[cell_index].vertices[v] = temp_vertex_indices[v];
917 std::vector<Edge<dim> > edge_list = build_edges(cells);
918 std::vector<Cell<dim> > cell_list = build_cells_and_connect_edges(cells, edge_list);
922 unsigned int next_cell_with_unoriented_edge = 0;
923 while ((next_cell_with_unoriented_edge = get_next_unoriented_cell(cell_list,
925 next_cell_with_unoriented_edge)) !=
939 for (
unsigned int l=0;
l<dim; ++
l)
940 if (edge_list[cell_list[next_cell_with_unoriented_edge].edge_indices[ParallelEdges<dim>::starter_edges[l]]].orientation_status
941 == Edge<dim>::not_oriented)
942 orient_one_set_of_parallel_edges (cell_list,
944 next_cell_with_unoriented_edge,
945 ParallelEdges<dim>::starter_edges[l]);
949 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
950 Assert (edge_list[cell_list[next_cell_with_unoriented_edge].edge_indices[l]].orientation_status
951 != Edge<dim>::not_oriented,
958 for (
unsigned int c=0; c<cells.size(); ++c)
959 rotate_cell (cell_list, edge_list, c, cells);
983 reorder_new_to_old_style (std::vector<
CellData<1> > &)
988 reorder_new_to_old_style (std::vector<
CellData<2> > &cells)
990 for (
unsigned int cell=0; cell<cells.size(); ++cell)
991 std::swap(cells[cell].vertices[2], cells[cell].vertices[3]);
996 reorder_new_to_old_style (std::vector<
CellData<3> > &cells)
999 for (
unsigned int cell=0; cell<cells.size(); ++cell)
1001 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1002 tmp[i] = cells[cell].vertices[i];
1003 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1013 reorder_old_to_new_style (std::vector<
CellData<1> > &)
1018 reorder_old_to_new_style (std::vector<
CellData<2> > &cells)
1021 reorder_new_to_old_style(cells);
1026 reorder_old_to_new_style (std::vector<
CellData<3> > &cells)
1030 for (
unsigned int cell=0; cell<cells.size(); ++cell)
1032 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1033 tmp[i] = cells[cell].vertices[i];
1034 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1042 template <
int dim,
int spacedim>
1045 const bool use_new_style_ordering)
1047 Assert (cells.size() != 0,
1048 ExcMessage(
"List of elements to orient must have at least one cell"));
1055 if (use_new_style_ordering ==
false)
1056 reorder_old_to_new_style(cells);
1060 if (!is_consistent (cells))
1075 if (use_new_style_ordering ==
false)
1076 reorder_new_to_old_style(cells);
1116 unsigned int n_negative_cells=0;
1117 for (
unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1122 for (
unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1124 if (GridTools::cell_measure<2>(all_vertices, vertices_lex) < 0)
1127 std::swap(cells[cell_no].vertices[1], cells[cell_no].vertices[3]);
1135 for (
unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1137 AssertThrow(GridTools::cell_measure<2>(all_vertices, vertices_lex) > 0,
1148 AssertThrow(n_negative_cells==0 || n_negative_cells==cells.size(),
1149 ExcMessage(std::string(
"This class assumes that either all cells have positive " 1150 "volume, or that all cells have been specified in an " 1151 "inverted vertex order so that their volume is negative. " 1152 "(In the latter case, this class automatically inverts " 1153 "every cell.) However, the mesh you have specified " 1154 "appears to have both cells with positive and cells with " 1155 "negative volume. You need to check your mesh which " 1156 "cells these are and how they got there.\n" 1157 "As a hint, of the total ")
1159 +
" cells in the mesh, " 1161 +
" appear to have a negative volume."));
1182 unsigned int n_negative_cells=0;
1183 for (
unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1188 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1190 if (GridTools::cell_measure<3>(all_vertices, vertices_lex) < 0)
1194 for (
unsigned int i=0; i<4; ++i)
1195 std::swap(cells[cell_no].vertices[i], cells[cell_no].vertices[i+4]);
1203 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1205 AssertThrow(GridTools::cell_measure<3>(all_vertices, vertices_lex) > 0,
1218 AssertThrow(n_negative_cells==0 || n_negative_cells==cells.size(),
1219 ExcMessage(
"While sorting the cells that will be passed for " 1220 "creating a Triangulation object, deal.II found that " 1221 "some but not all cells have a negative volume. (If " 1222 "all cells had a negative volume, they would simply " 1223 "all have been inverted.) This usually happens in " 1224 "hand-generated meshes if one accidentally uses an " 1225 "incorrect convention for ordering the vertices in " 1226 "one or more cells; in that case, you may want to " 1227 "double check that you specified the vertex indices " 1228 "in their correct order. If you are reading a mesh " 1229 "that was created by a mesh generator, then this " 1230 "exception indicates that some of the cells created " 1231 "are so badly distorted that their volume becomes " 1232 "negative; this commonly occurs at complex geometric " 1233 "features, and you may see if the problem can be " 1234 "fixed by playing with the parameters that control " 1235 "mesh properties in your mesh generator, such as " 1236 "the number of cells, the mesh density, etc."));
1249 DEAL_II_NAMESPACE_CLOSE
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)
static void invert_all_cells_of_negative_grid(const std::vector< Point< spacedim > > &all_vertices, std::vector< CellData< dim > > &original_cells)
static void reorder_cells(std::vector< CellData< dim > > &original_cells, const bool use_new_style_ordering=false)
void swap(Vector< Number > &u, Vector< Number > &v)
static ::ExceptionBase & ExcNotImplemented()
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()