]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Added copy_material_to_manifold
authorLuca Heltai <luca.heltai@sissa.it>
Fri, 9 Jan 2015 19:33:39 +0000 (20:33 +0100)
committerLuca Heltai <luca.heltai@sissa.it>
Sat, 10 Jan 2015 16:15:06 +0000 (17:15 +0100)
doc/news/changes.h
include/deal.II/grid/grid_tools.h
source/grid/grid_tools.cc
source/grid/grid_tools.inst.in
tests/manifold/copy_material_to_manifold_id_01.cc [new file with mode: 0644]
tests/manifold/copy_material_to_manifold_id_01.output [new file with mode: 0644]

index 52be64b6d50c1825d989d63b15eaaee9df1762f1..25eb7de2ca310f2e73e11d332137952f9b257695 100644 (file)
@@ -134,8 +134,10 @@ inconvenience this causes.
 <h3>Specific improvements</h3>
 
 <ol>
-  <li> New: GridTools::copy_boundary_to_manifold_id() copies 
-  boundary_ids to manifold_ids for faces on the boundary.
+  <li> New: GridTools::copy_boundary_to_manifold_id() and 
+  GridTools::copy_material_to_manifold_id() copy
+  boundary_ids and material_ids to manifold_ids for 
+  faces on the boundary and for cells respectively.
   <br>
   (Luca Heltai, 2015/01/09)
   </li>
index 6208c6480ba70a4719fc83f611e5f5884fe7730e..1b99809908a74bc3c6228a3fb33ea864be1f9dd7 100644 (file)
@@ -1180,6 +1180,43 @@ namespace GridTools
   void copy_boundary_to_manifold_id(Triangulation<dim, spacedim> &tria,
                                     bool reset_boundary_ids=false);
 
+  /**
+   * Copy material ids to manifold ids. The default manifold_id for
+   * new Triangulation objects is numbers::invalid_manifold_id. When
+   * refinements occurs, the Triangulation asks where to locate new
+   * points to the underlying manifold.
+   *
+   * When reading a Triangulation from a supported input format,
+   * tipical informations that can be stored in a file are boundary
+   * conditions (which we store in boundary_ids), material types
+   * (which we store in material_ids) and in some cases subdomain ids
+   * (which we store in subdomain_ids).
+   *
+   * If you read one of these grids into a Triangulation, you might
+   * still want to use the material_id specified in the input file as
+   * a manifold_id description. In this case you can associate a
+   * Manifold object to internal cells, and this object will be used
+   * by the Triangulation to query Manifold objects for new
+   * points. This function iterates over active cells and copies the
+   * material_ids to the manifold_ids.
+   *
+   * The optional parameter @p compute_face_ids, indicates wether this
+   * function should also set the manifold_ids of the faces (both for
+   * internal faces and for faces on the boundary). If set to true,
+   * then each face will get a manifold_id equal to the minimum of the
+   * surrounding manifold_ids, ensuring that a unique manifold id is
+   * selected for each face of the Triangulation. By default, face
+   * manifold_ids are not computed.
+   *
+   * @ingroup manifold
+   *
+   * @author Luca Heltai, 2015
+   */
+  template <int dim, int spacedim>
+  void copy_material_to_manifold_id(Triangulation<dim, spacedim> &tria,
+                                    bool compute_face_ids=false);
+
+
   /*@}*/
   /**
    * @name Exceptions
index 38b91b6504868dcffd7b853c3689f99575137a8b..9f0a5b9b915ae30f52a6590e5acaa6424ee15693 100644 (file)
@@ -2913,6 +2913,30 @@ next_cell:
           }
   }
 
+  template <int dim, int spacedim>
+  void copy_material_to_manifold_id(Triangulation<dim, spacedim> &tria,
+                                    bool compute_face_ids)
+  {
+    typename Triangulation<dim,spacedim>::active_cell_iterator
+    cell=tria.begin_active(), endc=tria.end();
+
+    for (; cell != endc; ++cell)
+      {
+        cell->set_manifold_id(cell->material_id());
+        if (compute_face_ids == true)
+          {
+            for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
+              {
+                cell->face(f)->set_manifold_id(cell->material_id());
+                if (cell->neighbor(f) != endc)
+                  cell->face(f)->set_manifold_id
+                  (std::min(cell->material_id(),
+                            cell->neighbor(f)->material_id()));
+              }
+          }
+      }
+  }
+
 } /* namespace GridTools */
 
 
