]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
A new test for GridTools::collect_periodic_face_pairs
authormaier <maier@0785d39b-7218-0410-832d-ea1e28bc413d>
Sat, 1 Dec 2012 17:17:02 +0000 (17:17 +0000)
committermaier <maier@0785d39b-7218-0410-832d-ea1e28bc413d>
Sat, 1 Dec 2012 17:17:02 +0000 (17:17 +0000)
git-svn-id: https://svn.dealii.org/trunk@27716 0785d39b-7218-0410-832d-ea1e28bc413d

tests/deal.II/grid_tools_05.cc
tests/deal.II/grid_tools_05/cmp/generic

index f403dc6691b11fbc8f1a201f7dff4a38063842dc..701eda7ebd2a7a41af0d5ad1646022d5dd6ffce9 100644 (file)
 //----------------------------  grid_tools.cc  ---------------------------
 
 
-// check GridTools::collect_periodic_cell_pairs for some simple cases
+//
+// check collect_periodic_face_pairs for correct return values
+//
 
 
 #include "../tests.h"
 #include <deal.II/base/logstream.h>
 
-#include <deal.II/grid/tria.h>
 #include <deal.II/grid/grid_generator.h>
 #include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
 
 #include <fstream>
-
 std::ofstream logfile("grid_tools_05/output");
 
 using namespace dealii;
 
-template<int dim>
-void test1 ()
+
+/*
+ * Generate a grid consisting of two disjoint cells, colorize the two
+ * outermost faces. They will be matched via collect_periodic_face_pairs
+ *
+ * The integer orientation determines the orientation of the second cell
+ * (to get something else than the boring default orientation).
+ */
+
+/* The 2D case */
+void generate_grid(Triangulation<2> &triangulation, int orientation)
 {
-  Triangulation<dim> triangulation;
+    Point<2> vertices_1[]
+      = {
+           Point<2> (-1.,-3.),
+           Point<2> (+1.,-3.),
+           Point<2> (-1.,-1.),
+           Point<2> (+1.,-1.),
+           Point<2> (-1.,+1.),
+           Point<2> (+1.,+1.),
+           Point<2> (-1.,+3.),
+           Point<2> (+1.,+3.),
+        };
+    std::vector<Point<2> > vertices (&vertices_1[0], &vertices_1[8]);
 
-  GridGenerator::hyper_cube(triangulation, 0., 1.);
-  triangulation.refine_global(2);
+    std::vector<CellData<2> > cells (2, CellData<2>());
 
-  deallog << std::endl << std::endl << "Hyper cube, dimension " << dim << ":" << std::endl;
+    /* cell 0 */
+    int cell_vertices_0[GeometryInfo<2>::vertices_per_cell] = {0, 1,  2,  3};
 
-  typedef typename Triangulation<dim>::cell_iterator CELL_IT;
-  typedef typename std::map<CELL_IT, CELL_IT> MAP;
+    /* cell 1 */
+    int cell_vertices_1[2][GeometryInfo<2>::vertices_per_cell]
+      = {
+        {4,5,6,7},
+        {7,6,5,4},
+      };
 
-  MAP matched_cells = GridTools::collect_periodic_cell_pairs(triangulation, 0, /* direction: */1);
+    for (unsigned int j=0;j<GeometryInfo<2>::vertices_per_cell;++j) {
+      cells[0].vertices[j] = cell_vertices_0[j];
+      cells[1].vertices[j] = cell_vertices_1[orientation][j];
+    }
+    cells[0].material_id = 0;
+    cells[1].material_id = 0;
 
-  deallog << "Matching:" << std::endl;
-  for (typename MAP::const_iterator it = matched_cells.begin(); it != matched_cells.end(); ++it) {
-    deallog << it->first << "  -  " << it->second << std::endl;
-  }
+    triangulation.create_triangulation(vertices, cells, SubCellData());
+
+    Triangulation<2>::cell_iterator cell_1 = triangulation.begin();
+    Triangulation<2>::cell_iterator cell_2 = cell_1++;
+    Triangulation<2>::face_iterator face_1;
+    Triangulation<2>::face_iterator face_2;
+
+    // Look for the two outermost faces:
+    for (unsigned int j=0;j<GeometryInfo<2>::faces_per_cell;++j) {
+      if (cell_1->face(j)->center()(1) > 2.9)
+        face_1 = cell_1->face(j);
+      if (cell_2->face(j)->center()(1) < -2.9)
+        face_2 = cell_2->face(j);
+    }
+    face_1->set_boundary_indicator(42);
+    face_2->set_boundary_indicator(43);
+
+    triangulation.refine_global(1);
 }
 
