Reference documentation for deal.II version 9.0.0
grid_reordering.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2017 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 
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>
21 
22 #include <algorithm>
23 #include <set>
24 #include <iostream>
25 #include <fstream>
26 #include <functional>
27 
28 DEAL_II_NAMESPACE_OPEN
29 
30 
31 namespace
32 {
38  struct CheapEdge
39  {
43  CheapEdge (const unsigned int v0,
44  const unsigned int v1)
45  :
46  v0(v0), v1(v1)
47  {}
48 
53  bool operator < (const CheapEdge &e) const
54  {
55  return ((v0 < e.v0) || ((v0 == e.v0) && (v1 < e.v1)));
56  }
57 
58  private:
62  const unsigned int v0, v1;
63  };
64 
65 
74  template <int dim>
75  bool
76  is_consistent (const std::vector<CellData<dim> > &cells)
77  {
78  std::set<CheapEdge> edges;
79 
80  for (typename std::vector<CellData<dim> >::const_iterator c = cells.begin();
81  c != cells.end(); ++c)
82  {
83  // construct the edges in reverse order. for each of them,
84  // ensure that the reverse edge is not yet in the list of
85  // edges (return false if the reverse edge already *is* in
86  // the list) and then add the actual edge to it; std::set
87  // eliminates duplicates automatically
88  for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
89  {
90  const CheapEdge reverse_edge (c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 1)],
91  c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 0)]);
92  if (edges.find (reverse_edge) != edges.end())
93  return false;
94 
95 
96  // ok, not. insert edge in correct order
97  const CheapEdge correct_edge (c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 0)],
98  c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 1)]);
99  edges.insert (correct_edge);
100  }
101  }
102 
103  // no conflicts found, so return true
104  return true;
105  }
106 
107 
114  template <int dim>
115  struct ParallelEdges
116  {
122  static const unsigned int starter_edges[dim];
123 
128  static const unsigned int n_other_parallel_edges = (1<<(dim-1)) - 1;
129  static const unsigned int parallel_edges[GeometryInfo<dim>::lines_per_cell][n_other_parallel_edges];
130  };
131 
132  template <>
133  const unsigned int ParallelEdges<2>::starter_edges[2] = { 0, 2 };
134 
135  template <>
136  const unsigned int ParallelEdges<2>::parallel_edges[4][1] = { {1}, {0}, {3}, {2} };
137 
138  template <>
139  const unsigned int ParallelEdges<3>::starter_edges[3] = { 0, 2, 8 };
140 
141  template <>
142  const unsigned int ParallelEdges<3>::parallel_edges[12][3] = { {1, 4, 5}, // line 0
143  {0, 4, 5}, // line 1
144  {3, 6, 7}, // line 2
145  {2, 6, 7}, // line 3
146  {0, 1, 5}, // line 4
147  {0, 1, 4}, // line 5
148  {2, 3, 7}, // line 6
149  {2, 3, 6}, // line 7
150  {9, 10, 11}, // line 8
151  {8, 10, 11}, // line 9
152  {8, 9, 11}, // line 10
153  {8, 9, 10} // line 11
154  };
155 
156 
161  struct AdjacentCell
162  {
166  AdjacentCell ()
167  :
168  cell_index (numbers::invalid_unsigned_int),
169  edge_within_cell (numbers::invalid_unsigned_int)
170  {}
171 
175  AdjacentCell (const unsigned int cell_index,
176  const unsigned int edge_within_cell)
177  :
178  cell_index (cell_index),
179  edge_within_cell (edge_within_cell)
180  {}
181 
182 
183  unsigned int cell_index;
184  unsigned int edge_within_cell;
185  };
186 
187 
188 
189  template <int dim> class AdjacentCells;
190 
196  template <>
197  class AdjacentCells<2>
198  {
199  public:
204  typedef const AdjacentCell *const_iterator;
205 
214  void push_back (const AdjacentCell &adjacent_cell)
215  {
216  if (adjacent_cells[0].cell_index == numbers::invalid_unsigned_int)
217  adjacent_cells[0] = adjacent_cell;
218  else
219  {
220  Assert (adjacent_cells[1].cell_index == numbers::invalid_unsigned_int,
221  ExcInternalError());
222  adjacent_cells[1] = adjacent_cell;
223  }
224  }
225 
226 
231  const_iterator begin () const
232  {
233  return adjacent_cells;
234  }
235 
236 
242  const_iterator end () const
243  {
244  // check whether the current object stores zero, one, or two
245  // adjacent cells, and use this to point to the element past the
246  // last valid one
247  if (adjacent_cells[0].cell_index == numbers::invalid_unsigned_int)
248  return adjacent_cells;
249  else if (adjacent_cells[1].cell_index == numbers::invalid_unsigned_int)
250  return adjacent_cells + 1;
251  else
252  return adjacent_cells + 2;
253  }
254 
255  private:
262  AdjacentCell adjacent_cells[2];
263  };
264 
265 
266 
274  template <>
275  class AdjacentCells<3> : public std::vector<AdjacentCell>
276  {};
277 
278 
288  template <int dim>
289  class Edge
290  {
291  public:
297  Edge (const CellData<dim> &cell,
298  const unsigned int edge_number)
299  :
300  orientation_status (not_oriented)
301  {
303 
304  // copy vertices for this particular line
305  vertex_indices[0] = cell.vertices[GeometryInfo<dim>::line_to_cell_vertices(edge_number, 0)];
306  vertex_indices[1] = cell.vertices[GeometryInfo<dim>::line_to_cell_vertices(edge_number, 1)];
307 
308  // bring them into standard orientation
309  if (vertex_indices[0] > vertex_indices[1])
310  std::swap (vertex_indices[0], vertex_indices[1]);
311  }
312 
317  bool operator< (const Edge<dim> &e) const
318  {
319  return ((vertex_indices[0] < e.vertex_indices[0])
320  ||
321  ((vertex_indices[0] == e.vertex_indices[0]) && (vertex_indices[1] < e.vertex_indices[1])));
322  }
323 
327  bool operator== (const Edge<dim> &e) const
328  {
329  return ((vertex_indices[0] == e.vertex_indices[0])
330  &&
331  (vertex_indices[1] == e.vertex_indices[1]));
332  }
333 
338  unsigned int vertex_indices[2];
339 
344  enum OrientationStatus
345  {
346  not_oriented,
347  forward,
348  backward
349  };
350 
351  OrientationStatus orientation_status;
352 
357  AdjacentCells<dim> adjacent_cells;
358  };
359 
360 
361 
366  template <int dim>
367  struct Cell
368  {
374  Cell (const CellData<dim> &c,
375  const std::vector<Edge<dim> > &edge_list)
376  {
377  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
378  vertex_indices[i] = c.vertices[i];
379 
380  // now for each of the edges of this cell, find the location inside the
381  // given edge_list array and store than index
382  for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
383  {
384  const Edge<dim> e (c, l);
385  edge_indices[l] = (std::lower_bound (edge_list.begin(), edge_list.end(), e)
386  -
387  edge_list.begin());
388  Assert (edge_indices[l] < edge_list.size(), ExcInternalError());
389  Assert (edge_list[edge_indices[l]] == e, ExcInternalError())
390  }
391  }
392 
396  unsigned int vertex_indices[GeometryInfo<dim>::vertices_per_cell];
397 
402  unsigned int edge_indices[GeometryInfo<dim>::lines_per_cell];
403  };
404 
405 
406 
407  template <int dim> class EdgeDeltaSet;
408 
418  template <>
419  class EdgeDeltaSet<2>
420  {
421  public:
425  typedef const unsigned int *const_iterator;
426 
431  EdgeDeltaSet ()
432  {
433  edge_indices[0] = edge_indices[1] = numbers::invalid_unsigned_int;
434  }
435 
436 
440  void clear ()
441  {
442  edge_indices[0] = edge_indices[1] = numbers::invalid_unsigned_int;
443  }
444 
449  void insert (const unsigned int edge_index)
450  {
451  if (edge_indices[0] == numbers::invalid_unsigned_int)
452  edge_indices[0] = edge_index;
453  else
454  {
455  Assert (edge_indices[1] == numbers::invalid_unsigned_int,
456  ExcInternalError());
457  edge_indices[1] = edge_index;
458  }
459  }
460 
461 
465  const_iterator begin () const
466  {
467  return edge_indices;
468  }
469 
470 
474  const_iterator end () const
475  {
476  // check whether the current object stores zero, one, or two
477  // indices, and use this to point to the element past the
478  // last valid one
479  if (edge_indices[0] == numbers::invalid_unsigned_int)
480  return edge_indices;
481  else if (edge_indices[1] == numbers::invalid_unsigned_int)
482  return edge_indices + 1;
483  else
484  return edge_indices + 2;
485  }
486 
487  private:
491  unsigned int edge_indices[2];
492  };
493 
494 
495 
507  template <>
508  class EdgeDeltaSet<3> : public std::set<unsigned int>
509  {};
510 
511 
512 
513 
518  template <int dim>
519  std::vector<Edge<dim> >
520  build_edges (const std::vector<CellData<dim> > &cells)
521  {
522  // build the edge list for all cells. because each cell has
523  // GeometryInfo<dim>::lines_per_cell edges, the total number
524  // of edges is this many times the number of cells. of course
525  // some of them will be duplicates, and we throw them out below
526  std::vector<Edge<dim> > edge_list;
527  edge_list.reserve(cells.size()*GeometryInfo<dim>::lines_per_cell);
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);
531 
532  // next sort the edge list and then remove duplicates
533  std::sort (edge_list.begin(), edge_list.end());
534  edge_list.erase(std::unique(edge_list.begin(),edge_list.end()),
535  edge_list.end());
536 
537  return edge_list;
538  }
539 
540 
541 
546  template <int dim>
547  std::vector<Cell<dim> >
548  build_cells_and_connect_edges (const std::vector<CellData<dim> > &cells,
549  std::vector<Edge<dim> > &edges)
550  {
551  std::vector<Cell<dim> > cell_list;
552  cell_list.reserve(cells.size());
553  for (unsigned int i=0; i<cells.size(); ++i)
554  {
555  // create our own data structure for the cells and let it
556  // connect to the edges array
557  cell_list.emplace_back (cells[i], edges);
558 
559  // then also inform the edges that they are adjacent
560  // to the current cell, and where within this cell
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));
563  }
564  Assert (cell_list.size() == cells.size(), ExcInternalError());
565 
566  return cell_list;
567  }
568 
569 
570 
575  template <int dim>
576  unsigned int
577  get_next_unoriented_cell(const std::vector<Cell<dim> > &cells,
578  const std::vector<Edge<dim> > &edges,
579  const unsigned int current_cell)
580  {
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)
584  return c;
585 
587  }
588 
589 
590 
596  template <int dim>
597  void
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)
602  {
603  // choose the direction of the first edge. we have free choice
604  // here and could simply choose "forward" if that's what pleases
605  // us. however, for backward compatibility with the previous
606  // implementation used till 2016, let us just choose the
607  // direction so that it matches what we have in the given cell.
608  //
609  // in fact, in what can only be assumed to be a bug in the
610  // original implementation, after orienting all edges, the code
611  // that rotates the cells so that they match edge orientations
612  // (see the rotate_cell() function below) rotated the cell two
613  // more times by 90 degrees. this is ok -- it simply flips all
614  // edge orientations, which leaves them valid. rather than do
615  // the same in the current implementation, we can achieve the
616  // same effect by modifying the rule above to choose the
617  // direction of the starting edge of this parallel set
618  // *opposite* to what it looks like in the current cell
619  //
620  // this bug only existed in the 2d implementation since there
621  // were different implementations for 2d and 3d. consequently,
622  // only replicate it for the 2d case and be "intuitive" in 3d.
623  if (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
624  ==
625  cells[cell].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices (local_edge, 0)])
626  // orient initial edge *opposite* to the way it is in the cell
627  // (see above for the reason)
628  edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
629  Edge<dim>::backward :
630  Edge<dim>::forward);
631  else
632  {
633  Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
634  ==
635  cells[cell].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices (local_edge, 1)],
636  ExcInternalError());
637  Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[1]
638  ==
639  cells[cell].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices (local_edge, 0)],
640  ExcInternalError());
641 
642  // orient initial edge *opposite* to the way it is in the cell
643  // (see above for the reason)
644  edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
645  Edge<dim>::forward :
646  Edge<dim>::backward);
647  }
648 
649  // walk outward from the given edge as described in
650  // the algorithm in the paper that documents all of
651  // this
652  //
653  // note that in 2d, each of the Deltas can at most
654  // contain two elements, whereas in 3d it can be arbitrarily many
655  EdgeDeltaSet<dim> Delta_k;
656  EdgeDeltaSet<dim> Delta_k_minus_1;
657  Delta_k_minus_1.insert (cells[cell].edge_indices[local_edge]);
658 
659  while (Delta_k_minus_1.begin() != Delta_k_minus_1.end()) // while set is not empty
660  {
661  Delta_k.clear ();
662 
663  for (typename EdgeDeltaSet<dim>::const_iterator delta = Delta_k_minus_1.begin();
664  delta != Delta_k_minus_1.end(); ++delta)
665  {
666  Assert (edges[*delta].orientation_status != Edge<dim>::not_oriented,
667  ExcInternalError());
668 
669  // now go through the cells adjacent to this edge
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)
673  {
674  const unsigned int K = adjacent_cell->cell_index;
675  const unsigned int delta_is_edge_in_K = adjacent_cell->edge_within_cell;
676 
677  // figure out the direction of delta with respect to the cell K
678  // (in the orientation in which the user has given it to us)
679  const unsigned int first_edge_vertex
680  = (edges[*delta].orientation_status == Edge<dim>::forward
681  ?
682  edges[*delta].vertex_indices[0]
683  :
684  edges[*delta].vertex_indices[1]);
685  const unsigned int first_edge_vertex_in_K
686  = cells[K].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices(delta_is_edge_in_K, 0)];
687  Assert (first_edge_vertex == first_edge_vertex_in_K
688  ||
689  first_edge_vertex == cells[K].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices(delta_is_edge_in_K, 1)],
690  ExcInternalError());
691 
692  // now figure out which direction the each of the "opposite" edges
693  // needs to be oriented into.
694  for (unsigned int o_e=0; o_e<ParallelEdges<dim>::n_other_parallel_edges; ++o_e)
695  {
696  // get the index of the opposite edge and select which its first
697  // vertex needs to be based on how the current edge is oriented
698  // in the current cell
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
702  = cells[K].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices(
703  ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K][o_e],
704  (first_edge_vertex == first_edge_vertex_in_K
705  ?
706  0
707  :
708  1))];
709 
710  // then determine the orientation of the edge based on
711  // whether the vertex we want to be the edge's first
712  // vertex is already the first vertex of the edge, or
713  // whether it points in the opposite direction
714  const typename Edge<dim>::OrientationStatus opposite_edge_orientation
715  = (edges[opposite_edge].vertex_indices[0]
716  ==
717  first_opposite_edge_vertex
718  ?
719  Edge<dim>::forward
720  :
721  Edge<dim>::backward);
722 
723  // see if the opposite edge (there is only one in 2d) has already been
724  // oriented.
725  if (edges[opposite_edge].orientation_status == Edge<dim>::not_oriented)
726  {
727  // the opposite edge is not yet oriented. do orient it and add it to
728  // Delta_k
729  edges[opposite_edge].orientation_status = opposite_edge_orientation;
730  Delta_k.insert (opposite_edge);
731  }
732  else
733  {
734  // this opposite edge has already been oriented. it should be
735  // consistent with the current one in 2d, while in 3d it may in fact
736  // be mis-oriented, and in that case the mesh will not be
737  // orientable. indicate this by throwing an exception that we can
738  // catch further up; this has the advantage that we can propagate
739  // through a couple of functions without having to do error
740  // checking and without modifying the 'cells' array that the
741  // user gave us
742  if (dim == 2)
743  {
744  Assert (edges[opposite_edge].orientation_status == opposite_edge_orientation,
746  }
747  else if (dim == 3)
748  {
749  if (edges[opposite_edge].orientation_status != opposite_edge_orientation)
750  throw ExcMeshNotOrientable ();
751  }
752  else
753  Assert (false, ExcNotImplemented());
754  }
755  }
756  }
757  }
758 
759  // finally copy the new set to the previous one
760  // (corresponding to increasing 'k' by one in the
761  // algorithm)
762  Delta_k_minus_1 = Delta_k;
763  }
764  }
765 
766 
774  template <int dim>
775  void
776  rotate_cell (const std::vector<Cell<dim> > &cell_list,
777  const std::vector<Edge<dim> > &edge_list,
778  const unsigned int cell_index,
779  std::vector<CellData<dim> > &raw_cells)
780  {
781  // find the first vertex of the cell. this is the vertex where dim edges
782  // originate, so for each of the edges record which the starting vertex is
783  unsigned int starting_vertex_of_edge[GeometryInfo<dim>::lines_per_cell];
784  for (unsigned int e=0; e<GeometryInfo<dim>::lines_per_cell; ++e)
785  {
786  Assert (edge_list[cell_list[cell_index].edge_indices[e]].orientation_status
787  != Edge<dim>::not_oriented,
788  ExcInternalError());
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];
791  else
792  starting_vertex_of_edge[e] = edge_list[cell_list[cell_index].edge_indices[e]].vertex_indices[1];
793  }
794 
795  // find the vertex number that appears dim times. this will then be
796  // the vertex at which we want to locate the origin of the cell's
797  // coordinate system (i.e., vertex 0)
798  unsigned int origin_vertex_of_cell = numbers::invalid_unsigned_int;
799  switch (dim)
800  {
801  case 2:
802  {
803  // in 2d, we can simply enumerate the possibilities where the
804  // origin may be located because edges zero and one don't share
805  // any vertices, and the same for edges two and three
806  if ((starting_vertex_of_edge[0] == starting_vertex_of_edge[2])
807  ||
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])
811  ||
812  (starting_vertex_of_edge[1] == starting_vertex_of_edge[3]))
813  origin_vertex_of_cell = starting_vertex_of_edge[1];
814  else
815  Assert (false, ExcInternalError());
816 
817  break;
818  }
819 
820  case 3:
821  {
822  // one could probably do something similar in 3d, but that seems
823  // more complicated than one wants to write down. just go
824  // through the list of possible starting vertices and check
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,
829  starting_vertex_of_edge + GeometryInfo<dim>::lines_per_cell,
830  cell_list[cell_index].vertex_indices[origin_vertex_of_cell])
831  == dim)
832  break;
833  Assert (origin_vertex_of_cell < GeometryInfo<dim>::vertices_per_cell,
834  ExcInternalError());
835 
836  break;
837  }
838 
839  default:
840  Assert (false, ExcNotImplemented());
841  }
842 
843  // now rotate raw_cells[cell_index] in such a way that its orientation
844  // matches that of cell_list[cell_index]
845  switch (dim)
846  {
847  case 2:
848  {
849  // in 2d, we can literally rotate the cell until its origin
850  // matches the one that we have determined above should be
851  // the origin vertex
852  //
853  // when doing a rotation, take into account the ordering of
854  // vertices (not in clockwise or counter-clockwise sense)
855  while (raw_cells[cell_index].vertices[0] != origin_vertex_of_cell)
856  {
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;
862  }
863  break;
864  }
865 
866  case 3:
867  {
868  // in 3d, the situation is a bit more complicated. from above, we
869  // now know which vertex is at the origin (because 3 edges originate
870  // from it), but that still leaves 3 possible rotations of the cube.
871  // the important realization is that we can choose any of them:
872  // in all 3 rotations, all edges originate from the one vertex,
873  // and that fixes the directions of all 12 edges in the cube because
874  // these 3 cover all 3 equivalence classes! consequently, we can
875  // select an arbitrary one among the permutations -- for
876  // example the following ones:
877  static const unsigned int cube_permutations[8][8] =
878  {
879  {0,1,2,3,4,5,6,7},
880  {1,5,3,7,0,4,2,6},
881  {2,6,0,4,3,7,1,5},
882  {3,2,1,0,7,6,5,4},
883  {4,0,6,2,5,1,7,3},
884  {5,4,7,6,1,0,3,2},
885  {6,7,4,5,2,3,0,1},
886  {7,3,5,1,6,2,4,0}
887  };
888 
889  unsigned int temp_vertex_indices[GeometryInfo<dim>::vertices_per_cell];
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];
895 
896  break;
897  }
898 
899  default:
900  {
901  Assert (false, ExcNotImplemented());
902  }
903  }
904  }
905 
906 
912  template <int dim>
913  void reorient (std::vector<CellData<dim> > &cells)
914  {
915  // first build the arrays that connect cells to edges and the other
916  // way around
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);
919 
920  // then loop over all cells and start orienting parallel edge sets
921  // of cells that still have non-oriented edges
922  unsigned int next_cell_with_unoriented_edge = 0;
923  while ((next_cell_with_unoriented_edge = get_next_unoriented_cell(cell_list,
924  edge_list,
925  next_cell_with_unoriented_edge)) !=
927  {
928  // see which edge sets are still not oriented
929  //
930  // we do not need to look at each edge because if we orient edge
931  // 0, we will end up with edge 1 also oriented (in 2d; in 3d, there
932  // will be 3 other edges that are also oriented). there are only
933  // dim independent sets of edges, so loop over these.
934  //
935  // we need to check whether each one of these starter edges may
936  // already be oriented because the line (sheet) that connects
937  // globally parallel edges may be self-intersecting in the
938  // current cell
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,
943  edge_list,
944  next_cell_with_unoriented_edge,
945  ParallelEdges<dim>::starter_edges[l]);
946 
947  // ensure that we have really oriented all edges now, not just
948  // the starter edges
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,
952  ExcInternalError());
953  }
954 
955  // now that we have oriented all edges, we need to rotate cells
956  // so that the edges point in the right direction with the now
957  // rotated coordinate system
958  for (unsigned int c=0; c<cells.size(); ++c)
959  rotate_cell (cell_list, edge_list, c, cells);
960  }
961 
962 
963  // overload of the function above for 1d -- there is nothing
964  // to orient in that case
965  void reorient (std::vector<CellData<1> > &)
966  {}
967 }
968 
969 
970 // anonymous namespace for internal helper functions
971 namespace
972 {
982  void
983  reorder_new_to_old_style (std::vector<CellData<1> > &)
984  {}
985 
986 
987  void
988  reorder_new_to_old_style (std::vector<CellData<2> > &cells)
989  {
990  for (unsigned int cell=0; cell<cells.size(); ++cell)
991  std::swap(cells[cell].vertices[2], cells[cell].vertices[3]);
992  }
993 
994 
995  void
996  reorder_new_to_old_style (std::vector<CellData<3> > &cells)
997  {
998  unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
999  for (unsigned int cell=0; cell<cells.size(); ++cell)
1000  {
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)
1004  cells[cell].vertices[i] = tmp[GeometryInfo<3>::ucd_to_deal[i]];
1005  }
1006  }
1007 
1008 
1012  void
1013  reorder_old_to_new_style (std::vector<CellData<1> > &)
1014  {}
1015 
1016 
1017  void
1018  reorder_old_to_new_style (std::vector<CellData<2> > &cells)
1019  {
1020  // just invert the permutation:
1021  reorder_new_to_old_style(cells);
1022  }
1023 
1024 
1025  void
1026  reorder_old_to_new_style (std::vector<CellData<3> > &cells)
1027  {
1028  // undo the ordering above
1029  unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
1030  for (unsigned int cell=0; cell<cells.size(); ++cell)
1031  {
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)
1035  cells[cell].vertices[GeometryInfo<3>::ucd_to_deal[i]] = tmp[i];
1036  }
1037  }
1038 }
1039 
1040 
1041 
1042 template <int dim, int spacedim>
1043 void
1045  const bool use_new_style_ordering)
1046 {
1047  Assert (cells.size() != 0,
1048  ExcMessage("List of elements to orient must have at least one cell"));
1049 
1050  // there is nothing for us to do in 1d
1051  if (dim == 1)
1052  return;
1053 
1054  // if necessary, convert to new-style format
1055  if (use_new_style_ordering == false)
1056  reorder_old_to_new_style(cells);
1057 
1058  // check if grids are already consistent. if so, do
1059  // nothing. if not, then do the reordering
1060  if (!is_consistent (cells))
1061  try
1062  {
1063  reorient(cells);
1064  }
1065  catch (const ExcMeshNotOrientable &)
1066  {
1067  // the mesh is not orientable. this is acceptable if we are in 3d,
1068  // as class Triangulation knows how to handle this, but it is
1069  // not in 2d; in that case, re-throw the exception
1070  if (dim < 3)
1071  throw;
1072  }
1073 
1074  // and convert back if necessary
1075  if (use_new_style_ordering == false)
1076  reorder_new_to_old_style(cells);
1077 }
1078 
1079 
1080 
1081 template <>
1082 void
1084  std::vector<CellData<1> > &)
1085 {
1086  // nothing to be done in 1d
1087 }
1088 
1089 
1090 
1091 template <>
1092 void
1094  std::vector<CellData<1> > &)
1095 {
1096  // nothing to be done in 1d
1097 }
1098 
1099 
1100 
1101 template <>
1102 void
1104  std::vector<CellData<1> > &)
1105 {
1106  // nothing to be done in 1d
1107 }
1108 
1109 
1110 template <>
1111 void
1112 GridReordering<2>::invert_all_cells_of_negative_grid(const std::vector<Point<2> > &all_vertices,
1113  std::vector<CellData<2> > &cells)
1114 {
1115  unsigned int vertices_lex[GeometryInfo<2>::vertices_per_cell];
1116  unsigned int n_negative_cells=0;
1117  for (unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1118  {
1119  // GridTools::cell_measure
1120  // requires the vertices to be
1121  // in lexicographic ordering
1122  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1123  vertices_lex[GeometryInfo<2>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1124  if (GridTools::cell_measure<2>(all_vertices, vertices_lex) < 0)
1125  {
1126  ++n_negative_cells;
1127  std::swap(cells[cell_no].vertices[1], cells[cell_no].vertices[3]);
1128 
1129  // check whether the
1130  // resulting cell is now ok.
1131  // if not, then the grid is
1132  // seriously broken and
1133  // should be sticked into the
1134  // bin
1135  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1136  vertices_lex[GeometryInfo<2>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1137  AssertThrow(GridTools::cell_measure<2>(all_vertices, vertices_lex) > 0,
1138  ExcInternalError());
1139  }
1140  }
1141 
1142  // We assume that all cells of a grid have
1143  // either positive or negative volumes but
1144  // not both mixed. Although above reordering
1145  // might work also on single cells, grids
1146  // with both kind of cells are very likely to
1147  // be broken. Check for this here.
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 ")
1158  + Utilities::to_string (cells.size())
1159  + " cells in the mesh, "
1160  + Utilities::to_string (n_negative_cells)
1161  + " appear to have a negative volume."));
1162 }
1163 
1164 
1165 
1166 template <>
1167 void
1169  std::vector<CellData<2> > &)
1170 {
1171  Assert(false, ExcNotImplemented());
1172 }
1173 
1174 
1175 
1176 template <>
1177 void
1178 GridReordering<3>::invert_all_cells_of_negative_grid(const std::vector<Point<3> > &all_vertices,
1179  std::vector<CellData<3> > &cells)
1180 {
1181  unsigned int vertices_lex[GeometryInfo<3>::vertices_per_cell];
1182  unsigned int n_negative_cells=0;
1183  for (unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1184  {
1185  // GridTools::cell_measure
1186  // requires the vertices to be
1187  // in lexicographic ordering
1188  for (unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1189  vertices_lex[GeometryInfo<3>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1190  if (GridTools::cell_measure<3>(all_vertices, vertices_lex) < 0)
1191  {
1192  ++n_negative_cells;
1193  // reorder vertices: swap front and back face
1194  for (unsigned int i=0; i<4; ++i)
1195  std::swap(cells[cell_no].vertices[i], cells[cell_no].vertices[i+4]);
1196 
1197  // check whether the
1198  // resulting cell is now ok.
1199  // if not, then the grid is
1200  // seriously broken and
1201  // should be sticked into the
1202  // bin
1203  for (unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1204  vertices_lex[GeometryInfo<3>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1205  AssertThrow(GridTools::cell_measure<3>(all_vertices, vertices_lex) > 0,
1206  ExcInternalError());
1207  }
1208  }
1209 
1210  // We assume that all cells of a
1211  // grid have either positive or
1212  // negative volumes but not both
1213  // mixed. Although above reordering
1214  // might work also on single cells,
1215  // grids with both kind of cells
1216  // are very likely to be
1217  // broken. Check for this here.
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."));
1237 }
1238 
1239 
1240 
1241 /* ------------------------ explicit instantiations ------------------- */
1242 template class GridReordering<1,1>;
1243 template class GridReordering<1,2>;
1244 template class GridReordering<1,3>;
1245 template class GridReordering<2,2>;
1246 template class GridReordering<2,3>;
1247 template class GridReordering<3,3>;
1248 
1249 DEAL_II_NAMESPACE_CLOSE
static const unsigned int invalid_unsigned_int
Definition: types.h:173
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)
Definition: exceptions.h:1221
Definition: point.h:104
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:107
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1142
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)
Definition: vector.h:1312
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]
Definition: tria.h:132
static ::ExceptionBase & ExcInternalError()