]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Fix computation of face normal vectors in the codim=1 case.
authorWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 20 Dec 2010 03:51:23 +0000 (03:51 +0000)
committerWolfgang Bangerth <bangerth@math.tamu.edu>
Mon, 20 Dec 2010 03:51:23 +0000 (03:51 +0000)
git-svn-id: https://svn.dealii.org/trunk@23015 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/source/fe/mapping_q1.cc
tests/codim_one/normal_vectors_01.cc
tests/codim_one/normal_vectors_01/cmp/generic

index 200ed752ac4bbcb54abe237b68e7f64ca9a86316..34f474b93b20675b6426462f1d7596c19174f47e 100644 (file)
@@ -1016,7 +1016,8 @@ namespace internal
            Assert (JxW_values.size() == n_q_points,
                    ExcDimensionMismatch(JxW_values.size(), n_q_points));
 
-                                          // map the unit tangentials to the real cell
+                                          // map the unit tangentials to the
+                                          // real cell
          for (unsigned int d=0; d<dim-1; ++d)
            {
              Assert (face_no+GeometryInfo<dim>::faces_per_cell*d <
@@ -1032,35 +1033,81 @@ namespace internal
                                 mapping_contravariant);
            }
 
-         switch (dim)
-           {
-             case 2:
-             {
-               for (unsigned int i=0; i<n_q_points; ++i)
+                                          // if dim==spacedim, we can use the
+                                          // unit tangentials to compute the
+                                          // boundary form by simply taking
+                                          // the cross product
+         if (dim == spacedim)
+           {
+             for (unsigned int i=0; i<n_q_points; ++i)
+               if (dim == 2)
                  cross_product (boundary_forms[i], data.aux[0][i]);
-               break;
-             }
-
-             case 3:
-             {
-               for (unsigned int i=0; i<n_q_points; ++i)
+               else if (dim == 3)
                  cross_product (boundary_forms[i], data.aux[0][i], data.aux[1][i]);
-               break;
-             }
-
-             default:
-                   Assert(false, ExcNotImplemented());
+               else
+                 Assert(false, ExcNotImplemented());
            }
+         else
+           {
+                                              // in the codim-one case, the
+                                              // boundary form results from
+                                              // the cross product of all the
+                                              // face tangential vectors and
+                                              // the cell normal vector
+                                              //
+                                              // to compute the cell normal,
+                                              // use the same method used in
+                                              // fill_fe_values for cells
+                                              // above
+             Assert (data.contravariant.size() == n_q_points,
+                     ExcInternalError());
+             for (unsigned int point=0; point<n_q_points; ++point)
+               {
+                 Tensor<1,spacedim> cell_normal;
+
+                 data.contravariant[point] = transpose(data.contravariant[point]);
+
+                                                  // compute the normal
+                                                  // vector to this cell
+                                                  // and put it into the
+                                                  // last row of
+                                                  // data.contravariant
+                 if ( (dim==1) && (spacedim==2) )
+                   cross_product(cell_normal,
+                                 -data.contravariant[point][0]);
+                 else if ( (dim==2) && (spacedim==3) )
+                   cross_product(cell_normal,
+                                 data.contravariant[point][0],
+                                 data.contravariant[point][1]);
+                 else
+                   Assert (false, ExcNotImplemented());
+
+                                                  // then compute the face
+                                                  // normal from the face
+                                                  // tangent and the cell
+                                                  // normal:
+                 if ( (dim==1) && (spacedim==2) )
+                                                    // need to think how to
+                                                    // figure out whether we
+                                                    // need to point to the
+                                                    // left or right
+                   Assert (false, ExcNotImplemented())
+                 else if ( (dim==2) && (spacedim==3) )
+                   cross_product (boundary_forms[point],
+                                  data.aux[0][point], cell_normal);
+                 else
+                   Assert (false, ExcNotImplemented());
+               }
+           }
+         
 
          if (update_flags & (update_normal_vectors
                              | update_JxW_values))
            for (unsigned int i=0;i<boundary_forms.size();++i)
              {
-               const double f = std::sqrt(contract(boundary_forms[i],
-                                                   boundary_forms[i]));
                if (update_flags & update_JxW_values)
                  {
-                   JxW_values[i] = f * weights[i];
+                   JxW_values[i] = boundary_forms[i].norm() * weights[i];
                    if (subface_no!=deal_II_numbers::invalid_unsigned_int)
                      {
                        const double area_ratio=GeometryInfo<dim>::subface_ratio(
@@ -1070,7 +1117,7 @@ namespace internal
                  }
 
                if (update_flags & update_normal_vectors)
-                 normal_vectors[i] = boundary_forms[i] / f;
+                 normal_vectors[i] = boundary_forms[i] / boundary_forms[i].norm();
              }
        }
     }
index 7b6a0160cf8f0475b63da5fb867fe12424fb0c36..84cda10e765eecfe1374bf24686861c617100065 100644 (file)
@@ -11,9 +11,8 @@
 //
 //----------------------------  normal_vectors_01.cc  ---------------------------
 
-/*
-  Asking for normal vectors in the codim-1 case led to aborts.
-*/
+//  Asking for face normal vectors in the codim-1 case led to aborts.
+
 
 #include "../tests.h"
 
 template <int dim>
 void test ()
 {
-  Triangulation<dim-1,dim> tria;
-  std::map<typename Triangulation<dim-1,dim>::cell_iterator,
-    typename Triangulation<dim,dim>::face_iterator> surface_to_volume_mapping;
-
-  HyperBallBoundary<dim> boundary_description;
-  Triangulation<dim> volume_mesh;
-  GridGenerator::half_hyper_ball(volume_mesh);
-  
-  volume_mesh.set_boundary (1, boundary_description);
-  volume_mesh.set_boundary (0, boundary_description);
-  volume_mesh.refine_global (1);
-  
-  static HyperBallBoundary<dim-1,dim> surface_description;
-  tria.set_boundary (1, surface_description);
-  tria.set_boundary (0, surface_description);
-  
-  std::set<unsigned char> boundary_ids;
-  boundary_ids.insert(0);
+  Triangulation<dim,dim> volume_mesh;
+  GridGenerator::hyper_cube(volume_mesh);
   
-  surface_to_volume_mapping
-    = GridTools::extract_boundary_mesh (volume_mesh, tria,
-                                       boundary_ids);
+  Triangulation<dim-1,dim> tria;
+  GridTools::extract_boundary_mesh (volume_mesh, tria);
   
   FE_Q<dim-1,dim> fe (1);
   DoFHandler<dim-1,dim> dh(tria);
@@ -63,14 +45,27 @@ void test ()
   
   dh.distribute_dofs (fe);
   
-  FEFaceValues<dim-1,dim> fe_face_values (mapping, fe, QGauss<dim-2>(2), 
-                                         update_values         | update_quadrature_points  |
-                                         update_normal_vectors | update_JxW_values);
-  
-  typename DoFHandler<dim-1,dim>::active_cell_iterator
-    cell = dh.begin_active();
-  
-  fe_face_values.reinit (cell,0);
+  FEFaceValues<dim-1,dim> fe_face_values (mapping, fe, QTrapez<dim-2>(),
+                                         update_normal_vectors);
+
+  for (typename DoFHandler<dim-1,dim>::active_cell_iterator
+        cell = dh.begin_active();
+       cell != dh.end(); ++cell)
+    {
+      deallog << "Face centered at " << cell->center()
+             << std::endl;
+
+      for (unsigned int f=0; f<GeometryInfo<dim-1>::faces_per_cell; ++f)
+       {
+         deallog << "  Edge centered at " << cell->face(f)->center()
+                 << std::endl;
+      
+         fe_face_values.reinit (cell,f);
+         for (unsigned int q=0; q<fe_face_values.n_quadrature_points; ++q)
+           deallog << "    normal_vector=" << fe_face_values.normal_vector(q)
+                   << std::endl;
+       }
+    }
 }
 
 
index d8b5693c0e6fe8399e1e6496c2ca3f3e21f84573..f65a59f7561abccaf4f41f97baa3f2f9e0058c5e 100644 (file)
 
-DEAL::Testing hyper_cube in dim: 3...
-DEAL::The last column should be 0
-DEAL::Vol faces                Surf cell               Distance
-DEAL::
-DEAL::-0.577350 -0.577350 -0.577350            -0.577350 -0.577350 -0.577350                   0
-DEAL::0.00000 -0.577350 -0.577350              0.00000 -0.577350 -0.577350                     0
-DEAL::-0.577350 0.00000 -0.577350              -0.577350 0.00000 -0.577350                     0
-DEAL::0.00000 0.00000 -0.577350                0.00000 0.00000 -0.577350                       0
-DEAL::
-DEAL::0.00000 -0.577350 -0.577350              0.00000 -0.577350 -0.577350                     0
-DEAL::0.577350 -0.577350 -0.577350             0.577350 -0.577350 -0.577350                    0
-DEAL::0.00000 0.00000 -0.577350                0.00000 0.00000 -0.577350                       0
-DEAL::0.577350 0.00000 -0.577350               0.577350 0.00000 -0.577350                      0
-DEAL::
-DEAL::-0.577350 0.00000 -0.577350              -0.577350 0.00000 -0.577350                     0
-DEAL::0.00000 0.00000 -0.577350                0.00000 0.00000 -0.577350                       0
-DEAL::-0.577350 0.577350 -0.577350             -0.577350 0.577350 -0.577350                    0
-DEAL::0.00000 0.577350 -0.577350               0.00000 0.577350 -0.577350                      0
-DEAL::
-DEAL::0.00000 0.00000 -0.577350                0.00000 0.00000 -0.577350                       0
-DEAL::0.577350 0.00000 -0.577350               0.577350 0.00000 -0.577350                      0
-DEAL::0.00000 0.577350 -0.577350               0.00000 0.577350 -0.577350                      0
-DEAL::0.577350 0.577350 -0.577350              0.577350 0.577350 -0.577350                     0
-DEAL::
-DEAL::0.577350 -0.577350 -0.577350             0.577350 -0.577350 -0.577350                    0
-DEAL::0.577350 -0.577350 0.00000               0.577350 -0.577350 0.00000                      0
-DEAL::0.577350 0.00000 -0.577350               0.577350 0.00000 -0.577350                      0
-DEAL::0.577350 0.00000 0.00000         0.577350 0.00000 0.00000                        0
-DEAL::
-DEAL::0.577350 -0.577350 0.00000               0.577350 -0.577350 0.00000                      0
-DEAL::0.577350 -0.577350 0.577350              0.577350 -0.577350 0.577350                     0
-DEAL::0.577350 0.00000 0.00000         0.577350 0.00000 0.00000                        0
-DEAL::0.577350 0.00000 0.577350                0.577350 0.00000 0.577350                       0
-DEAL::
-DEAL::0.577350 0.00000 -0.577350               0.577350 0.00000 -0.577350                      0
-DEAL::0.577350 0.00000 0.00000         0.577350 0.00000 0.00000                        0
-DEAL::0.577350 0.577350 -0.577350              0.577350 0.577350 -0.577350                     0
-DEAL::0.577350 0.577350 0.00000                0.577350 0.577350 0.00000                       0
-DEAL::
-DEAL::0.577350 0.00000 0.00000         0.577350 0.00000 0.00000                        0
-DEAL::0.577350 0.00000 0.577350                0.577350 0.00000 0.577350                       0
-DEAL::0.577350 0.577350 0.00000                0.577350 0.577350 0.00000                       0
-DEAL::0.577350 0.577350 0.577350               0.577350 0.577350 0.577350                      0
-DEAL::
-DEAL::-0.577350 -0.577350 0.577350             -0.577350 -0.577350 0.577350                    0
-DEAL::-0.577350 0.00000 0.577350               -0.577350 0.00000 0.577350                      0
-DEAL::0.00000 -0.577350 0.577350               0.00000 -0.577350 0.577350                      0
-DEAL::0.00000 0.00000 0.577350         0.00000 0.00000 0.577350                        0
-DEAL::
-DEAL::-0.577350 0.00000 0.577350               -0.577350 0.00000 0.577350                      0
-DEAL::-0.577350 0.577350 0.577350              -0.577350 0.577350 0.577350                     0
-DEAL::0.00000 0.00000 0.577350         0.00000 0.00000 0.577350                        0
-DEAL::0.00000 0.577350 0.577350                0.00000 0.577350 0.577350                       0
-DEAL::
-DEAL::0.00000 -0.577350 0.577350               0.00000 -0.577350 0.577350                      0
-DEAL::0.00000 0.00000 0.577350         0.00000 0.00000 0.577350                        0
-DEAL::0.577350 -0.577350 0.577350              0.577350 -0.577350 0.577350                     0
-DEAL::0.577350 0.00000 0.577350                0.577350 0.00000 0.577350                       0
-DEAL::
-DEAL::0.00000 0.00000 0.577350         0.00000 0.00000 0.577350                        0
-DEAL::0.00000 0.577350 0.577350                0.00000 0.577350 0.577350                       0
-DEAL::0.577350 0.00000 0.577350                0.577350 0.00000 0.577350                       0
-DEAL::0.577350 0.577350 0.577350               0.577350 0.577350 0.577350                      0
-DEAL::
-DEAL::-0.577350 -0.577350 -0.577350            -0.577350 -0.577350 -0.577350                   0
-DEAL::-0.577350 0.00000 -0.577350              -0.577350 0.00000 -0.577350                     0
-DEAL::-0.577350 -0.577350 0.00000              -0.577350 -0.577350 0.00000                     0
-DEAL::-0.577350 0.00000 0.00000                -0.577350 0.00000 0.00000                       0
-DEAL::
-DEAL::-0.577350 0.00000 -0.577350              -0.577350 0.00000 -0.577350                     0
-DEAL::-0.577350 0.577350 -0.577350             -0.577350 0.577350 -0.577350                    0
-DEAL::-0.577350 0.00000 0.00000                -0.577350 0.00000 0.00000                       0
-DEAL::-0.577350 0.577350 0.00000               -0.577350 0.577350 0.00000                      0
-DEAL::
-DEAL::-0.577350 -0.577350 0.00000              -0.577350 -0.577350 0.00000                     0
-DEAL::-0.577350 0.00000 0.00000                -0.577350 0.00000 0.00000                       0
-DEAL::-0.577350 -0.577350 0.577350             -0.577350 -0.577350 0.577350                    0
-DEAL::-0.577350 0.00000 0.577350               -0.577350 0.00000 0.577350                      0
-DEAL::
-DEAL::-0.577350 0.00000 0.00000                -0.577350 0.00000 0.00000                       0
-DEAL::-0.577350 0.577350 0.00000               -0.577350 0.577350 0.00000                      0
-DEAL::-0.577350 0.00000 0.577350               -0.577350 0.00000 0.577350                      0
-DEAL::-0.577350 0.577350 0.577350              -0.577350 0.577350 0.577350                     0
-DEAL::
-DEAL::-0.577350 -0.577350 -0.577350            -0.577350 -0.577350 -0.577350                   0
-DEAL::-0.577350 -0.577350 0.00000              -0.577350 -0.577350 0.00000                     0
-DEAL::0.00000 -0.577350 -0.577350              0.00000 -0.577350 -0.577350                     0
-DEAL::0.00000 -0.577350 0.00000                0.00000 -0.577350 0.00000                       0
-DEAL::
-DEAL::-0.577350 -0.577350 0.00000              -0.577350 -0.577350 0.00000                     0
-DEAL::-0.577350 -0.577350 0.577350             -0.577350 -0.577350 0.577350                    0
-DEAL::0.00000 -0.577350 0.00000                0.00000 -0.577350 0.00000                       0
-DEAL::0.00000 -0.577350 0.577350               0.00000 -0.577350 0.577350                      0
-DEAL::
-DEAL::0.00000 -0.577350 -0.577350              0.00000 -0.577350 -0.577350                     0
-DEAL::0.00000 -0.577350 0.00000                0.00000 -0.577350 0.00000                       0
-DEAL::0.577350 -0.577350 -0.577350             0.577350 -0.577350 -0.577350                    0
-DEAL::0.577350 -0.577350 0.00000               0.577350 -0.577350 0.00000                      0
-DEAL::
-DEAL::0.00000 -0.577350 0.00000                0.00000 -0.577350 0.00000                       0
-DEAL::0.00000 -0.577350 0.577350               0.00000 -0.577350 0.577350                      0
-DEAL::0.577350 -0.577350 0.00000               0.577350 -0.577350 0.00000                      0
-DEAL::0.577350 -0.577350 0.577350              0.577350 -0.577350 0.577350                     0
-DEAL::
-DEAL::-0.577350 0.577350 -0.577350             -0.577350 0.577350 -0.577350                    0
-DEAL::0.00000 0.577350 -0.577350               0.00000 0.577350 -0.577350                      0
-DEAL::-0.577350 0.577350 0.00000               -0.577350 0.577350 0.00000                      0
-DEAL::0.00000 0.577350 0.00000         0.00000 0.577350 0.00000                        0
-DEAL::
-DEAL::0.00000 0.577350 -0.577350               0.00000 0.577350 -0.577350                      0
-DEAL::0.577350 0.577350 -0.577350              0.577350 0.577350 -0.577350                     0
-DEAL::0.00000 0.577350 0.00000         0.00000 0.577350 0.00000                        0
-DEAL::0.577350 0.577350 0.00000                0.577350 0.577350 0.00000                       0
-DEAL::
-DEAL::-0.577350 0.577350 0.00000               -0.577350 0.577350 0.00000                      0
-DEAL::0.00000 0.577350 0.00000         0.00000 0.577350 0.00000                        0
-DEAL::-0.577350 0.577350 0.577350              -0.577350 0.577350 0.577350                     0
-DEAL::0.00000 0.577350 0.577350                0.00000 0.577350 0.577350                       0
-DEAL::
-DEAL::0.00000 0.577350 0.00000         0.00000 0.577350 0.00000                        0
-DEAL::0.577350 0.577350 0.00000                0.577350 0.577350 0.00000                       0
-DEAL::0.00000 0.577350 0.577350                0.00000 0.577350 0.577350                       0
-DEAL::0.577350 0.577350 0.577350               0.577350 0.577350 0.577350                      0
-DEAL::
--0.577350 -0.577350 -0.577350 1 0
-0.00000 -0.577350 -0.577350 1 0
-0.00000 0.00000 -0.577350 1 0
--0.577350 0.00000 -0.577350 1 0
--0.577350 -0.577350 -0.577350 1 0
-
-
-0.00000 -0.577350 -0.577350 1 0
-0.577350 -0.577350 -0.577350 1 0
-0.577350 0.00000 -0.577350 1 0
-0.00000 0.00000 -0.577350 1 0
-0.00000 -0.577350 -0.577350 1 0
-
-
--0.577350 0.00000 -0.577350 1 0
-0.00000 0.00000 -0.577350 1 0
-0.00000 0.577350 -0.577350 1 0
--0.577350 0.577350 -0.577350 1 0
--0.577350 0.00000 -0.577350 1 0
-
-
-0.00000 0.00000 -0.577350 1 0
-0.577350 0.00000 -0.577350 1 0
-0.577350 0.577350 -0.577350 1 0
-0.00000 0.577350 -0.577350 1 0
-0.00000 0.00000 -0.577350 1 0
-
-
-0.577350 -0.577350 -0.577350 1 0
-0.577350 -0.577350 0.00000 1 0
-0.577350 0.00000 0.00000 1 0
-0.577350 0.00000 -0.577350 1 0
-0.577350 -0.577350 -0.577350 1 0
-
-
-0.577350 -0.577350 0.00000 1 0
-0.577350 -0.577350 0.577350 1 0
-0.577350 0.00000 0.577350 1 0
-0.577350 0.00000 0.00000 1 0
-0.577350 -0.577350 0.00000 1 0
-
-
-0.577350 0.00000 -0.577350 1 0
-0.577350 0.00000 0.00000 1 0
-0.577350 0.577350 0.00000 1 0
-0.577350 0.577350 -0.577350 1 0
-0.577350 0.00000 -0.577350 1 0
-
-
-0.577350 0.00000 0.00000 1 0
-0.577350 0.00000 0.577350 1 0
-0.577350 0.577350 0.577350 1 0
-0.577350 0.577350 0.00000 1 0
-0.577350 0.00000 0.00000 1 0
-
-
--0.577350 -0.577350 0.577350 1 0
--0.577350 0.00000 0.577350 1 0
-0.00000 0.00000 0.577350 1 0
-0.00000 -0.577350 0.577350 1 0
--0.577350 -0.577350 0.577350 1 0
-
-
--0.577350 0.00000 0.577350 1 0
--0.577350 0.577350 0.577350 1 0
-0.00000 0.577350 0.577350 1 0
-0.00000 0.00000 0.577350 1 0
--0.577350 0.00000 0.577350 1 0
-
-
-0.00000 -0.577350 0.577350 1 0
-0.00000 0.00000 0.577350 1 0
-0.577350 0.00000 0.577350 1 0
-0.577350 -0.577350 0.577350 1 0
-0.00000 -0.577350 0.577350 1 0
-
-
-0.00000 0.00000 0.577350 1 0
-0.00000 0.577350 0.577350 1 0
-0.577350 0.577350 0.577350 1 0
-0.577350 0.00000 0.577350 1 0
-0.00000 0.00000 0.577350 1 0
-
-
--0.577350 -0.577350 -0.577350 1 0
--0.577350 0.00000 -0.577350 1 0
--0.577350 0.00000 0.00000 1 0
--0.577350 -0.577350 0.00000 1 0
--0.577350 -0.577350 -0.577350 1 0
-
-
--0.577350 0.00000 -0.577350 1 0
--0.577350 0.577350 -0.577350 1 0
--0.577350 0.577350 0.00000 1 0
--0.577350 0.00000 0.00000 1 0
--0.577350 0.00000 -0.577350 1 0
-
-
--0.577350 -0.577350 0.00000 1 0
--0.577350 0.00000 0.00000 1 0
--0.577350 0.00000 0.577350 1 0
--0.577350 -0.577350 0.577350 1 0
--0.577350 -0.577350 0.00000 1 0
-
-
--0.577350 0.00000 0.00000 1 0
--0.577350 0.577350 0.00000 1 0
--0.577350 0.577350 0.577350 1 0
--0.577350 0.00000 0.577350 1 0
--0.577350 0.00000 0.00000 1 0
-
-
--0.577350 -0.577350 -0.577350 1 0
--0.577350 -0.577350 0.00000 1 0
-0.00000 -0.577350 0.00000 1 0
-0.00000 -0.577350 -0.577350 1 0
--0.577350 -0.577350 -0.577350 1 0
-
-
--0.577350 -0.577350 0.00000 1 0
--0.577350 -0.577350 0.577350 1 0
-0.00000 -0.577350 0.577350 1 0
-0.00000 -0.577350 0.00000 1 0
--0.577350 -0.577350 0.00000 1 0
-
-
-0.00000 -0.577350 -0.577350 1 0
-0.00000 -0.577350 0.00000 1 0
-0.577350 -0.577350 0.00000 1 0
-0.577350 -0.577350 -0.577350 1 0
-0.00000 -0.577350 -0.577350 1 0
-
-
-0.00000 -0.577350 0.00000 1 0
-0.00000 -0.577350 0.577350 1 0
-0.577350 -0.577350 0.577350 1 0
-0.577350 -0.577350 0.00000 1 0
-0.00000 -0.577350 0.00000 1 0
-
-
--0.577350 0.577350 -0.577350 1 0
-0.00000 0.577350 -0.577350 1 0
-0.00000 0.577350 0.00000 1 0
--0.577350 0.577350 0.00000 1 0
--0.577350 0.577350 -0.577350 1 0
-
-
-0.00000 0.577350 -0.577350 1 0
-0.577350 0.577350 -0.577350 1 0
-0.577350 0.577350 0.00000 1 0
-0.00000 0.577350 0.00000 1 0
-0.00000 0.577350 -0.577350 1 0
-
-
--0.577350 0.577350 0.00000 1 0
-0.00000 0.577350 0.00000 1 0
-0.00000 0.577350 0.577350 1 0
--0.577350 0.577350 0.577350 1 0
--0.577350 0.577350 0.00000 1 0
-
-
-0.00000 0.577350 0.00000 1 0
-0.577350 0.577350 0.00000 1 0
-0.577350 0.577350 0.577350 1 0
-0.00000 0.577350 0.577350 1 0
-0.00000 0.577350 0.00000 1 0
-
-
+DEAL::Face centered at 0.00000 0.500000 0.500000
+DEAL::  Edge centered at 0.00000 0.00000 0.500000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::  Edge centered at 0.00000 1.00000 0.500000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::  Edge centered at 0.00000 0.500000 0.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::  Edge centered at 0.00000 0.500000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::Face centered at 1.00000 0.500000 0.500000
+DEAL::  Edge centered at 1.00000 0.00000 0.500000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::  Edge centered at 1.00000 1.00000 0.500000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::  Edge centered at 1.00000 0.500000 0.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::  Edge centered at 1.00000 0.500000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::Face centered at 0.500000 0.00000 0.500000
+DEAL::  Edge centered at 0.500000 0.00000 0.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::  Edge centered at 0.500000 0.00000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::  Edge centered at 0.00000 0.00000 0.500000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::  Edge centered at 1.00000 0.00000 0.500000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::Face centered at 0.500000 1.00000 0.500000
+DEAL::  Edge centered at 0.500000 1.00000 0.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::    normal_vector=0.00000 0.00000 -1.00000
+DEAL::  Edge centered at 0.500000 1.00000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::    normal_vector=0.00000 0.00000 1.00000
+DEAL::  Edge centered at 0.00000 1.00000 0.500000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::  Edge centered at 1.00000 1.00000 0.500000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::Face centered at 0.500000 0.500000 0.00000
+DEAL::  Edge centered at 0.00000 0.500000 0.00000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::  Edge centered at 1.00000 0.500000 0.00000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::  Edge centered at 0.500000 0.00000 0.00000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::  Edge centered at 0.500000 1.00000 0.00000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::Face centered at 0.500000 0.500000 1.00000
+DEAL::  Edge centered at 0.00000 0.500000 1.00000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::    normal_vector=-1.00000 0.00000 0.00000
+DEAL::  Edge centered at 1.00000 0.500000 1.00000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::    normal_vector=1.00000 0.00000 0.00000
+DEAL::  Edge centered at 0.500000 0.00000 1.00000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::    normal_vector=0.00000 -1.00000 0.00000
+DEAL::  Edge centered at 0.500000 1.00000 1.00000
+DEAL::    normal_vector=0.00000 1.00000 0.00000
+DEAL::    normal_vector=0.00000 1.00000 0.00000

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.