-template<int dim>
-void test2 ()
+
+/* The 3D case */
+void generate_grid(Triangulation<3> &triangulation, int orientation)
 {
-  /* Unfortunately, there is no parallelogram for 3D.. */
-  Assert(dim == 2,
-         ExcNotImplemented());
+    Point<3> vertices_1[]
+      = {
+           Point<3> (-1.,-1.,-3.),
+           Point<3> (+1.,-1.,-3.),
+           Point<3> (-1.,+1.,-3.),
+           Point<3> (+1.,+1.,-3.),
+           Point<3> (-1.,-1.,-1.),
+           Point<3> (+1.,-1.,-1.),
+           Point<3> (-1.,+1.,-1.),
+           Point<3> (+1.,+1.,-1.),
+           Point<3> (-1.,-1.,+1.),
+           Point<3> (+1.,-1.,+1.),
+           Point<3> (-1.,+1.,+1.),
+           Point<3> (+1.,+1.,+1.),
+           Point<3> (-1.,-1.,+3.),
+           Point<3> (+1.,-1.,+3.),
+           Point<3> (-1.,+1.,+3.),
+           Point<3> (+1.,+1.,+3.)
+        };
+    std::vector<Point<3> > vertices (&vertices_1[0], &vertices_1[16]);
 
-  Triangulation<dim> triangulation;
+    std::vector<CellData<3> > cells (2, CellData<3>());
 
-  Tensor<2,dim> vectors;
-  vectors[0][0] = 1.;
-  vectors[0][1] = 0.;
-  vectors[1][0] = 0.5;
-  vectors[1][1] = 1.;
+    /* cell 0 */
+    int cell_vertices_0[GeometryInfo<3>::vertices_per_cell] = {0, 1,  2,  3,  4,  5,  6,  7};
 
-  Tensor<1,dim> offset;
-  offset[0] = 0.5;
-  offset[1] = 1.;
+    /* cell 1 */
+    int cell_vertices_1[8][GeometryInfo<3>::vertices_per_cell]
+      = {
+        {8,9,10,11,12,13,14,15},
+        {9,11,8,10,13,15,12,14},
+        {11,10,9,8,15,14,13,12},
+        {10,8,11,9,14,12,15,13},
+        {13,12,15,14,9,8,11,10},
+        {12,14,13,15,8,10,9,11},
+        {14,15,12,13,10,11,8,9},
+        {15,13,14,12,11,9,10,8},
+      };
 
-  GridGenerator::parallelogram(triangulation, vectors);
-  triangulation.refine_global(3);
+    for (unsigned int j=0;j<GeometryInfo<3>::vertices_per_cell;++j) {
+      cells[0].vertices[j] = cell_vertices_0[j];
+      cells[1].vertices[j] = cell_vertices_1[orientation][j];
+    }
+    cells[0].material_id = 0;
+    cells[1].material_id = 0;
 
-  deallog << std::endl << std::endl << "Parallelogram, dimension " << dim << ":" << std::endl;
 
+    triangulation.create_triangulation(vertices, cells, SubCellData());
 
-  typedef typename Triangulation<dim>::cell_iterator CELL_IT;
-  typedef typename std::map<CELL_IT, CELL_IT> MAP;
+    Triangulation<3>::cell_iterator cell_1 = triangulation.begin();
+    Triangulation<3>::cell_iterator cell_2 = cell_1++;
+    Triangulation<3>::face_iterator face_1;
+    Triangulation<3>::face_iterator face_2;
 
-  MAP matched_cells = GridTools::collect_periodic_cell_pairs(triangulation, 0, /* direction: */0);
-  deallog << "Matching, direction 0:" << std::endl;
-  for (typename MAP::const_iterator it = matched_cells.begin(); it != matched_cells.end(); ++it) {
-    deallog << it->first << "  -  " << it->second << std::endl;
-  }
+    // Look for the two outermost faces:
+    for (unsigned int j=0;j<GeometryInfo<3>::faces_per_cell;++j) {
+      if (cell_1->face(j)->center()(2) > 2.9)
+        face_1 = cell_1->face(j);
+      if (cell_2->face(j)->center()(2) < -2.9)
+        face_2 = cell_2->face(j);
+    }
+    face_1->set_boundary_indicator(42);
+    face_2->set_boundary_indicator(43);
 
-  MAP matched_cells_2 = GridTools::collect_periodic_cell_pairs(triangulation, 0, /* direction: */1, offset);
-  deallog << "Matching, direction 1:" << std::endl;
-  for (typename MAP::const_iterator it = matched_cells_2.begin(); it != matched_cells_2.end(); ++it) {
-    deallog << it->first << "  -  " << it->second << std::endl;
-  }
+    triangulation.refine_global(1);
+}
+
+
+
+/*
+ * Print out the face vertices as well as the orientation of a match:
+ */
+template<typename FaceIterator>
+void print_match(const FaceIterator &face_1,
+                 const FaceIterator &face_2,
+                 const std::bitset<3> &orientation)
+{
+  static const int dim = FaceIterator::AccessorType::dimension;
+
+  deallog << "face 1";
+  for (unsigned int j=0;j<GeometryInfo<dim>::vertices_per_face;++j)
+    deallog << " :: " << face_1->vertex(j);
+  deallog << std::endl;
+
+  deallog << "face 2";
+  for (unsigned int j=0;j<GeometryInfo<dim>::vertices_per_face;++j)
+    deallog << " :: " << face_2->vertex(j);
+  deallog << std::endl;
+
+  deallog << "orientation: " << orientation[0]
+            << "  flip: " << orientation[1]
+            << "  rotation: " << orientation[2]
+            << std::endl
+            << std::endl;
 }
 
 int main()
