Reference documentation for deal.II version 8.5.1
grid_reordering.cc
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2000 - 2016 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/std_cxx11/bind.h>
21 #include <deal.II/base/timer.h>
22 
23 #include <algorithm>
24 #include <set>
25 #include <iostream>
26 #include <fstream>
27 #include <functional>
28 
29 DEAL_II_NAMESPACE_OPEN
30 
31 
32 namespace
33 {
39  struct CheapEdge
40  {
44  CheapEdge (const unsigned int v0,
45  const unsigned int v1)
46  :
47  v0(v0), v1(v1)
48  {}
49 
54  bool operator < (const CheapEdge &e) const
55  {
56  return ((v0 < e.v0) || ((v0 == e.v0) && (v1 < e.v1)));
57  }
58 
59  private:
63  const unsigned int v0, v1;
64  };
65 
66 
75  template <int dim>
76  bool
77  is_consistent (const std::vector<CellData<dim> > &cells)
78  {
79  std::set<CheapEdge> edges;
80 
81  for (typename std::vector<CellData<dim> >::const_iterator c = cells.begin();
82  c != cells.end(); ++c)
83  {
84  // construct the edges in reverse order. for each of them,
85  // ensure that the reverse edge is not yet in the list of
86  // edges (return false if the reverse edge already *is* in
87  // the list) and then add the actual edge to it; std::set
88  // eliminates duplicates automatically
89  for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
90  {
91  const CheapEdge reverse_edge (c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 1)],
92  c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 0)]);
93  if (edges.find (reverse_edge) != edges.end())
94  return false;
95 
96 
97  // ok, not. insert edge in correct order
98  const CheapEdge correct_edge (c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 0)],
99  c->vertices[GeometryInfo<dim>::line_to_cell_vertices(l, 1)]);
100  edges.insert (correct_edge);
101  }
102  }
103 
104  // no conflicts found, so return true
105  return true;
106  }
107 
108 
115  template <int dim>
116  struct ParallelEdges
117  {
123  static const unsigned int starter_edges[dim];
124 
129  static const unsigned int n_other_parallel_edges = (1<<(dim-1)) - 1;
130  static const unsigned int parallel_edges[GeometryInfo<dim>::lines_per_cell][n_other_parallel_edges];
131  };
132 
133  template <>
134  const unsigned int ParallelEdges<2>::starter_edges[2] = { 0, 2 };
135 
136  template <>
137  const unsigned int ParallelEdges<2>::parallel_edges[4][1] = { {1}, {0}, {3}, {2} };
138 
139  template <>
140  const unsigned int ParallelEdges<3>::starter_edges[3] = { 0, 2, 8 };
141 
142  template <>
143  const unsigned int ParallelEdges<3>::parallel_edges[12][3] = { {1, 4, 5}, // line 0
144  {0, 4, 5}, // line 1
145  {3, 6, 7}, // line 2
146  {2, 6, 7}, // line 3
147  {0, 1, 5}, // line 4
148  {0, 1, 4}, // line 5
149  {2, 3, 7}, // line 6
150  {2, 3, 6}, // line 7
151  {9, 10, 11}, // line 8
152  {8, 10, 11}, // line 9
153  {8, 9, 11}, // line 10
154  {8, 9, 10} // line 11
155  };
156 
157 
162  struct AdjacentCell
163  {
167  AdjacentCell ()
168  :
169  cell_index (numbers::invalid_unsigned_int),
170  edge_within_cell (numbers::invalid_unsigned_int)
171  {}
172 
176  AdjacentCell (const unsigned int cell_index,
177  const unsigned int edge_within_cell)
178  :
179  cell_index (cell_index),
180  edge_within_cell (edge_within_cell)
181  {}
182 
183 
184  unsigned int cell_index;
185  unsigned int edge_within_cell;
186  };
187 
188 
189 
190  template <int dim> class AdjacentCells;
191 
197  template <>
198  class AdjacentCells<2>
199  {
200  public:
205  typedef const AdjacentCell *const_iterator;
206 
215  void push_back (const AdjacentCell &adjacent_cell)
216  {
217  if (adjacent_cells[0].cell_index == numbers::invalid_unsigned_int)
218  adjacent_cells[0] = adjacent_cell;
219  else
220  {
221  Assert (adjacent_cells[1].cell_index == numbers::invalid_unsigned_int,
222  ExcInternalError());
223  adjacent_cells[1] = adjacent_cell;
224  }
225  }
226 
227 
232  const_iterator begin () const
233  {
234  return &adjacent_cells[0];
235  }
236 
237 
243  const_iterator end () const
244  {
245  // check whether the current object stores zero, one, or two
246  // adjacent cells, and use this to point to the element past the
247  // last valid one
248  if (adjacent_cells[0].cell_index == numbers::invalid_unsigned_int)
249  return &adjacent_cells[0];
250  else if (adjacent_cells[1].cell_index == numbers::invalid_unsigned_int)
251  return &adjacent_cells[0]+1;
252  else
253  return &adjacent_cells[0]+2;
254  }
255 
256  private:
263  AdjacentCell adjacent_cells[2];
264  };
265 
266 
267 
275  template <>
276  class AdjacentCells<3> : public std::vector<AdjacentCell>
277  {};
278 
279 
289  template <int dim>
290  class Edge
291  {
292  public:
296  Edge ()
297  :
298  orientation_status (not_oriented)
299  {
300  for (unsigned int i=0; i<2; ++i)
301  vertex_indices[i] = numbers::invalid_unsigned_int;
302  }
303 
309  Edge (const CellData<dim> &cell,
310  const unsigned int edge_number)
311  :
312  orientation_status (not_oriented)
313  {
315 
316  // copy vertices for this particular line
317  vertex_indices[0] = cell.vertices[GeometryInfo<dim>::line_to_cell_vertices(edge_number, 0)];
318  vertex_indices[1] = cell.vertices[GeometryInfo<dim>::line_to_cell_vertices(edge_number, 1)];
319 
320  // bring them into standard orientation
321  if (vertex_indices[0] > vertex_indices[1])
322  std::swap (vertex_indices[0], vertex_indices[1]);
323  }
324 
329  bool operator< (const Edge<dim> &e) const
330  {
331  return ((vertex_indices[0] < e.vertex_indices[0])
332  ||
333  ((vertex_indices[0] == e.vertex_indices[0]) && (vertex_indices[1] < e.vertex_indices[1])));
334  }
335 
339  bool operator== (const Edge<dim> &e) const
340  {
341  return ((vertex_indices[0] == e.vertex_indices[0])
342  &&
343  (vertex_indices[1] == e.vertex_indices[1]));
344  }
345 
350  unsigned int vertex_indices[2];
351 
356  enum OrientationStatus
357  {
358  not_oriented,
359  forward,
360  backward
361  };
362 
363  OrientationStatus orientation_status;
364 
369  AdjacentCells<dim> adjacent_cells;
370  };
371 
372 
373 
378  template <int dim>
379  struct Cell
380  {
384  Cell ()
385  {
386  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
387  vertex_indices[i] = numbers::invalid_unsigned_int;
388  for (unsigned int i=0; i<GeometryInfo<dim>::lines_per_cell; ++i)
389  edge_indices[i] = numbers::invalid_unsigned_int;
390  }
391 
397  Cell (const CellData<dim> &c,
398  const std::vector<Edge<dim> > &edge_list)
399  {
400  for (unsigned int i=0; i<GeometryInfo<dim>::vertices_per_cell; ++i)
401  vertex_indices[i] = c.vertices[i];
402 
403  // now for each of the edges of this cell, find the location inside the
404  // given edge_list array and store than index
405  for (unsigned int l=0; l<GeometryInfo<dim>::lines_per_cell; ++l)
406  {
407  const Edge<dim> e (c, l);
408  edge_indices[l] = (std::lower_bound (edge_list.begin(), edge_list.end(), e)
409  -
410  edge_list.begin());
411  Assert (edge_indices[l] < edge_list.size(), ExcInternalError());
412  Assert (edge_list[edge_indices[l]] == e, ExcInternalError())
413  }
414  }
415 
419  unsigned int vertex_indices[GeometryInfo<dim>::vertices_per_cell];
420 
425  unsigned int edge_indices[GeometryInfo<dim>::lines_per_cell];
426  };
427 
428 
429 
430  template <int dim> class EdgeDeltaSet;
431 
441  template <>
442  class EdgeDeltaSet<2>
443  {
444  public:
448  typedef const unsigned int *const_iterator;
449 
454  EdgeDeltaSet ()
455  {
456  edge_indices[0] = edge_indices[1] = numbers::invalid_unsigned_int;
457  }
458 
459 
463  void clear ()
464  {
465  edge_indices[0] = edge_indices[1] = numbers::invalid_unsigned_int;
466  }
467 
472  void insert (const unsigned int edge_index)
473  {
474  if (edge_indices[0] == numbers::invalid_unsigned_int)
475  edge_indices[0] = edge_index;
476  else
477  {
478  Assert (edge_indices[1] == numbers::invalid_unsigned_int,
479  ExcInternalError());
480  edge_indices[1] = edge_index;
481  }
482  }
483 
484 
488  const_iterator begin () const
489  {
490  return &edge_indices[0];
491  }
492 
493 
497  const_iterator end () const
498  {
499  // check whether the current object stores zero, one, or two
500  // indices, and use this to point to the element past the
501  // last valid one
502  if (edge_indices[0] == numbers::invalid_unsigned_int)
503  return &edge_indices[0];
504  else if (edge_indices[1] == numbers::invalid_unsigned_int)
505  return &edge_indices[0]+1;
506  else
507  return &edge_indices[0]+2;
508  }
509 
510  private:
514  unsigned int edge_indices[2];
515  };
516 
517 
518 
530  template <>
531  class EdgeDeltaSet<3> : public std::set<unsigned int>
532  {};
533 
534 
535 
536 
541  template <int dim>
542  std::vector<Edge<dim> >
543  build_edges (const std::vector<CellData<dim> > &cells)
544  {
545  // build the edge list for all cells. because each cell has
546  // GeometryInfo<dim>::lines_per_cell edges, the total number
547  // of edges is this many times the number of cells. of course
548  // some of them will be duplicates, and we throw them out below
549  std::vector<Edge<dim> > edge_list;
550  edge_list.reserve(cells.size()*GeometryInfo<dim>::lines_per_cell);
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));
554 
555  // next sort the edge list and then remove duplicates
556  std::sort (edge_list.begin(), edge_list.end());
557  edge_list.erase(std::unique(edge_list.begin(),edge_list.end()),
558  edge_list.end());
559 
560  return edge_list;
561  }
562 
563 
564 
569  template <int dim>
570  std::vector<Cell<dim> >
571  build_cells_and_connect_edges (const std::vector<CellData<dim> > &cells,
572  std::vector<Edge<dim> > &edges)
573  {
574  std::vector<Cell<dim> > cell_list;
575  cell_list.reserve(cells.size());
576  for (unsigned int i=0; i<cells.size(); ++i)
577  {
578  // create our own data structure for the cells and let it
579  // connect to the edges array
580  cell_list.push_back (Cell<dim>(cells[i], edges));
581 
582  // then also inform the edges that they are adjacent
583  // to the current cell, and where within this cell
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));
586  }
587  Assert (cell_list.size() == cells.size(), ExcInternalError());
588 
589  return cell_list;
590  }
591 
592 
593 
598  template <int dim>
599  unsigned int
600  get_next_unoriented_cell(const std::vector<Cell<dim> > &cells,
601  const std::vector<Edge<dim> > &edges,
602  const unsigned int current_cell)
603  {
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)
607  return c;
608 
610  }
611 
612 
613 
619  template <int dim>
620  void
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)
625  {
626  // choose the direction of the first edge. we have free choice
627  // here and could simply choose "forward" if that's what pleases
628  // us. however, for backward compatibility with the previous
629  // implementation used till 2016, let us just choose the
630  // direction so that it matches what we have in the given cell.
631  //
632  // in fact, in what can only be assumed to be a bug in the
633  // original implementation, after orienting all edges, the code
634  // that rotates the cells so that they match edge orientations
635  // (see the rotate_cell() function below) rotated the cell two
636  // more times by 90 degrees. this is ok -- it simply flips all
637  // edge orientations, which leaves them valid. rather than do
638  // the same in the current implementation, we can achieve the
639  // same effect by modifying the rule above to choose the
640  // direction of the starting edge of this parallel set
641  // *opposite* to what it looks like in the current cell
642  //
643  // this bug only existed in the 2d implementation since there
644  // were different implementations for 2d and 3d. consequently,
645  // only replicate it for the 2d case and be "intuitive" in 3d.
646  if (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
647  ==
648  cells[cell].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices (local_edge, 0)])
649  // orient initial edge *opposite* to the way it is in the cell
650  // (see above for the reason)
651  edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
652  Edge<dim>::backward :
653  Edge<dim>::forward);
654  else
655  {
656  Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[0]
657  ==
658  cells[cell].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices (local_edge, 1)],
659  ExcInternalError());
660  Assert (edges[cells[cell].edge_indices[local_edge]].vertex_indices[1]
661  ==
662  cells[cell].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices (local_edge, 0)],
663  ExcInternalError());
664 
665  // orient initial edge *opposite* to the way it is in the cell
666  // (see above for the reason)
667  edges[cells[cell].edge_indices[local_edge]].orientation_status = (dim == 2 ?
668  Edge<dim>::forward :
669  Edge<dim>::backward);
670  }
671 
672  // walk outward from the given edge as described in
673  // the algorithm in the paper that documents all of
674  // this
675  //
676  // note that in 2d, each of the Deltas can at most
677  // contain two elements, whereas in 3d it can be arbitrarily many
678  EdgeDeltaSet<dim> Delta_k;
679  EdgeDeltaSet<dim> Delta_k_minus_1;
680  Delta_k_minus_1.insert (cells[cell].edge_indices[local_edge]);
681 
682  while (Delta_k_minus_1.begin() != Delta_k_minus_1.end()) // while set is not empty
683  {
684  Delta_k.clear ();
685 
686  for (typename EdgeDeltaSet<dim>::const_iterator delta = Delta_k_minus_1.begin();
687  delta != Delta_k_minus_1.end(); ++delta)
688  {
689  Assert (edges[*delta].orientation_status != Edge<dim>::not_oriented,
690  ExcInternalError());
691 
692  // now go through the cells adjacent to this edge
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)
696  {
697  const unsigned int K = adjacent_cell->cell_index;
698  const unsigned int delta_is_edge_in_K = adjacent_cell->edge_within_cell;
699 
700  // figure out the direction of delta with respect to the cell K
701  // (in the orientation in which the user has given it to us)
702  const unsigned int first_edge_vertex
703  = (edges[*delta].orientation_status == Edge<dim>::forward
704  ?
705  edges[*delta].vertex_indices[0]
706  :
707  edges[*delta].vertex_indices[1]);
708  const unsigned int first_edge_vertex_in_K
709  = cells[K].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices(delta_is_edge_in_K, 0)];
710  Assert (first_edge_vertex == first_edge_vertex_in_K
711  ||
712  first_edge_vertex == cells[K].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices(delta_is_edge_in_K, 1)],
713  ExcInternalError());
714 
715  // now figure out which direction the each of the "opposite" edges
716  // needs to be oriented into.
717  for (unsigned int o_e=0; o_e<ParallelEdges<dim>::n_other_parallel_edges; ++o_e)
718  {
719  // get the index of the opposite edge and select which its first
720  // vertex needs to be based on how the current edge is oriented
721  // in the current cell
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
725  = cells[K].vertex_indices[GeometryInfo<dim>::line_to_cell_vertices(
726  ParallelEdges<dim>::parallel_edges[delta_is_edge_in_K][o_e],
727  (first_edge_vertex == first_edge_vertex_in_K
728  ?
729  0
730  :
731  1))];
732 
733  // then determine the orientation of the edge based on
734  // whether the vertex we want to be the edge's first
735  // vertex is already the first vertex of the edge, or
736  // whether it points in the opposite direction
737  const typename Edge<dim>::OrientationStatus opposite_edge_orientation
738  = (edges[opposite_edge].vertex_indices[0]
739  ==
740  first_opposite_edge_vertex
741  ?
742  Edge<dim>::forward
743  :
744  Edge<dim>::backward);
745 
746  // see if the opposite edge (there is only one in 2d) has already been
747  // oriented.
748  if (edges[opposite_edge].orientation_status == Edge<dim>::not_oriented)
749  {
750  // the opposite edge is not yet oriented. do orient it and add it to
751  // Delta_k
752  edges[opposite_edge].orientation_status = opposite_edge_orientation;
753  Delta_k.insert (opposite_edge);
754  }
755  else
756  {
757  // this opposite edge has already been oriented. it should be
758  // consistent with the current one in 2d, while in 3d it may in fact
759  // be mis-oriented, and in that case the mesh will not be
760  // orientable. indicate this by throwing an exception that we can
761  // catch further up; this has the advantage that we can propagate
762  // through a couple of functions without having to do error
763  // checking and without modifying the 'cells' array that the
764  // user gave us
765  if (dim == 2)
766  {
767  Assert (edges[opposite_edge].orientation_status == opposite_edge_orientation,
769  }
770  else if (dim == 3)
771  {
772  if (edges[opposite_edge].orientation_status != opposite_edge_orientation)
773  throw ExcMeshNotOrientable ();
774  }
775  else
776  Assert (false, ExcNotImplemented());
777  }
778  }
779  }
780  }
781 
782  // finally copy the new set to the previous one
783  // (corresponding to increasing 'k' by one in the
784  // algorithm)
785  Delta_k_minus_1 = Delta_k;
786  }
787  }
788 
789 
797  template <int dim>
798  void
799  rotate_cell (const std::vector<Cell<dim> > &cell_list,
800  const std::vector<Edge<dim> > &edge_list,
801  const unsigned int cell_index,
802  std::vector<CellData<dim> > &raw_cells)
803  {
804  // find the first vertex of the cell. this is the vertex where dim edges
805  // originate, so for each of the edges record which the starting vertex is
806  unsigned int starting_vertex_of_edge[GeometryInfo<dim>::lines_per_cell];
807  for (unsigned int e=0; e<GeometryInfo<dim>::lines_per_cell; ++e)
808  {
809  Assert (edge_list[cell_list[cell_index].edge_indices[e]].orientation_status
810  != Edge<dim>::not_oriented,
811  ExcInternalError());
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];
814  else
815  starting_vertex_of_edge[e] = edge_list[cell_list[cell_index].edge_indices[e]].vertex_indices[1];
816  }
817 
818  // find the vertex number that appears dim times. this will then be
819  // the vertex at which we want to locate the origin of the cell's
820  // coordinate system (i.e., vertex 0)
821  unsigned int origin_vertex_of_cell = numbers::invalid_unsigned_int;
822  switch (dim)
823  {
824  case 2:
825  {
826  // in 2d, we can simply enumerate the possibilities where the
827  // origin may be located because edges zero and one don't share
828  // any vertices, and the same for edges two and three
829  if ((starting_vertex_of_edge[0] == starting_vertex_of_edge[2])
830  ||
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])
834  ||
835  (starting_vertex_of_edge[1] == starting_vertex_of_edge[3]))
836  origin_vertex_of_cell = starting_vertex_of_edge[1];
837  else
838  Assert (false, ExcInternalError());
839 
840  break;
841  }
842 
843  case 3:
844  {
845  // one could probably do something similar in 3d, but that seems
846  // more complicated than one wants to write down. just go
847  // through the list of possible starting vertices and check
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],
852  &starting_vertex_of_edge[0]+GeometryInfo<dim>::lines_per_cell,
853  cell_list[cell_index].vertex_indices[origin_vertex_of_cell])
854  == dim)
855  break;
856  Assert (origin_vertex_of_cell < GeometryInfo<dim>::vertices_per_cell,
857  ExcInternalError());
858 
859  break;
860  }
861 
862  default:
863  Assert (false, ExcNotImplemented());
864  }
865 
866  // now rotate raw_cells[cell_index] in such a way that its orientation
867  // matches that of cell_list[cell_index]
868  switch (dim)
869  {
870  case 2:
871  {
872  // in 2d, we can literally rotate the cell until its origin
873  // matches the one that we have determined above should be
874  // the origin vertex
875  //
876  // when doing a rotation, take into account the ordering of
877  // vertices (not in clockwise or counter-clockwise sense)
878  while (raw_cells[cell_index].vertices[0] != origin_vertex_of_cell)
879  {
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;
885  }
886  break;
887  }
888 
889  case 3:
890  {
891  // in 3d, the situation is a bit more complicated. from above, we
892  // now know which vertex is at the origin (because 3 edges originate
893  // from it), but that still leaves 3 possible rotations of the cube.
894  // the important realization is that we can choose any of them:
895  // in all 3 rotations, all edges originate from the one vertex,
896  // and that fixes the directions of all 12 edges in the cube because
897  // these 3 cover all 3 equivalence classes! consequently, we can
898  // select an arbitrary one among the permutations -- for
899  // example the following ones:
900  static const unsigned int cube_permutations[8][8] =
901  {
902  {0,1,2,3,4,5,6,7},
903  {1,5,3,7,0,4,2,6},
904  {2,6,0,4,3,7,1,5},
905  {3,2,1,0,7,6,5,4},
906  {4,0,6,2,5,1,7,3},
907  {5,4,7,6,1,0,3,2},
908  {6,7,4,5,2,3,0,1},
909  {7,3,5,1,6,2,4,0}
910  };
911 
912  unsigned int temp_vertex_indices[GeometryInfo<dim>::vertices_per_cell];
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];
918 
919  break;
920  }
921 
922  default:
923  {
924  Assert (false, ExcNotImplemented());
925  }
926  }
927  }
928 
929 
935  template <int dim>
936  void reorient (std::vector<CellData<dim> > &cells)
937  {
938  // first build the arrays that connect cells to edges and the other
939  // way around
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);
942 
943  // then loop over all cells and start orienting parallel edge sets
944  // of cells that still have non-oriented edges
945  unsigned int next_cell_with_unoriented_edge = 0;
946  while ((next_cell_with_unoriented_edge = get_next_unoriented_cell(cell_list,
947  edge_list,
948  next_cell_with_unoriented_edge)) !=
950  {
951  // see which edge sets are still not oriented
952  //
953  // we do not need to look at each edge because if we orient edge
954  // 0, we will end up with edge 1 also oriented (in 2d; in 3d, there
955  // will be 3 other edges that are also oriented). there are only
956  // dim independent sets of edges, so loop over these.
957  //
958  // we need to check whether each one of these starter edges may
959  // already be oriented because the line (sheet) that connects
960  // globally parallel edges may be self-intersecting in the
961  // current cell
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,
966  edge_list,
967  next_cell_with_unoriented_edge,
968  ParallelEdges<dim>::starter_edges[l]);
969 
970  // ensure that we have really oriented all edges now, not just
971  // the starter edges
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,
975  ExcInternalError());
976  }
977 
978  // now that we have oriented all edges, we need to rotate cells
979  // so that the edges point in the right direction with the now
980  // rotated coordinate system
981  for (unsigned int c=0; c<cells.size(); ++c)
982  rotate_cell (cell_list, edge_list, c, cells);
983  }
984 
985 
986  // overload of the function above for 1d -- there is nothing
987  // to orient in that case
988  void reorient (std::vector<CellData<1> > &)
989  {}
990 }
991 
992 
993 // anonymous namespace for internal helper functions
994 namespace
995 {
1005  void
1006  reorder_new_to_old_style (std::vector<CellData<1> > &)
1007  {}
1008 
1009 
1010  void
1011  reorder_new_to_old_style (std::vector<CellData<2> > &cells)
1012  {
1013  for (unsigned int cell=0; cell<cells.size(); ++cell)
1014  std::swap(cells[cell].vertices[2], cells[cell].vertices[3]);
1015  }
1016 
1017 
1018  void
1019  reorder_new_to_old_style (std::vector<CellData<3> > &cells)
1020  {
1021  unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
1022  for (unsigned int cell=0; cell<cells.size(); ++cell)
1023  {
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)
1027  cells[cell].vertices[i] = tmp[GeometryInfo<3>::ucd_to_deal[i]];
1028  }
1029  }
1030 
1031 
1035  void
1036  reorder_old_to_new_style (std::vector<CellData<1> > &)
1037  {}
1038 
1039 
1040  void
1041  reorder_old_to_new_style (std::vector<CellData<2> > &cells)
1042  {
1043  // just invert the permutation:
1044  reorder_new_to_old_style(cells);
1045  }
1046 
1047 
1048  void
1049  reorder_old_to_new_style (std::vector<CellData<3> > &cells)
1050  {
1051  // undo the ordering above
1052  unsigned int tmp[GeometryInfo<3>::vertices_per_cell];
1053  for (unsigned int cell=0; cell<cells.size(); ++cell)
1054  {
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)
1058  cells[cell].vertices[GeometryInfo<3>::ucd_to_deal[i]] = tmp[i];
1059  }
1060  }
1061 }
1062 
1063 
1064 
1065 template <int dim, int spacedim>
1066 void
1068  const bool use_new_style_ordering)
1069 {
1070  Assert (cells.size() != 0,
1071  ExcMessage("List of elements to orient must have at least one cell"));
1072 
1073  // there is nothing for us to do in 1d
1074  if (dim == 1)
1075  return;
1076 
1077  // if necessary, convert to new-style format
1078  if (use_new_style_ordering == false)
1079  reorder_old_to_new_style(cells);
1080 
1081  // check if grids are already consistent. if so, do
1082  // nothing. if not, then do the reordering
1083  if (!is_consistent (cells))
1084  try
1085  {
1086  reorient(cells);
1087  }
1088  catch (const ExcMeshNotOrientable &)
1089  {
1090  // the mesh is not orientable. this is acceptable if we are in 3d,
1091  // as class Triangulation knows how to handle this, but it is
1092  // not in 2d; in that case, re-throw the exception
1093  if (dim < 3)
1094  throw;
1095  }
1096 
1097  // and convert back if necessary
1098  if (use_new_style_ordering == false)
1099  reorder_new_to_old_style(cells);
1100 }
1101 
1102 
1103 
1104 template<>
1105 void
1107  std::vector<CellData<1> > &)
1108 {
1109  // nothing to be done in 1d
1110 }
1111 
1112 
1113 
1114 template<>
1115 void
1117  std::vector<CellData<1> > &)
1118 {
1119  // nothing to be done in 1d
1120 }
1121 
1122 
1123 
1124 template<>
1125 void
1127  std::vector<CellData<1> > &)
1128 {
1129  // nothing to be done in 1d
1130 }
1131 
1132 
1133 template<>
1134 void
1135 GridReordering<2>::invert_all_cells_of_negative_grid(const std::vector<Point<2> > &all_vertices,
1136  std::vector<CellData<2> > &cells)
1137 {
1138  unsigned int vertices_lex[GeometryInfo<2>::vertices_per_cell];
1139  unsigned int n_negative_cells=0;
1140  for (unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1141  {
1142  // GridTools::cell_measure
1143  // requires the vertices to be
1144  // in lexicographic ordering
1145  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1146  vertices_lex[GeometryInfo<2>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1147  if (GridTools::cell_measure<2>(all_vertices, vertices_lex) < 0)
1148  {
1149  ++n_negative_cells;
1150  std::swap(cells[cell_no].vertices[1], cells[cell_no].vertices[3]);
1151 
1152  // check whether the
1153  // resulting cell is now ok.
1154  // if not, then the grid is
1155  // seriously broken and
1156  // should be sticked into the
1157  // bin
1158  for (unsigned int i=0; i<GeometryInfo<2>::vertices_per_cell; ++i)
1159  vertices_lex[GeometryInfo<2>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1160  AssertThrow(GridTools::cell_measure<2>(all_vertices, vertices_lex) > 0,
1161  ExcInternalError());
1162  }
1163  }
1164 
1165  // We assume that all cells of a grid have
1166  // either positive or negative volumes but
1167  // not both mixed. Although above reordering
1168  // might work also on single cells, grids
1169  // with both kind of cells are very likely to
1170  // be broken. Check for this here.
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 ")
1181  + Utilities::to_string (cells.size())
1182  + " cells in the mesh, "
1183  + Utilities::to_string (n_negative_cells)
1184  + " appear to have a negative volume."));
1185 }
1186 
1187 
1188 
1189 template<>
1190 void
1192  std::vector<CellData<2> > &)
1193 {
1194  Assert(false, ExcNotImplemented());
1195 }
1196 
1197 
1198 
1199 template<>
1200 void
1201 GridReordering<3>::invert_all_cells_of_negative_grid(const std::vector<Point<3> > &all_vertices,
1202  std::vector<CellData<3> > &cells)
1203 {
1204  unsigned int vertices_lex[GeometryInfo<3>::vertices_per_cell];
1205  unsigned int n_negative_cells=0;
1206  for (unsigned int cell_no=0; cell_no<cells.size(); ++cell_no)
1207  {
1208  // GridTools::cell_measure
1209  // requires the vertices to be
1210  // in lexicographic ordering
1211  for (unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1212  vertices_lex[GeometryInfo<3>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1213  if (GridTools::cell_measure<3>(all_vertices, vertices_lex) < 0)
1214  {
1215  ++n_negative_cells;
1216  // reorder vertices: swap front and back face
1217  for (unsigned int i=0; i<4; ++i)
1218  std::swap(cells[cell_no].vertices[i], cells[cell_no].vertices[i+4]);
1219 
1220  // check whether the
1221  // resulting cell is now ok.
1222  // if not, then the grid is
1223  // seriously broken and
1224  // should be sticked into the
1225  // bin
1226  for (unsigned int i=0; i<GeometryInfo<3>::vertices_per_cell; ++i)
1227  vertices_lex[GeometryInfo<3>::ucd_to_deal[i]]=cells[cell_no].vertices[i];
1228  AssertThrow(GridTools::cell_measure<3>(all_vertices, vertices_lex) > 0,
1229  ExcInternalError());
1230  }
1231  }
1232 
1233  // We assume that all cells of a
1234  // grid have either positive or
1235  // negative volumes but not both
1236  // mixed. Although above reordering
1237  // might work also on single cells,
1238  // grids with both kind of cells
1239  // are very likely to be
1240  // broken. Check for this here.
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."));
1260 }
1261 
1262 
1263 
1264 /* ------------------------ explicit instantiations ------------------- */
1265 template class GridReordering<1,1>;
1266 template class GridReordering<1,2>;
1267 template class GridReordering<1,3>;
1268 template class GridReordering<2,2>;
1269 template class GridReordering<2,3>;
1270 template class GridReordering<3,3>;
1271 
1272 DEAL_II_NAMESPACE_CLOSE
1273 
static const unsigned int invalid_unsigned_int
Definition: types.h:170
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:369
Definition: point.h:89
std::string to_string(const number value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:92
static ::ExceptionBase & ExcMessage(std::string arg1)
#define Assert(cond, exc)
Definition: exceptions.h:313
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]
Definition: tria.h:139
static ::ExceptionBase & ExcInternalError()