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/std_cxx11/bind.h> 21 #include <deal.II/base/timer.h> 29 DEAL_II_NAMESPACE_OPEN
44 CheapEdge (
const unsigned int v0,
45 const unsigned int v1)
54 bool 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();
82 c != cells.end(); ++c)
89 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
93 if (edges.find (reverse_edge) != edges.end())
100 edges.insert (correct_edge);
123 static const unsigned int starter_edges[dim];
129 static const unsigned int n_other_parallel_edges = (1<<(dim-1)) - 1;
134 const unsigned int ParallelEdges<2>::starter_edges[2] = { 0, 2 };
137 const unsigned int ParallelEdges<2>::parallel_edges[4][1] = { {1}, {0}, {3}, {2} };
140 const unsigned int ParallelEdges<3>::starter_edges[3] = { 0, 2, 8 };
143 const unsigned int ParallelEdges<3>::parallel_edges[12][3] = { {1, 4, 5},
169 cell_index (
numbers::invalid_unsigned_int),
170 edge_within_cell (
numbers::invalid_unsigned_int)
176 AdjacentCell (
const unsigned int cell_index,
177 const unsigned int edge_within_cell)
179 cell_index (cell_index),
180 edge_within_cell (edge_within_cell)
184 unsigned int cell_index;
185 unsigned int edge_within_cell;
190 template <
int dim>
class AdjacentCells;
198 class AdjacentCells<2>
205 typedef const AdjacentCell *const_iterator;
215 void push_back (
const AdjacentCell &adjacent_cell)
218 adjacent_cells[0] = adjacent_cell;
223 adjacent_cells[1] = adjacent_cell;
232 const_iterator begin ()
const 234 return &adjacent_cells[0];
243 const_iterator end ()
const 249 return &adjacent_cells[0];
251 return &adjacent_cells[0]+1;
253 return &adjacent_cells[0]+2;
263 AdjacentCell adjacent_cells[2];
276 class AdjacentCells<3> :
public std::vector<AdjacentCell>
298 orientation_status (not_oriented)
300 for (
unsigned int i=0; i<2; ++i)
310 const unsigned int edge_number)
312 orientation_status (not_oriented)
321 if (vertex_indices[0] > vertex_indices[1])
322 std::swap (vertex_indices[0], vertex_indices[1]);
329 bool operator< (const Edge<dim> &
e)
const 331 return ((vertex_indices[0] <
e.vertex_indices[0])
333 ((vertex_indices[0] ==
e.vertex_indices[0]) && (vertex_indices[1] <
e.vertex_indices[1])));
339 bool operator== (
const Edge<dim> &e)
const 341 return ((vertex_indices[0] ==
e.vertex_indices[0])
343 (vertex_indices[1] ==
e.vertex_indices[1]));
350 unsigned int vertex_indices[2];
356 enum OrientationStatus
363 OrientationStatus orientation_status;
369 AdjacentCells<dim> adjacent_cells;
386 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
388 for (
unsigned int i=0; i<GeometryInfo<dim>::lines_per_cell; ++i)
398 const std::vector<Edge<dim> > &edge_list)
400 for (
unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
405 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
407 const Edge<dim>
e (c, l);
408 edge_indices[
l] = (std::lower_bound (edge_list.begin(), edge_list.end(),
e)
430 template <
int dim>
class EdgeDeltaSet;
442 class EdgeDeltaSet<2>
448 typedef const unsigned int *const_iterator;
472 void insert (
const unsigned int edge_index)
475 edge_indices[0] = edge_index;
480 edge_indices[1] = edge_index;
488 const_iterator begin ()
const 490 return &edge_indices[0];
497 const_iterator end ()
const 503 return &edge_indices[0];
505 return &edge_indices[0]+1;
507 return &edge_indices[0]+2;
514 unsigned int edge_indices[2];
531 class EdgeDeltaSet<3> :
public std::set<unsigned int>
542 std::vector<Edge<dim> >
549 std::vector<Edge<dim> > edge_list;
551 for (
unsigned int i=0; i<cells.size(); ++i)
552 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
553 edge_list.push_back (Edge<dim>(cells[i], l));
556 std::sort (edge_list.begin(), edge_list.end());
557 edge_list.erase(std::unique(edge_list.begin(),edge_list.end()),
570 std::vector<Cell<dim> >
571 build_cells_and_connect_edges (
const std::vector<
CellData<dim> > &cells,
572 std::vector<Edge<dim> > &edges)
574 std::vector<Cell<dim> > cell_list;
575 cell_list.reserve(cells.size());
576 for (
unsigned int i=0; i<cells.size(); ++i)
580 cell_list.push_back (Cell<dim>(cells[i], edges));
584 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
585 edges[cell_list.back().edge_indices[
l]].adjacent_cells.push_back (AdjacentCell (i, l));
600 get_next_unoriented_cell(
const std::vector<Cell<dim> > &cells,
601 const std::vector<Edge<dim> > &edges,
602 const unsigned int current_cell)
604 for (
unsigned int c=current_cell; c<cells.size(); ++c)
605 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
606 if (edges[cells[c].edge_indices[l]].orientation_status == Edge<dim>::not_oriented)
621 orient_one_set_of_parallel_edges (
const std::vector<Cell<dim> > &cells,
622 std::vector<Edge<dim> > &edges,
623 const unsigned int cell,
624 const unsigned int local_edge)
646 if (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
651 edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
652 Edge<dim>::backward :
656 Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
660 Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[1]
667 edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
669 Edge<dim>::backward);
678 EdgeDeltaSet<dim> Delta_k;
679 EdgeDeltaSet<dim> Delta_k_minus_1;
680 Delta_k_minus_1.insert (cells[cell].edge_indices[local_edge]);
682 while (Delta_k_minus_1.begin() != Delta_k_minus_1.end())
686 for (
typename EdgeDeltaSet<dim>::const_iterator delta = Delta_k_minus_1.begin();
687 delta != Delta_k_minus_1.end(); ++delta)
689 Assert (edges[*delta].orientation_status != Edge<dim>::not_oriented,
693 for (
typename AdjacentCells<dim>::const_iterator
694 adjacent_cell = edges[*delta].adjacent_cells.begin();
695 adjacent_cell != edges[*delta].adjacent_cells.end(); ++adjacent_cell)
697 const unsigned int K = adjacent_cell->cell_index;
698 const unsigned int delta_is_edge_in_K = adjacent_cell->edge_within_cell;
702 const unsigned int first_edge_vertex
703 = (edges[*delta].orientation_status == Edge<dim>::forward
705 edges[*delta].vertex_indices[0]
707 edges[*delta].vertex_indices[1]);
708 const unsigned int first_edge_vertex_in_K
710 Assert (first_edge_vertex == first_edge_vertex_in_K
717 for (
unsigned int o_e=0; o_e<ParallelEdges<dim>::n_other_parallel_edges; ++o_e)
722 const unsigned int opposite_edge
723 = cells[K].edge_indices[ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K][o_e]];
724 const unsigned int first_opposite_edge_vertex
726 ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K][o_e],
727 (first_edge_vertex == first_edge_vertex_in_K
737 const typename Edge<dim>::OrientationStatus opposite_edge_orientation
738 = (edges[opposite_edge].vertex_indices[0]
740 first_opposite_edge_vertex
744 Edge<dim>::backward);
748 if (edges[opposite_edge].orientation_status == Edge<dim>::not_oriented)
752 edges[opposite_edge].orientation_status = opposite_edge_orientation;
753 Delta_k.insert (opposite_edge);
767 Assert (edges[opposite_edge].orientation_status == opposite_edge_orientation,
772 if (edges[opposite_edge].orientation_status != opposite_edge_orientation)
785 Delta_k_minus_1 = Delta_k;
799 rotate_cell (
const std::vector<Cell<dim> > &cell_list,
800 const std::vector<Edge<dim> > &edge_list,
801 const unsigned int cell_index,
807 for (
unsigned int e=0; e<GeometryInfo<dim>::lines_per_cell; ++
e)
809 Assert (edge_list[cell_list[cell_index].edge_indices[e]].orientation_status
810 != Edge<dim>::not_oriented,
812 if (edge_list[cell_list[cell_index].edge_indices[e]].orientation_status == Edge<dim>::forward)
813 starting_vertex_of_edge[
e] = edge_list[cell_list[cell_index].edge_indices[
e]].vertex_indices[0];
815 starting_vertex_of_edge[
e] = edge_list[cell_list[cell_index].edge_indices[
e]].vertex_indices[1];
829 if ((starting_vertex_of_edge[0] == starting_vertex_of_edge[2])
831 (starting_vertex_of_edge[0] == starting_vertex_of_edge[3]))
832 origin_vertex_of_cell = starting_vertex_of_edge[0];
833 else if ((starting_vertex_of_edge[1] == starting_vertex_of_edge[2])
835 (starting_vertex_of_edge[1] == starting_vertex_of_edge[3]))
836 origin_vertex_of_cell = starting_vertex_of_edge[1];
848 for (origin_vertex_of_cell=0;
849 origin_vertex_of_cell<GeometryInfo<dim>::vertices_per_cell;
850 ++origin_vertex_of_cell)
851 if (std::count (&starting_vertex_of_edge[0],
853 cell_list[cell_index].vertex_indices[origin_vertex_of_cell])
878 while (raw_cells[cell_index].vertices[0] != origin_vertex_of_cell)
880 const unsigned int tmp = raw_cells[cell_index].vertices[0];
881 raw_cells[cell_index].vertices[0] = raw_cells[cell_index].vertices[1];
882 raw_cells[cell_index].vertices[1] = raw_cells[cell_index].vertices[3];
883 raw_cells[cell_index].vertices[3] = raw_cells[cell_index].vertices[2];
884 raw_cells[cell_index].vertices[2] = tmp;
900 static const unsigned int cube_permutations[8][8] =
913 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
914 temp_vertex_indices[v]
915 = raw_cells[cell_index].vertices[cube_permutations[origin_vertex_of_cell][v]];
916 for (
unsigned int v=0; v<GeometryInfo<dim>::vertices_per_cell; ++v)
917 raw_cells[cell_index].vertices[v] = temp_vertex_indices[v];
940 std::vector<Edge<dim> > edge_list = build_edges(cells);
941 std::vector<Cell<dim> > cell_list = build_cells_and_connect_edges(cells, edge_list);
945 unsigned int next_cell_with_unoriented_edge = 0;
946 while ((next_cell_with_unoriented_edge = get_next_unoriented_cell(cell_list,
948 next_cell_with_unoriented_edge)) !=
962 for (
unsigned int l=0;
l<dim; ++
l)
963 if (edge_list[cell_list[next_cell_with_unoriented_edge].edge_indices[ParallelEdges<dim>::starter_edges[l]]].orientation_status
964 == Edge<dim>::not_oriented)
965 orient_one_set_of_parallel_edges (cell_list,
967 next_cell_with_unoriented_edge,
968 ParallelEdges<dim>::starter_edges[l]);
972 for (
unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++
l)
973 Assert (edge_list[cell_list[next_cell_with_unoriented_edge].edge_indices[l]].orientation_status
974 != Edge<dim>::not_oriented,
981 for (
unsigned int c=0; c<cells.size(); ++c)
982 rotate_cell (cell_list, edge_list, c, cells);
1006 reorder_new_to_old_style (std::vector<
CellData<1> > &)
1011 reorder_new_to_old_style (std::vector<
CellData<2> > &cells)
1013 for (
unsigned int cell=0; cell<cells.size(); ++cell)
1014 std::swap(cells[cell].vertices[2], cells[cell].vertices[3]);
1019 reorder_new_to_old_style (std::vector<
CellData<3> > &cells)
1022 for (
unsigned int cell=0; cell<cells.size(); ++cell)
1024 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1025 tmp[i] = cells[cell].vertices[i];
1026 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1036 reorder_old_to_new_style (std::vector<
CellData<1> > &)
1041 reorder_old_to_new_style (std::vector<
CellData<2> > &cells)
1044 reorder_new_to_old_style(cells);
1049 reorder_old_to_new_style (std::vector<
CellData<3> > &cells)
1053 for (
unsigned int cell=0; cell<cells.size(); ++cell)
1055 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1056 tmp[i] = cells[cell].vertices[i];
1057 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1065 template <
int dim,
int spacedim>
1068 const bool use_new_style_ordering)
1070 Assert (cells.size() != 0,
1071 ExcMessage(
"List of elements to orient must have at least one cell"));
1078 if (use_new_style_ordering ==
false)
1079 reorder_old_to_new_style(cells);
1083 if (!is_consistent (cells))
1098 if (use_new_style_ordering ==
false)
1099 reorder_new_to_old_style(cells);
1139 unsigned int n_negative_cells=0;
1140 for (
unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1145 for (
unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1147 if (GridTools::cell_measure<2>(all_vertices, vertices_lex) < 0)
1150 std::swap(cells[cell_no].vertices[1], cells[cell_no].vertices[3]);
1158 for (
unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1160 AssertThrow(GridTools::cell_measure<2>(all_vertices, vertices_lex) > 0,
1171 AssertThrow(n_negative_cells==0 || n_negative_cells==cells.size(),
1172 ExcMessage(std::string(
"This class assumes that either all cells have positive " 1173 "volume, or that all cells have been specified in an " 1174 "inverted vertex order so that their volume is negative. " 1175 "(In the latter case, this class automatically inverts " 1176 "every cell.) However, the mesh you have specified " 1177 "appears to have both cells with positive and cells with " 1178 "negative volume. You need to check your mesh which " 1179 "cells these are and how they got there.\n" 1180 "As a hint, of the total ")
1182 +
" cells in the mesh, " 1184 +
" appear to have a negative volume."));
1205 unsigned int n_negative_cells=0;
1206 for (
unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1211 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1213 if (GridTools::cell_measure<3>(all_vertices, vertices_lex) < 0)
1217 for (
unsigned int i=0; i<4; ++i)
1218 std::swap(cells[cell_no].vertices[i], cells[cell_no].vertices[i+4]);
1226 for (
unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1228 AssertThrow(GridTools::cell_measure<3>(all_vertices, vertices_lex) > 0,
1241 AssertThrow(n_negative_cells==0 || n_negative_cells==cells.size(),
1242 ExcMessage(
"While sorting the cells that will be passed for " 1243 "creating a Triangulation object, deal.II found that " 1244 "some but not all cells have a negative volume. (If " 1245 "all cells had a negative volume, they would simply " 1246 "all have been inverted.) This usually happens in " 1247 "hand-generated meshes if one accidentally uses an " 1248 "incorrect convention for ordering the vertices in " 1249 "one or more cells; in that case, you may want to " 1250 "double check that you specified the vertex indices " 1251 "in their correct order. If you are reading a mesh " 1252 "that was created by a mesh generator, then this " 1253 "exception indicates that some of the cells created " 1254 "are so badly distorted that their volume becomes " 1255 "negative; this commonly occurs at complex geometric " 1256 "features, and you may see if the problem can be " 1257 "fixed by playing with the parameters that control " 1258 "mesh properties in your mesh generator, such as " 1259 "the number of cells, the mesh density, etc."));
1272 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)
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()