@@ -98,11 +201,73 @@ int main()
   deallog.depth_console(0);
   deallog.threshold_double(1.e-10);
 
-  /* Let's try to match the hyper_cube in 3D:*/
-  test1<3> ();
+  deallog << "Test for 2D: Hypercube" << std::endl << std::endl;
+
+  for(int i = 0; i < 2; ++i) {
+    // Generate a triangulation and match:
+    Triangulation<2> triangulation;
+
+    generate_grid(triangulation, i);
+
+    typedef Triangulation<2>::face_iterator FaceIterator;
+    typedef std::map<FaceIterator, std::pair<FaceIterator, std::bitset<3> > > FaceMap;
+    FaceMap test = GridTools::collect_periodic_face_pairs (triangulation, 42, 43, 1,
+                                                           dealii::Tensor<1,2>());
+
+    deallog << "Triangulation: " << i << std::endl;
+    deallog << "Coarse match: " << test.size() << std::endl;
+
+    for (FaceMap::iterator facepair = test.begin(); facepair != test.end(); ++facepair)
+      print_match(facepair->first, facepair->second.first, facepair->second.second);
+
+
+    typedef Triangulation<2>::active_face_iterator FaceIterator2;
+    typedef std::map<FaceIterator2, std::pair<FaceIterator2, std::bitset<3> > > FaceMap2;
+    FaceMap2 test2 = GridTools::collect_periodic_face_pairs (triangulation.begin_active_face(),
+                                                            triangulation.end_face(),
+                                                            42, 43,
+                                                            1,
+                                                            dealii::Tensor<1,2>());
+
+    deallog << "Fine match: " << test2.size() << std::endl;
+
+    for (FaceMap2::iterator facepair = test2.begin(); facepair != test2.end(); ++facepair)
+      print_match(facepair->first, facepair->second.first, facepair->second.second);
+  }
 