index cef43b58c7fdaf4c92f0ade4d2e537d6daed659e..246dcd76fd23060a4af49fd9539e248dde8cbeed 100644 (file)
@@ -394,6 +394,9 @@ for (deal_II_dimension : DIMENSIONS ; deal_II_space_dimension : SPACE_DIMENSIONS
  namespace GridTools \{
       template void copy_boundary_to_manifold_id<deal_II_dimension, deal_II_space_dimension>
       (Triangulation<deal_II_dimension, deal_II_space_dimension> &, bool);
+      template void copy_material_to_manifold_id<deal_II_dimension, deal_II_space_dimension>
+      (Triangulation<deal_II_dimension, deal_II_space_dimension> &, bool);
+
 \}
 #endif
 }
diff --git a/tests/manifold/copy_material_to_manifold_id_01.cc b/tests/manifold/copy_material_to_manifold_id_01.cc
new file mode 100644 (file)
index 0000000..ec654bc
--- /dev/null
@@ -0,0 +1,84 @@
+//------------------------------------------------------------
+//    Copyright (C) 2015, the deal.II authors
+//
+//    This file is subject to LGPL and may not be  distributed
+//    without copyright and license information. Please refer
+//    to the file deal.II/doc/license.html for the  text  and
+//    further information on this license.
+//
+//------------------------------------------------------------
+
+
+// Copy from boundary ids to manifold ids
+
+#include "../tests.h"
+#include <fstream>
+#include <base/logstream.h>
+
+
+// all include files you need here
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/tria_boundary_lib.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/grid_out.h>
+
+template<int dim, int spacedim>
+void print_info(Triangulation<dim,spacedim> &tria) {
+  typename Triangulation<dim,spacedim>::active_cell_iterator cell;
+
+  for (cell = tria.begin_active(); cell != tria.end(); ++cell)
+    {
+      deallog << "cell: " << cell 
+             << ", material_id: " 
+             << (int)cell->material_id()
+             << ", manifold_id: " 
+             << (int)cell->manifold_id() << std::endl;
+
+      for (unsigned int f=0; f<GeometryInfo<dim>::faces_per_cell; ++f)
+         deallog << "face: " << cell->face(f)
+                 << ", boundary_id: " 
+                 << (int)cell->face(f)->boundary_indicator()
+                 << ", manifold_id: " 
+                 << (int)cell->face(f)->manifold_id() << std::endl;
+    }
+}
+
+
+// Helper function
+template <int dim, int spacedim>
+void test()
+{
+  deallog << "Testing dim=" << dim
+          << ", spacedim="<< spacedim << std::endl;
+
+  Triangulation<dim,spacedim> tria;
+  GridGenerator::hyper_cube (tria, 0., 1.);
+  tria.refine_global(1);
+  tria.begin_active()->set_material_id(1);
+
+  deallog << "Original mesh ==============================" << std::endl;
+  print_info(tria);
+  GridTools::copy_material_to_manifold_id(tria);
+  deallog << "Copied mesh ================================" << std::endl;
+  print_info(tria);
+  GridTools::copy_material_to_manifold_id(tria, true);
+  deallog << "Copied mesh with boundary  =================" << std::endl;
+  print_info(tria);
+}
+
+int main ()
+{
+  initlog(true);
+  
+  test<1,1>();
+  test<1,2>();
+  test<2,2>();
+  test<2,3>();
+  test<3,3>();
+
+  return 0;
+}
+
diff --git a/tests/manifold/copy_material_to_manifold_id_01.output b/tests/manifold/copy_material_to_manifold_id_01.output
new file mode 100644 (file)
index 0000000..eb0b4c5
--- /dev/null
@@ -0,0 +1,345 @@
+
+DEAL::Testing dim=1, spacedim=1
+DEAL::Original mesh ==============================
+DEAL::cell: 1.0, material_id: 1, manifold_id: -1
+DEAL::face: 0, boundary_id: 0, manifold_id: -1
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: -1
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::face: 1, boundary_id: 1, manifold_id: -1
+DEAL::Copied mesh ================================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 0, boundary_id: 0, manifold_id: -1
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::face: 1, boundary_id: 1, manifold_id: -1
+DEAL::Copied mesh with boundary  =================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 0, boundary_id: 0, manifold_id: 1
+DEAL::face: 2, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 2, boundary_id: 255, manifold_id: 0
+DEAL::face: 1, boundary_id: 1, manifold_id: 0
+DEAL::Testing dim=1, spacedim=2
+DEAL::Original mesh ==============================
+DEAL::cell: 1.0, material_id: 1, manifold_id: -1
+DEAL::face: 0, boundary_id: 0, manifold_id: -1
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: -1
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::face: 1, boundary_id: 1, manifold_id: -1
+DEAL::Copied mesh ================================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 0, boundary_id: 0, manifold_id: -1
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 2, boundary_id: 255, manifold_id: -1
+DEAL::face: 1, boundary_id: 1, manifold_id: -1
+DEAL::Copied mesh with boundary  =================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 0, boundary_id: 0, manifold_id: 1
+DEAL::face: 2, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 2, boundary_id: 255, manifold_id: 0
+DEAL::face: 1, boundary_id: 1, manifold_id: 0
+DEAL::Testing dim=2, spacedim=2
+DEAL::Original mesh ==============================
+DEAL::cell: 1.0, material_id: 1, manifold_id: -1
+DEAL::face: 6, boundary_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 4, boundary_id: 0, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 8, boundary_id: 0, manifold_id: -1
+DEAL::face: 5, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.2, material_id: 0, manifold_id: -1
+DEAL::face: 7, boundary_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::face: 10, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.3, material_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 9, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::face: 11, boundary_id: 0, manifold_id: -1
+DEAL::Copied mesh ================================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 6, boundary_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 4, boundary_id: 0, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 8, boundary_id: 0, manifold_id: -1
+DEAL::face: 5, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.2, material_id: 0, manifold_id: 0
+DEAL::face: 7, boundary_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::face: 10, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.3, material_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 9, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::face: 11, boundary_id: 0, manifold_id: -1
+DEAL::Copied mesh with boundary  =================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 6, boundary_id: 0, manifold_id: 1
+DEAL::face: 12, boundary_id: 255, manifold_id: 0
+DEAL::face: 4, boundary_id: 0, manifold_id: 1
+DEAL::face: 14, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 12, boundary_id: 255, manifold_id: 0
+DEAL::face: 8, boundary_id: 0, manifold_id: 0
+DEAL::face: 5, boundary_id: 0, manifold_id: 0
+DEAL::face: 15, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.2, material_id: 0, manifold_id: 0
+DEAL::face: 7, boundary_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 255, manifold_id: 0
+DEAL::face: 14, boundary_id: 255, manifold_id: 0
+DEAL::face: 10, boundary_id: 0, manifold_id: 0
+DEAL::cell: 1.3, material_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 255, manifold_id: 0
+DEAL::face: 9, boundary_id: 0, manifold_id: 0
+DEAL::face: 15, boundary_id: 255, manifold_id: 0
+DEAL::face: 11, boundary_id: 0, manifold_id: 0
+DEAL::Testing dim=2, spacedim=3
+DEAL::Original mesh ==============================
+DEAL::cell: 1.0, material_id: 1, manifold_id: -1
+DEAL::face: 6, boundary_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 4, boundary_id: 0, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 8, boundary_id: 0, manifold_id: -1
+DEAL::face: 5, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.2, material_id: 0, manifold_id: -1
+DEAL::face: 7, boundary_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::face: 10, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.3, material_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 9, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::face: 11, boundary_id: 0, manifold_id: -1
+DEAL::Copied mesh ================================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 6, boundary_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 4, boundary_id: 0, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 12, boundary_id: 255, manifold_id: -1
+DEAL::face: 8, boundary_id: 0, manifold_id: -1
+DEAL::face: 5, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.2, material_id: 0, manifold_id: 0
+DEAL::face: 7, boundary_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 14, boundary_id: 255, manifold_id: -1
+DEAL::face: 10, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.3, material_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 255, manifold_id: -1
+DEAL::face: 9, boundary_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 255, manifold_id: -1
+DEAL::face: 11, boundary_id: 0, manifold_id: -1
+DEAL::Copied mesh with boundary  =================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 6, boundary_id: 0, manifold_id: 1
+DEAL::face: 12, boundary_id: 255, manifold_id: 0
+DEAL::face: 4, boundary_id: 0, manifold_id: 1
+DEAL::face: 14, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 12, boundary_id: 255, manifold_id: 0
+DEAL::face: 8, boundary_id: 0, manifold_id: 0
+DEAL::face: 5, boundary_id: 0, manifold_id: 0
+DEAL::face: 15, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.2, material_id: 0, manifold_id: 0
+DEAL::face: 7, boundary_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 255, manifold_id: 0
+DEAL::face: 14, boundary_id: 255, manifold_id: 0
+DEAL::face: 10, boundary_id: 0, manifold_id: 0
+DEAL::cell: 1.3, material_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 255, manifold_id: 0
+DEAL::face: 9, boundary_id: 0, manifold_id: 0
+DEAL::face: 15, boundary_id: 255, manifold_id: 0
+DEAL::face: 11, boundary_id: 0, manifold_id: 0
+DEAL::Testing dim=3, spacedim=3
+DEAL::Original mesh ==============================
+DEAL::cell: 1.0, material_id: 1, manifold_id: -1
+DEAL::face: 14, boundary_id: 0, manifold_id: -1
+DEAL::face: 41, boundary_id: 255, manifold_id: -1
+DEAL::face: 6, boundary_id: 0, manifold_id: -1
+DEAL::face: 37, boundary_id: 255, manifold_id: -1
+DEAL::face: 10, boundary_id: 0, manifold_id: -1
+DEAL::face: 33, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: -1
+DEAL::face: 41, boundary_id: 255, manifold_id: -1
+DEAL::face: 18, boundary_id: 0, manifold_id: -1
+DEAL::face: 8, boundary_id: 0, manifold_id: -1
+DEAL::face: 35, boundary_id: 255, manifold_id: -1
+DEAL::face: 11, boundary_id: 0, manifold_id: -1
+DEAL::face: 32, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.2, material_id: 0, manifold_id: -1
+DEAL::face: 15, boundary_id: 0, manifold_id: -1
+DEAL::face: 40, boundary_id: 255, manifold_id: -1
+DEAL::face: 37, boundary_id: 255, manifold_id: -1
+DEAL::face: 22, boundary_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 0, manifold_id: -1
+DEAL::face: 31, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.3, material_id: 0, manifold_id: -1
+DEAL::face: 40, boundary_id: 255, manifold_id: -1
+DEAL::face: 19, boundary_id: 0, manifold_id: -1
+DEAL::face: 35, boundary_id: 255, manifold_id: -1
+DEAL::face: 24, boundary_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 0, manifold_id: -1
+DEAL::face: 30, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.4, material_id: 0, manifold_id: -1
+DEAL::face: 16, boundary_id: 0, manifold_id: -1
+DEAL::face: 39, boundary_id: 255, manifold_id: -1
+DEAL::face: 7, boundary_id: 0, manifold_id: -1
+DEAL::face: 36, boundary_id: 255, manifold_id: -1
+DEAL::face: 33, boundary_id: 255, manifold_id: -1
+DEAL::face: 26, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.5, material_id: 0, manifold_id: -1
+DEAL::face: 39, boundary_id: 255, manifold_id: -1
+DEAL::face: 20, boundary_id: 0, manifold_id: -1
+DEAL::face: 9, boundary_id: 0, manifold_id: -1
+DEAL::face: 34, boundary_id: 255, manifold_id: -1
+DEAL::face: 32, boundary_id: 255, manifold_id: -1
+DEAL::face: 27, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.6, material_id: 0, manifold_id: -1
+DEAL::face: 17, boundary_id: 0, manifold_id: -1
+DEAL::face: 38, boundary_id: 255, manifold_id: -1
+DEAL::face: 36, boundary_id: 255, manifold_id: -1
+DEAL::face: 23, boundary_id: 0, manifold_id: -1
+DEAL::face: 31, boundary_id: 255, manifold_id: -1
+DEAL::face: 28, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.7, material_id: 0, manifold_id: -1
+DEAL::face: 38, boundary_id: 255, manifold_id: -1
+DEAL::face: 21, boundary_id: 0, manifold_id: -1
+DEAL::face: 34, boundary_id: 255, manifold_id: -1
+DEAL::face: 25, boundary_id: 0, manifold_id: -1
+DEAL::face: 30, boundary_id: 255, manifold_id: -1
+DEAL::face: 29, boundary_id: 0, manifold_id: -1
+DEAL::Copied mesh ================================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 14, boundary_id: 0, manifold_id: -1
+DEAL::face: 41, boundary_id: 255, manifold_id: -1
+DEAL::face: 6, boundary_id: 0, manifold_id: -1
+DEAL::face: 37, boundary_id: 255, manifold_id: -1
+DEAL::face: 10, boundary_id: 0, manifold_id: -1
+DEAL::face: 33, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 41, boundary_id: 255, manifold_id: -1
+DEAL::face: 18, boundary_id: 0, manifold_id: -1
+DEAL::face: 8, boundary_id: 0, manifold_id: -1
+DEAL::face: 35, boundary_id: 255, manifold_id: -1
+DEAL::face: 11, boundary_id: 0, manifold_id: -1
+DEAL::face: 32, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.2, material_id: 0, manifold_id: 0
+DEAL::face: 15, boundary_id: 0, manifold_id: -1
+DEAL::face: 40, boundary_id: 255, manifold_id: -1
+DEAL::face: 37, boundary_id: 255, manifold_id: -1
+DEAL::face: 22, boundary_id: 0, manifold_id: -1
+DEAL::face: 12, boundary_id: 0, manifold_id: -1
+DEAL::face: 31, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.3, material_id: 0, manifold_id: 0
+DEAL::face: 40, boundary_id: 255, manifold_id: -1
+DEAL::face: 19, boundary_id: 0, manifold_id: -1
+DEAL::face: 35, boundary_id: 255, manifold_id: -1
+DEAL::face: 24, boundary_id: 0, manifold_id: -1
+DEAL::face: 13, boundary_id: 0, manifold_id: -1
+DEAL::face: 30, boundary_id: 255, manifold_id: -1
+DEAL::cell: 1.4, material_id: 0, manifold_id: 0
+DEAL::face: 16, boundary_id: 0, manifold_id: -1
+DEAL::face: 39, boundary_id: 255, manifold_id: -1
+DEAL::face: 7, boundary_id: 0, manifold_id: -1
+DEAL::face: 36, boundary_id: 255, manifold_id: -1
+DEAL::face: 33, boundary_id: 255, manifold_id: -1
+DEAL::face: 26, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.5, material_id: 0, manifold_id: 0
+DEAL::face: 39, boundary_id: 255, manifold_id: -1
+DEAL::face: 20, boundary_id: 0, manifold_id: -1
+DEAL::face: 9, boundary_id: 0, manifold_id: -1
+DEAL::face: 34, boundary_id: 255, manifold_id: -1
+DEAL::face: 32, boundary_id: 255, manifold_id: -1
+DEAL::face: 27, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.6, material_id: 0, manifold_id: 0
+DEAL::face: 17, boundary_id: 0, manifold_id: -1
+DEAL::face: 38, boundary_id: 255, manifold_id: -1
+DEAL::face: 36, boundary_id: 255, manifold_id: -1
+DEAL::face: 23, boundary_id: 0, manifold_id: -1
+DEAL::face: 31, boundary_id: 255, manifold_id: -1
+DEAL::face: 28, boundary_id: 0, manifold_id: -1
+DEAL::cell: 1.7, material_id: 0, manifold_id: 0
+DEAL::face: 38, boundary_id: 255, manifold_id: -1
+DEAL::face: 21, boundary_id: 0, manifold_id: -1
+DEAL::face: 34, boundary_id: 255, manifold_id: -1
+DEAL::face: 25, boundary_id: 0, manifold_id: -1
+DEAL::face: 30, boundary_id: 255, manifold_id: -1
+DEAL::face: 29, boundary_id: 0, manifold_id: -1
+DEAL::Copied mesh with boundary  =================
+DEAL::cell: 1.0, material_id: 1, manifold_id: 1
+DEAL::face: 14, boundary_id: 0, manifold_id: 1
+DEAL::face: 41, boundary_id: 255, manifold_id: 0
+DEAL::face: 6, boundary_id: 0, manifold_id: 1
+DEAL::face: 37, boundary_id: 255, manifold_id: 0
+DEAL::face: 10, boundary_id: 0, manifold_id: 1
+DEAL::face: 33, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.1, material_id: 0, manifold_id: 0
+DEAL::face: 41, boundary_id: 255, manifold_id: 0
+DEAL::face: 18, boundary_id: 0, manifold_id: 0
+DEAL::face: 8, boundary_id: 0, manifold_id: 0
+DEAL::face: 35, boundary_id: 255, manifold_id: 0
+DEAL::face: 11, boundary_id: 0, manifold_id: 0
+DEAL::face: 32, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.2, material_id: 0, manifold_id: 0
+DEAL::face: 15, boundary_id: 0, manifold_id: 0
+DEAL::face: 40, boundary_id: 255, manifold_id: 0
+DEAL::face: 37, boundary_id: 255, manifold_id: 0
+DEAL::face: 22, boundary_id: 0, manifold_id: 0
+DEAL::face: 12, boundary_id: 0, manifold_id: 0
+DEAL::face: 31, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.3, material_id: 0, manifold_id: 0
+DEAL::face: 40, boundary_id: 255, manifold_id: 0
+DEAL::face: 19, boundary_id: 0, manifold_id: 0
+DEAL::face: 35, boundary_id: 255, manifold_id: 0
+DEAL::face: 24, boundary_id: 0, manifold_id: 0
+DEAL::face: 13, boundary_id: 0, manifold_id: 0
+DEAL::face: 30, boundary_id: 255, manifold_id: 0
+DEAL::cell: 1.4, material_id: 0, manifold_id: 0
+DEAL::face: 16, boundary_id: 0, manifold_id: 0
+DEAL::face: 39, boundary_id: 255, manifold_id: 0
+DEAL::face: 7, boundary_id: 0, manifold_id: 0
+DEAL::face: 36, boundary_id: 255, manifold_id: 0
+DEAL::face: 33, boundary_id: 255, manifold_id: 0
+DEAL::face: 26, boundary_id: 0, manifold_id: 0
+DEAL::cell: 1.5, material_id: 0, manifold_id: 0
+DEAL::face: 39, boundary_id: 255, manifold_id: 0
+DEAL::face: 20, boundary_id: 0, manifold_id: 0
+DEAL::face: 9, boundary_id: 0, manifold_id: 0
+DEAL::face: 34, boundary_id: 255, manifold_id: 0
+DEAL::face: 32, boundary_id: 255, manifold_id: 0
+DEAL::face: 27, boundary_id: 0, manifold_id: 0
+DEAL::cell: 1.6, material_id: 0, manifold_id: 0
+DEAL::face: 17, boundary_id: 0, manifold_id: 0
+DEAL::face: 38, boundary_id: 255, manifold_id: 0
+DEAL::face: 36, boundary_id: 255, manifold_id: 0
+DEAL::face: 23, boundary_id: 0, manifold_id: 0
+DEAL::face: 31, boundary_id: 255, manifold_id: 0
+DEAL::face: 28, boundary_id: 0, manifold_id: 0
+DEAL::cell: 1.7, material_id: 0, manifold_id: 0
+DEAL::face: 38, boundary_id: 255, manifold_id: 0
+DEAL::face: 21, boundary_id: 0, manifold_id: 0
+DEAL::face: 34, boundary_id: 255, manifold_id: 0
+DEAL::face: 25, boundary_id: 0, manifold_id: 0
+DEAL::face: 30, boundary_id: 255, manifold_id: 0
+DEAL::face: 29, boundary_id: 0, manifold_id: 0

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.