-  /* Let's try to match a parallelogramm in 2D */
-  test2<2> ();
+  deallog << "Test for 3D: Hypercube" << std::endl << std::endl;
+
+  for(int i = 0; i < 8; ++i) {
+    // Generate a triangulation and match:
+    Triangulation<3> triangulation;
+
+    generate_grid(triangulation, i);
+
+    typedef Triangulation<3>::face_iterator FaceIterator;
+    typedef std::map<FaceIterator, std::pair<FaceIterator, std::bitset<3> > > FaceMap;
+    FaceMap test = GridTools::collect_periodic_face_pairs (triangulation, 42, 43, 2,
+                                                           dealii::Tensor<1,3>());
+
+    deallog << "Triangulation: " << i << std::endl;
+    deallog << "Coarse match: " << test.size() << std::endl;
+
+    for (FaceMap::iterator facepair = test.begin(); facepair != test.end(); ++facepair)
+      print_match(facepair->first, facepair->second.first, facepair->second.second);
+
+
+    typedef Triangulation<3>::active_face_iterator FaceIterator2;
+    typedef std::map<FaceIterator2, std::pair<FaceIterator2, std::bitset<3> > > FaceMap2;
+    FaceMap2 test2 = GridTools::collect_periodic_face_pairs (triangulation.begin_active_face(),
+                                                            triangulation.end_face(),
+                                                            42, 43,
+                                                            2,
+                                                            dealii::Tensor<1,3>());
+
+    deallog << "Fine match: " << test2.size() << std::endl;
+
+    for (FaceMap2::iterator facepair = test2.begin(); facepair != test2.end(); ++facepair)
+      print_match(facepair->first, facepair->second.first, facepair->second.second);
+  }
 
   return 0;
 }
index 677327fbb2cb7c51fa11e5c6f12a3ae41f4cde1e..36b7388a9cd9ca71fee02150fc5e444762d28e7c 100644 (file)
 
+DEAL::Test for 2D: Hypercube
 DEAL::
+DEAL::Triangulation: 0
+DEAL::Coarse match: 1
+DEAL::face 1 :: -1.000 3.000 :: 1.000 3.000
+DEAL::face 2 :: -1.000 -3.000 :: 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::Fine match: 2
+DEAL::face 1 :: -1.000 3.000 :: 0.000 3.000
+DEAL::face 2 :: -1.000 -3.000 :: 0.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 3.000 :: 1.000 3.000
+DEAL::face 2 :: 0.000 -3.000 :: 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::Triangulation: 1
+DEAL::Coarse match: 1
+DEAL::face 1 :: 1.000 3.000 :: -1.000 3.000
+DEAL::face 2 :: -1.000 -3.000 :: 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::Fine match: 2
+DEAL::face 1 :: 1.000 3.000 :: 0.000 3.000
+DEAL::face 2 :: 0.000 -3.000 :: 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 3.000 :: -1.000 3.000
+DEAL::face 2 :: -1.000 -3.000 :: 0.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::Test for 3D: Hypercube
+DEAL::
+DEAL::Triangulation: 0
+DEAL::Coarse match: 1
+DEAL::face 1 :: -1.000 -1.000 3.000 :: 1.000 -1.000 3.000 :: -1.000 1.000 3.000 :: 1.000 1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: -1.000 -1.000 3.000 :: 0.000 -1.000 3.000 :: -1.000 0.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 -1.000 3.000 :: 1.000 -1.000 3.000 :: 0.000 0.000 3.000 :: 1.000 0.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: -1.000 0.000 3.000 :: 0.000 0.000 3.000 :: -1.000 1.000 3.000 :: 0.000 1.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: 1.000 0.000 3.000 :: 0.000 1.000 3.000 :: 1.000 1.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 0
+DEAL::
+DEAL::Triangulation: 1
+DEAL::Coarse match: 1
+DEAL::face 1 :: 1.000 -1.000 3.000 :: 1.000 1.000 3.000 :: -1.000 -1.000 3.000 :: -1.000 1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 1
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: 1.000 -1.000 3.000 :: 1.000 0.000 3.000 :: 0.000 -1.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 1
+DEAL::
+DEAL::face 1 :: 1.000 0.000 3.000 :: 1.000 1.000 3.000 :: 0.000 0.000 3.000 :: 0.000 1.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 -1.000 3.000 :: 0.000 0.000 3.000 :: -1.000 -1.000 3.000 :: -1.000 0.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: 0.000 1.000 3.000 :: -1.000 0.000 3.000 :: -1.000 1.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 1  flip: 0  rotation: 1
+DEAL::
+DEAL::Triangulation: 2
+DEAL::Coarse match: 1
+DEAL::face 1 :: 1.000 1.000 3.000 :: -1.000 1.000 3.000 :: 1.000 -1.000 3.000 :: -1.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: 1.000 1.000 3.000 :: 0.000 1.000 3.000 :: 1.000 0.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 1.000 3.000 :: -1.000 1.000 3.000 :: 0.000 0.000 3.000 :: -1.000 0.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 1.000 0.000 3.000 :: 0.000 0.000 3.000 :: 1.000 -1.000 3.000 :: 0.000 -1.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: -1.000 0.000 3.000 :: 0.000 -1.000 3.000 :: -1.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 0
+DEAL::
+DEAL::Triangulation: 3
+DEAL::Coarse match: 1
+DEAL::face 1 :: -1.000 1.000 3.000 :: -1.000 -1.000 3.000 :: 1.000 1.000 3.000 :: 1.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 1
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: -1.000 1.000 3.000 :: -1.000 0.000 3.000 :: 0.000 1.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 1
+DEAL::
+DEAL::face 1 :: -1.000 0.000 3.000 :: -1.000 -1.000 3.000 :: 0.000 0.000 3.000 :: 0.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 1.000 3.000 :: 0.000 0.000 3.000 :: 1.000 1.000 3.000 :: 1.000 0.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: 0.000 -1.000 3.000 :: 1.000 0.000 3.000 :: 1.000 -1.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 1  flip: 1  rotation: 1
+DEAL::
+DEAL::Triangulation: 4
+DEAL::Coarse match: 1
+DEAL::face 1 :: 1.000 -1.000 3.000 :: -1.000 -1.000 3.000 :: 1.000 1.000 3.000 :: -1.000 1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 1
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: 1.000 -1.000 3.000 :: 0.000 -1.000 3.000 :: 1.000 0.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 -1.000 3.000 :: -1.000 -1.000 3.000 :: 0.000 0.000 3.000 :: -1.000 0.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 1
+DEAL::
+DEAL::face 1 :: 1.000 0.000 3.000 :: 0.000 0.000 3.000 :: 1.000 1.000 3.000 :: 0.000 1.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: -1.000 0.000 3.000 :: 0.000 1.000 3.000 :: -1.000 1.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 1
+DEAL::
+DEAL::Triangulation: 5
+DEAL::Coarse match: 1
+DEAL::face 1 :: -1.000 -1.000 3.000 :: -1.000 1.000 3.000 :: 1.000 -1.000 3.000 :: 1.000 1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 0
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: -1.000 -1.000 3.000 :: -1.000 0.000 3.000 :: 0.000 -1.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: -1.000 0.000 3.000 :: -1.000 1.000 3.000 :: 0.000 0.000 3.000 :: 0.000 1.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 -1.000 3.000 :: 0.000 0.000 3.000 :: 1.000 -1.000 3.000 :: 1.000 0.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: 0.000 1.000 3.000 :: 1.000 0.000 3.000 :: 1.000 1.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 0
+DEAL::
+DEAL::Triangulation: 6
+DEAL::Coarse match: 1
+DEAL::face 1 :: -1.000 1.000 3.000 :: 1.000 1.000 3.000 :: -1.000 -1.000 3.000 :: 1.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 1
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: -1.000 1.000 3.000 :: 0.000 1.000 3.000 :: -1.000 0.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 1.000 3.000 :: 1.000 1.000 3.000 :: 0.000 0.000 3.000 :: 1.000 0.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 1
+DEAL::
+DEAL::face 1 :: -1.000 0.000 3.000 :: 0.000 0.000 3.000 :: -1.000 -1.000 3.000 :: 0.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 1
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: 1.000 0.000 3.000 :: 0.000 -1.000 3.000 :: 1.000 -1.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 0  flip: 0  rotation: 1
+DEAL::
+DEAL::Triangulation: 7
+DEAL::Coarse match: 1
+DEAL::face 1 :: 1.000 1.000 3.000 :: 1.000 -1.000 3.000 :: -1.000 1.000 3.000 :: -1.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: -1.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 0
+DEAL::
+DEAL::Fine match: 4
+DEAL::face 1 :: 1.000 1.000 3.000 :: 1.000 0.000 3.000 :: 0.000 1.000 3.000 :: 0.000 0.000 3.000
+DEAL::face 2 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000 :: 0.000 1.000 -3.000 :: 1.000 1.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 1.000 0.000 3.000 :: 1.000 -1.000 3.000 :: 0.000 0.000 3.000 :: 0.000 -1.000 3.000
+DEAL::face 2 :: 0.000 -1.000 -3.000 :: 1.000 -1.000 -3.000 :: 0.000 0.000 -3.000 :: 1.000 0.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 1.000 3.000 :: 0.000 0.000 3.000 :: -1.000 1.000 3.000 :: -1.000 0.000 3.000
+DEAL::face 2 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000 :: -1.000 1.000 -3.000 :: 0.000 1.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 0
+DEAL::
+DEAL::face 1 :: 0.000 0.000 3.000 :: 0.000 -1.000 3.000 :: -1.000 0.000 3.000 :: -1.000 -1.000 3.000
+DEAL::face 2 :: -1.000 -1.000 -3.000 :: 0.000 -1.000 -3.000 :: -1.000 0.000 -3.000 :: 0.000 0.000 -3.000
+DEAL::orientation: 0  flip: 1  rotation: 0
 DEAL::
-DEAL::Hyper cube, dimension 3:
-DEAL::Matching:
-DEAL::2.0  -  2.18
-DEAL::2.1  -  2.19
-DEAL::2.4  -  2.22
-DEAL::2.5  -  2.23
-DEAL::2.8  -  2.26
-DEAL::2.9  -  2.27
-DEAL::2.12  -  2.30
-DEAL::2.13  -  2.31
-DEAL::2.32  -  2.50
-DEAL::2.33  -  2.51
-DEAL::2.36  -  2.54
-DEAL::2.37  -  2.55
-DEAL::2.40  -  2.58
-DEAL::2.41  -  2.59
-DEAL::2.44  -  2.62
-DEAL::2.45  -  2.63
-DEAL::
-DEAL::
-DEAL::Parallelogram, dimension 2:
-DEAL::Matching, direction 0:
-DEAL::3.0  -  3.21
-DEAL::3.2  -  3.23
-DEAL::3.8  -  3.29
-DEAL::3.10  -  3.31
-DEAL::3.32  -  3.53
-DEAL::3.34  -  3.55
-DEAL::3.40  -  3.61
-DEAL::3.42  -  3.63
-DEAL::Matching, direction 1:
-DEAL::3.0  -  3.42
-DEAL::3.1  -  3.43
-DEAL::3.4  -  3.46
-DEAL::3.5  -  3.47
-DEAL::3.16  -  3.58
-DEAL::3.17  -  3.59
-DEAL::3.20  -  3.62
-DEAL::3.21  -  3.63

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.