]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add an option to copy manifold ids in extrude_triangulation.
authorDavid Wells <drwells@email.unc.edu>
Sat, 4 Aug 2018 02:48:21 +0000 (22:48 -0400)
committerDavid Wells <drwells@email.unc.edu>
Mon, 13 Aug 2018 14:44:41 +0000 (10:44 -0400)
doc/news/changes/major/20180811DavidWells [new file with mode: 0644]
include/deal.II/grid/grid_generator.h
source/grid/grid_generator.cc
tests/grid/extrude_copy_manifold.cc [new file with mode: 0644]
tests/grid/extrude_copy_manifold.output [new file with mode: 0644]

diff --git a/doc/news/changes/major/20180811DavidWells b/doc/news/changes/major/20180811DavidWells
new file mode 100644 (file)
index 0000000..7a7ca6e
--- /dev/null
@@ -0,0 +1,4 @@
+Improved: GridGenerator::extrude_triangulation can now optionally extrude manifold
+ids.
+<br>
+(2018/08/11, David Wells)
index 46796992d7c4a03b1734387c206c6ad1d75b1613..bd277e3592cff21f4b349be3aadac69ca98bae95 100644 (file)
@@ -1280,18 +1280,56 @@ namespace GridGenerator
    * $0$, $1$, and $42$, then the $z = 0$ boundary id of @p result will be $43$
    * and the $z = \text{height}$ boundary id will be $44$.
    *
+   * This function does not, by default, copy manifold ids. The reason for
+   * this is that there is no way to set the manifold ids on the lines of the
+   * resulting Triangulation without more information: for example, if two
+   * faces of @p input with different manifold ids meet at a shared vertex then
+   * there is no <em>a priori</em> reason to pick one manifold id or another
+   * for the lines created in @p result that are parallel to the $z$-axis and
+   * pass through that point. If @p copy_manifold_ids is <code>true</code>
+   * then this function sets line manifold ids by picking the one that appears
+   * <em>first</em> in @p manifold_priorities. For example: if @p
+   * manifold_priorities is <code>{0, 42, numbers::flat_manifold_id}</code>
+   * and the line under consideration is adjacent to faces with manifold ids of
+   * <code>0</code> and <code>42</code>, then that line will have a manifold id
+   * of <code>0</code>. The correct ordering is almost always
+   * <ol>
+   *   <li>manifold ids set on the boundary,</li>
+   *   <li>manifold ids that describe most of the cells in the Triangulation
+   *   (e.g., numbers::flat_manifold_id), and</li>
+   *   <li>any manifold ids corresponding to TransfiniteInterpolationManifold
+   *   manifolds.</li>
+   * </ol>
+   *
+   * In particular, since TransfiniteInterpolationManifold interpolates
+   * between surrounding manifolds, its manifold id should usually not be set
+   * on lines or faces that are adjacent to cells with different manifold
+   * ids. The default value for @p manifold_priorities follows this ranking
+   * (where each category is sorted in ascending order):
+   * <ol>
+   *   <li>manifold ids associated with manifolds that are not
+   *   TransfiniteInterpolationManifold, and</li>
+   *   <li>manifold ids associated with any TransfiniteInterpolationManifold
+   *   objects.</li>
+   * </ol>
+   * Note that numbers::flat_manifold_id (should it be a manifold id of @p
+   * input) will always be the last entry in the first category.
+   *
    * @note The 2d input triangulation @p input must be a coarse mesh that has
    * no refined cells.
    *
    * @note Since @p input and @p output have different spatial dimensions no
-   * manifold objects are copied (nor are any manifold ids set) by this
-   * function.
+   * manifold objects are copied by this function regardless of the value of
+   * @p copy_manifold_ids.
    */
   void
-  extrude_triangulation(const Triangulation<2, 2> &input,
-                        const unsigned int         n_slices,
-                        const double               height,
-                        Triangulation<3, 3> &      result);
+  extrude_triangulation(
+    const Triangulation<2, 2> &            input,
+    const unsigned int                     n_slices,
+    const double                           height,
+    Triangulation<3, 3> &                  result,
+    const bool                             copy_manifold_ids   = false,
+    const std::vector<types::manifold_id> &manifold_priorities = {});
 
   /**
    * Overload of the previous function. Take a 2d Triangulation that is being
@@ -1312,9 +1350,12 @@ namespace GridGenerator
    * @author Weixiong Zheng, 2018
    */
   void
-  extrude_triangulation(const Triangulation<2, 2> &input,
-                        const std::vector<double> &slice_coordinates,
-                        Triangulation<3, 3> &      result);
+  extrude_triangulation(
+    const Triangulation<2, 2> &            input,
+    const std::vector<double> &            slice_coordinates,
+    Triangulation<3, 3> &                  result,
+    const bool                             copy_manifold_ids   = false,
+    const std::vector<types::manifold_id> &manifold_priorities = {});
 
 
   /**
index 672638ff3866468513a8e466c35f50887e3bc56e..fa87331b0297db5a440b8781f8e452cd039b4f62 100644 (file)
@@ -4476,10 +4476,13 @@ namespace GridGenerator
 
 
   void
-  extrude_triangulation(const Triangulation<2, 2> &input,
-                        const unsigned int         n_slices,
-                        const double               height,
-                        Triangulation<3, 3> &      result)
+  extrude_triangulation(
+    const Triangulation<2, 2> &            input,
+    const unsigned int                     n_slices,
+    const double                           height,
+    Triangulation<3, 3> &                  result,
+    const bool                             copy_manifold_ids,
+    const std::vector<types::manifold_id> &manifold_priorities)
   {
     Assert(input.n_levels() == 1,
            ExcMessage(
@@ -4497,13 +4500,19 @@ namespace GridGenerator
     std::vector<double> slices_z_values;
     for (unsigned int i = 0; i < n_slices; ++i)
       slices_z_values.push_back(i * delta_h);
-    extrude_triangulation(input, slices_z_values, result);
+    extrude_triangulation(
+      input, slices_z_values, result, copy_manifold_ids, manifold_priorities);
   }
 
+
+
   void
-  extrude_triangulation(const Triangulation<2, 2> &input,
-                        const std::vector<double> &slice_coordinates,
-                        Triangulation<3, 3> &      result)
+  extrude_triangulation(
+    const Triangulation<2, 2> &            input,
+    const std::vector<double> &            slice_coordinates,
+    Triangulation<3, 3> &                  result,
+    const bool                             copy_manifold_ids,
+    const std::vector<types::manifold_id> &manifold_priorities)
   {
     Assert(input.n_levels() == 1,
            ExcMessage(
@@ -4516,6 +4525,66 @@ namespace GridGenerator
              "The number of slices for extrusion must be at least 2."));
     Assert(std::is_sorted(slice_coordinates.begin(), slice_coordinates.end()),
            ExcMessage("Slice z-coordinates should be in ascending order"));
+
+    const auto priorities = [&]() -> std::vector<types::manifold_id> {
+      // if a non-empty (i.e., not the default) vector is given for
+      // manifold_priorities then use it (but check its validity in debug
+      // mode)
+      if (0 < manifold_priorities.size())
+        {
+#ifdef DEBUG
+          // check that the provided manifold_priorities is valid
+          std::vector<types::manifold_id> sorted_manifold_priorities =
+            manifold_priorities;
+          std::sort(sorted_manifold_priorities.begin(),
+                    sorted_manifold_priorities.end());
+          Assert(std::unique(sorted_manifold_priorities.begin(),
+                             sorted_manifold_priorities.end()) ==
+                   sorted_manifold_priorities.end(),
+                 ExcMessage(
+                   "The given vector of manifold ids may not contain any "
+                   "duplicated entries."));
+          std::vector<types::manifold_id> sorted_manifold_ids =
+            input.get_manifold_ids();
+          std::sort(sorted_manifold_ids.begin(), sorted_manifold_ids.end());
+          if (sorted_manifold_priorities != sorted_manifold_ids)
+            {
+              std::ostringstream message;
+              message << "The given triangulation has manifold ids {";
+              for (const types::manifold_id manifold_id : sorted_manifold_ids)
+                if (manifold_id != sorted_manifold_ids.back())
+                  message << manifold_id << ", ";
+              message << sorted_manifold_ids.back() << "}, but \n"
+                      << "    the given vector of manifold ids is {";
+              for (const types::manifold_id manifold_id : manifold_priorities)
+                if (manifold_id != manifold_priorities.back())
+                  message << manifold_id << ", ";
+              message
+                << manifold_priorities.back() << "}.\n"
+                << "    These vectors should contain the same elements.\n";
+              const std::string m = message.str();
+              Assert(false, ExcMessage(m));
+            }
+#endif
+          return manifold_priorities;
+        }
+      // otherwise use the default ranking: ascending order, but TFI manifolds
+      // are at the end.
+      std::vector<types::manifold_id> default_priorities =
+        input.get_manifold_ids();
+      const auto first_tfi_it = std::partition(
+        default_priorities.begin(),
+        default_priorities.end(),
+        [&input](const types::manifold_id &id) {
+          return dynamic_cast<const TransfiniteInterpolationManifold<2, 2> *>(
+                   &input.get_manifold(id)) == nullptr;
+        });
+      std::sort(default_priorities.begin(), first_tfi_it);
+      std::sort(first_tfi_it, default_priorities.end());
+
+      return default_priorities;
+    }();
+
     const std::size_t        n_slices = slice_coordinates.size();
     std::vector<Point<3>>    points(n_slices * input.n_vertices());
     std::vector<CellData<3>> cells;
@@ -4523,9 +4592,9 @@ namespace GridGenerator
 
     // copy the array of points as many times as there will be slices,
     // one slice at a time. The z-axis value are defined in slices_coordinates
-    for (unsigned int slice_n = 0; slice_n < n_slices; ++slice_n)
+    for (std::size_t slice_n = 0; slice_n < n_slices; ++slice_n)
       {
-        for (unsigned int vertex_n = 0; vertex_n < input.n_vertices();
+        for (std::size_t vertex_n = 0; vertex_n < input.n_vertices();
              ++vertex_n)
           {
             const Point<2> vertex = input.get_vertices()[vertex_n];
@@ -4538,7 +4607,7 @@ namespace GridGenerator
     // time
     for (const auto &cell : input.active_cell_iterators())
       {
-        for (unsigned int slice_n = 0; slice_n < n_slices - 1; ++slice_n)
+        for (std::size_t slice_n = 0; slice_n < n_slices - 1; ++slice_n)
           {
             CellData<3> this_cell;
             for (unsigned int vertex_n = 0;
@@ -4554,14 +4623,14 @@ namespace GridGenerator
               }
 
             this_cell.material_id = cell->material_id();
+            if (copy_manifold_ids)
+              this_cell.manifold_id = cell->manifold_id();
             cells.push_back(this_cell);
           }
       }
 
-    // next, create face data for all of the outer faces for which the
-    // boundary indicator will not be equal to zero (where we would
-    // explicitly set it to something that is already the default --
-    // no need to do that)
+    // Next, create face data for all faces that are orthogonal to the x-y
+    // plane
     SubCellData               subcell_data;
     std::vector<CellData<2>> &quads           = subcell_data.boundary_quads;
     types::boundary_id        max_boundary_id = 0;
@@ -4573,20 +4642,44 @@ namespace GridGenerator
         quad.boundary_id = face->boundary_id();
         if (face->at_boundary())
           max_boundary_id = std::max(max_boundary_id, quad.boundary_id);
-        for (unsigned int slice = 0; slice < n_slices - 1; ++slice)
+        if (copy_manifold_ids)
+          quad.manifold_id = face->manifold_id();
+        for (std::size_t slice_n = 0; slice_n < n_slices - 1; ++slice_n)
           {
             quad.vertices[0] =
-              face->vertex_index(0) + slice * input.n_vertices();
+              face->vertex_index(0) + slice_n * input.n_vertices();
             quad.vertices[1] =
-              face->vertex_index(1) + slice * input.n_vertices();
+              face->vertex_index(1) + slice_n * input.n_vertices();
             quad.vertices[2] =
-              face->vertex_index(0) + (slice + 1) * input.n_vertices();
+              face->vertex_index(0) + (slice_n + 1) * input.n_vertices();
             quad.vertices[3] =
-              face->vertex_index(1) + (slice + 1) * input.n_vertices();
+              face->vertex_index(1) + (slice_n + 1) * input.n_vertices();
             quads.push_back(quad);
           }
       }
 
+    // if necessary, create face data for faces parallel to the x-y
+    // plane. This is only necessary if we need to set manifolds.
+    if (copy_manifold_ids)
+      for (const auto &cell : input.active_cell_iterators())
+        {
+          CellData<2> quad;
+          quad.boundary_id = numbers::internal_face_boundary_id;
+          quad.manifold_id = cell->manifold_id(); // check is outside loop
+          for (std::size_t slice_n = 1; slice_n < n_slices - 1; ++slice_n)
+            {
+              quad.vertices[0] =
+                cell->vertex_index(0) + slice_n * input.n_vertices();
+              quad.vertices[1] =
+                cell->vertex_index(1) + slice_n * input.n_vertices();
+              quad.vertices[2] =
+                cell->vertex_index(2) + slice_n * input.n_vertices();
+              quad.vertices[3] =
+                cell->vertex_index(3) + slice_n * input.n_vertices();
+              quads.push_back(quad);
+            }
+        }
+
     // then mark the bottom and top boundaries of the extruded mesh
     // with max_boundary_id+1 and max_boundary_id+2. check that this
     // remains valid
@@ -4609,6 +4702,8 @@ namespace GridGenerator
         quad.vertices[1] = cell->vertex_index(1);
         quad.vertices[2] = cell->vertex_index(2);
         quad.vertices[3] = cell->vertex_index(3);
+        if (copy_manifold_ids)
+          quad.manifold_id = cell->manifold_id();
         quads.push_back(quad);
 
         quad.boundary_id = top_boundary_id;
@@ -4616,6 +4711,8 @@ namespace GridGenerator
              vertex_n < GeometryInfo<3>::vertices_per_face;
              ++vertex_n)
           quad.vertices[vertex_n] += (n_slices - 1) * input.n_vertices();
+        if (copy_manifold_ids)
+          quad.manifold_id = cell->manifold_id();
         quads.push_back(quad);
       }
 
@@ -4628,6 +4725,16 @@ namespace GridGenerator
     // as discussed in the edge orientation paper mentioned in the
     // introduction to the GridReordering class.
     result.create_triangulation(points, cells, subcell_data);
+
+    for (auto manifold_id_it = priorities.rbegin();
+         manifold_id_it != priorities.rend();
+         ++manifold_id_it)
+      for (auto &face : result.active_face_iterators())
+        if (face->manifold_id() == *manifold_id_it)
+          for (unsigned int line_n = 0;
+               line_n < GeometryInfo<3>::lines_per_face;
+               ++line_n)
+            face->line(line_n)->set_manifold_id(*manifold_id_it);
   }
 
 
diff --git a/tests/grid/extrude_copy_manifold.cc b/tests/grid/extrude_copy_manifold.cc
new file mode 100644 (file)
index 0000000..979b2ee
--- /dev/null
@@ -0,0 +1,141 @@
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2018 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE.md at
+// the top level directory of deal.II.
+//
+// ---------------------------------------------------------------------
+
+// Check that we can correctly copy manifold ids during extrusion. The output
+// was manually checked by plotting points in paraview.
+
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/manifold_lib.h>
+#include <deal.II/grid/tria.h>
+
+#include "../tests.h"
+
+void
+test()
+{
+  Triangulation<2> triangulation_2;
+  GridGenerator::hyper_ball(triangulation_2);
+  for (auto &cell : triangulation_2.active_cell_iterators())
+    if (cell->at_boundary())
+      {
+        cell->set_manifold_id(1);
+        for (unsigned int face_n = 0; face_n < GeometryInfo<2>::faces_per_cell;
+             ++face_n)
+          if (!cell->face(face_n)->at_boundary())
+            cell->face(face_n)->set_manifold_id(1);
+      }
+  for (auto &cell : triangulation_2.active_cell_iterators())
+    if (!cell->at_boundary())
+      cell->set_all_manifold_ids(numbers::flat_manifold_id);
+  TransfiniteInterpolationManifold<2> tfi_manifold_2;
+  tfi_manifold_2.initialize(triangulation_2);
+  triangulation_2.set_manifold(1, tfi_manifold_2);
+
+  Triangulation<3> triangulation_3;
+  GridGenerator::extrude_triangulation(
+    triangulation_2, 3, 1.0, triangulation_3, true);
+  TransfiniteInterpolationManifold<3> tfi_manifold;
+  tfi_manifold.initialize(triangulation_3);
+  CylindricalManifold<3> cylinder_manifold(2);
+  triangulation_3.set_manifold(0, cylinder_manifold);
+  triangulation_3.set_manifold(42, tfi_manifold);
+  triangulation_3.refine_global(1);
+
+  std::map<types::manifold_id, std::vector<Point<3>>> line_centers;
+  std::map<types::manifold_id, std::vector<Point<3>>> face_centers;
+  std::map<types::manifold_id, std::vector<Point<3>>> cell_centers;
+
+  for (const auto &cell : triangulation_3.active_cell_iterators())
+    {
+      cell_centers[cell->manifold_id()].push_back(cell->center());
+      for (unsigned int face_n = 0; face_n < GeometryInfo<3>::faces_per_cell;
+           ++face_n)
+        face_centers[cell->face(face_n)->manifold_id()].push_back(
+          cell->face(face_n)->center());
+      for (unsigned int line_n = 0; line_n < GeometryInfo<3>::lines_per_cell;
+           ++line_n)
+        line_centers[cell->line(line_n)->manifold_id()].push_back(
+          cell->line(line_n)->center());
+    }
+
+  deallog << "cell centers:" << std::endl;
+  for (const std::pair<const types::manifold_id, std::vector<Point<3>>> &key :
+       cell_centers)
+    {
+      // cell centers are already unique: we only visited each cell once
+      deallog << "manifold id: " << key.first << std::endl;
+      for (const Point<3> &point : key.second)
+        deallog << point[0] << ", " << point[1] << ", " << point[2]
+                << std::endl;
+    }
+  deallog << std::endl;
+
+  // The generated arrays of face and line centers contain duplicates: get rid
+  // of them by sorting and then std::unique-ing
+  auto point_comparator = [](const Point<3> &a, const Point<3> &b) {
+    // just use std::tuple's lexical ordering so that we don't have to think
+    // too hard
+    return std::tie(a[0], a[1], a[2]) < std::tie(b[0], b[1], b[2]);
+  };
+
+  deallog << "face centers:" << std::endl;
+  for (std::pair<const types::manifold_id, std::vector<Point<3>>> &key :
+       face_centers)
+    {
+      deallog << "manifold id: " << key.first << std::endl;
+      std::sort(key.second.begin(), key.second.end(), point_comparator);
+      key.second.erase(std::unique(key.second.begin(), key.second.end()),
+                       key.second.end());
+      for (const Point<3> &point : key.second)
+        deallog << point[0] << ", " << point[1] << ", " << point[2]
+                << std::endl;
+    }
+
+  deallog << "line centers:" << std::endl;
+  for (std::pair<const types::manifold_id, std::vector<Point<3>>> &key :
+       line_centers)
+    {
+      deallog << "manifold id: " << key.first << std::endl;
+      std::sort(key.second.begin(), key.second.end(), point_comparator);
+      key.second.erase(std::unique(key.second.begin(), key.second.end()),
+                       key.second.end());
+      for (const Point<3> &point : key.second)
+        deallog << point[0] << ", " << point[1] << ", " << point[2]
+                << std::endl;
+    }
+
+  // reenable if we want to look at pictures
+  if (false)
+    {
+      std::ofstream out("out-" + std::to_string(3) + ".vtu");
+      GridOut       grid_out;
+      grid_out.write_vtu(triangulation_3, out);
+
+      Triangulation<2, 3> triangulation_23;
+      GridGenerator::extract_boundary_mesh(triangulation_3, triangulation_23);
+      std::ofstream out23("out-23.vtu");
+      grid_out.write_vtu(triangulation_23, out23);
+    }
+}
+
+
+int
+main()
+{
+  initlog();
+  test();
+}
diff --git a/tests/grid/extrude_copy_manifold.output b/tests/grid/extrude_copy_manifold.output
new file mode 100644 (file)
index 0000000..d7bf109
--- /dev/null
@@ -0,0 +1,689 @@
+
+DEAL::cell centers:
+DEAL::manifold id: 1
+DEAL::-0.301777, -0.713388, 0.125000
+DEAL::0.301777, -0.713388, 0.125000
+DEAL::-0.198223, -0.433058, 0.125000
+DEAL::0.198223, -0.433058, 0.125000
+DEAL::-0.301777, -0.713388, 0.375000
+DEAL::0.301777, -0.713388, 0.375000
+DEAL::-0.198223, -0.433058, 0.375000
+DEAL::0.198223, -0.433058, 0.375000
+DEAL::-0.301777, -0.713388, 0.625000
+DEAL::0.301777, -0.713388, 0.625000
+DEAL::-0.198223, -0.433058, 0.625000
+DEAL::0.198223, -0.433058, 0.625000
+DEAL::-0.301777, -0.713388, 0.875000
+DEAL::0.301777, -0.713388, 0.875000
+DEAL::-0.198223, -0.433058, 0.875000
+DEAL::0.198223, -0.433058, 0.875000
+DEAL::-0.713388, -0.301777, 0.125000
+DEAL::-0.433058, -0.198223, 0.125000
+DEAL::-0.713388, 0.301777, 0.125000
+DEAL::-0.433058, 0.198223, 0.125000
+DEAL::-0.713388, -0.301777, 0.375000
+DEAL::-0.433058, -0.198223, 0.375000
+DEAL::-0.713388, 0.301777, 0.375000
+DEAL::-0.433058, 0.198223, 0.375000
+DEAL::-0.713388, -0.301777, 0.625000
+DEAL::-0.433058, -0.198223, 0.625000
+DEAL::-0.713388, 0.301777, 0.625000
+DEAL::-0.433058, 0.198223, 0.625000
+DEAL::-0.713388, -0.301777, 0.875000
+DEAL::-0.433058, -0.198223, 0.875000
+DEAL::-0.713388, 0.301777, 0.875000
+DEAL::-0.433058, 0.198223, 0.875000
+DEAL::0.713388, -0.301777, 0.125000
+DEAL::0.713388, 0.301777, 0.125000
+DEAL::0.433058, -0.198223, 0.125000
+DEAL::0.433058, 0.198223, 0.125000
+DEAL::0.713388, -0.301777, 0.375000
+DEAL::0.713388, 0.301777, 0.375000
+DEAL::0.433058, -0.198223, 0.375000
+DEAL::0.433058, 0.198223, 0.375000
+DEAL::0.713388, -0.301777, 0.625000
+DEAL::0.713388, 0.301777, 0.625000
+DEAL::0.433058, -0.198223, 0.625000
+DEAL::0.433058, 0.198223, 0.625000
+DEAL::0.713388, -0.301777, 0.875000
+DEAL::0.713388, 0.301777, 0.875000
+DEAL::0.433058, -0.198223, 0.875000
+DEAL::0.433058, 0.198223, 0.875000
+DEAL::-0.301777, 0.713388, 0.125000
+DEAL::-0.198223, 0.433058, 0.125000
+DEAL::0.301777, 0.713388, 0.125000
+DEAL::0.198223, 0.433058, 0.125000
+DEAL::-0.301777, 0.713388, 0.375000
+DEAL::-0.198223, 0.433058, 0.375000
+DEAL::0.301777, 0.713388, 0.375000
+DEAL::0.198223, 0.433058, 0.375000
+DEAL::-0.301777, 0.713388, 0.625000
+DEAL::-0.198223, 0.433058, 0.625000
+DEAL::0.301777, 0.713388, 0.625000
+DEAL::0.198223, 0.433058, 0.625000
+DEAL::-0.301777, 0.713388, 0.875000
+DEAL::-0.198223, 0.433058, 0.875000
+DEAL::0.301777, 0.713388, 0.875000
+DEAL::0.198223, 0.433058, 0.875000
+DEAL::manifold id: 4294967295
+DEAL::-0.146447, -0.146447, 0.125000
+DEAL::0.146447, -0.146447, 0.125000
+DEAL::-0.146447, 0.146447, 0.125000
+DEAL::0.146447, 0.146447, 0.125000
+DEAL::-0.146447, -0.146447, 0.375000
+DEAL::0.146447, -0.146447, 0.375000
+DEAL::-0.146447, 0.146447, 0.375000
+DEAL::0.146447, 0.146447, 0.375000
+DEAL::-0.146447, -0.146447, 0.625000
+DEAL::0.146447, -0.146447, 0.625000
+DEAL::-0.146447, 0.146447, 0.625000
+DEAL::0.146447, 0.146447, 0.625000
+DEAL::-0.146447, -0.146447, 0.875000
+DEAL::0.146447, -0.146447, 0.875000
+DEAL::-0.146447, 0.146447, 0.875000
+DEAL::0.146447, 0.146447, 0.875000
+DEAL::
+DEAL::face centers:
+DEAL::manifold id: 0
+DEAL::-0.853553, -0.353553, 0.125000
+DEAL::-0.853553, -0.353553, 0.375000
+DEAL::-0.853553, -0.353553, 0.625000
+DEAL::-0.853553, -0.353553, 0.875000
+DEAL::-0.853553, 0.353553, 0.125000
+DEAL::-0.853553, 0.353553, 0.375000
+DEAL::-0.853553, 0.353553, 0.625000
+DEAL::-0.853553, 0.353553, 0.875000
+DEAL::-0.353553, -0.853553, 0.125000
+DEAL::-0.353553, -0.853553, 0.625000
+DEAL::-0.353553, -0.853553, 0.375000
+DEAL::-0.353553, -0.853553, 0.875000
+DEAL::-0.353553, 0.853553, 0.375000
+DEAL::-0.353553, 0.853553, 0.875000
+DEAL::-0.353553, 0.853553, 0.125000
+DEAL::-0.353553, 0.853553, 0.625000
+DEAL::0.353553, -0.853553, 0.125000
+DEAL::0.353553, -0.853553, 0.625000
+DEAL::0.353553, -0.853553, 0.375000
+DEAL::0.353553, -0.853553, 0.875000
+DEAL::0.353553, 0.853553, 0.125000
+DEAL::0.353553, 0.853553, 0.375000
+DEAL::0.353553, 0.853553, 0.625000
+DEAL::0.353553, 0.853553, 0.875000
+DEAL::0.853553, -0.353553, 0.375000
+DEAL::0.853553, -0.353553, 0.875000
+DEAL::0.853553, -0.353553, 0.125000
+DEAL::0.853553, -0.353553, 0.625000
+DEAL::0.853553, 0.353553, 0.125000
+DEAL::0.853553, 0.353553, 0.375000
+DEAL::0.853553, 0.353553, 0.625000
+DEAL::0.853553, 0.353553, 0.875000
+DEAL::manifold id: 1
+DEAL::-0.823223, 1.66533e-16, 0.125000
+DEAL::-0.823223, 1.66533e-16, 0.625000
+DEAL::-0.823223, 1.66533e-16, 0.375000
+DEAL::-0.823223, 1.66533e-16, 0.875000
+DEAL::-0.713388, 0.301777, 0.00000
+DEAL::-0.713388, 0.301777, 0.500000
+DEAL::-0.713388, 0.301777, 1.00000
+DEAL::-0.713388, -0.301777, 0.00000
+DEAL::-0.713388, -0.301777, 0.250000
+DEAL::-0.713388, -0.301777, 0.500000
+DEAL::-0.713388, -0.301777, 0.750000
+DEAL::-0.713388, -0.301777, 1.00000
+DEAL::-0.713388, 0.301777, 0.250000
+DEAL::-0.713388, 0.301777, 0.750000
+DEAL::-0.603553, -0.603553, 0.125000
+DEAL::-0.603553, -0.603553, 0.375000
+DEAL::-0.603553, -0.603553, 0.625000
+DEAL::-0.603553, -0.603553, 0.875000
+DEAL::-0.603553, 0.603553, 0.125000
+DEAL::-0.603553, 0.603553, 0.375000
+DEAL::-0.603553, 0.603553, 0.625000
+DEAL::-0.603553, 0.603553, 0.875000
+DEAL::-0.573223, -0.250000, 0.125000
+DEAL::-0.573223, -0.250000, 0.375000
+DEAL::-0.573223, -0.250000, 0.625000
+DEAL::-0.573223, -0.250000, 0.875000
+DEAL::-0.573223, 0.250000, 0.125000
+DEAL::-0.573223, 0.250000, 0.375000
+DEAL::-0.573223, 0.250000, 0.625000
+DEAL::-0.573223, 0.250000, 0.875000
+DEAL::-0.469670, 5.55112e-17, 0.125000
+DEAL::-0.469670, 5.55112e-17, 0.375000
+DEAL::-0.469670, 5.55112e-17, 0.625000
+DEAL::-0.469670, 5.55112e-17, 0.875000
+DEAL::-0.433058, -0.198223, 0.00000
+DEAL::-0.433058, -0.198223, 0.500000
+DEAL::-0.433058, -0.198223, 1.00000
+DEAL::-0.433058, 0.198223, 0.00000
+DEAL::-0.433058, 0.198223, 0.250000
+DEAL::-0.433058, 0.198223, 0.500000
+DEAL::-0.433058, 0.198223, 0.750000
+DEAL::-0.433058, 0.198223, 1.00000
+DEAL::-0.433058, -0.198223, 0.250000
+DEAL::-0.433058, -0.198223, 0.750000
+DEAL::-0.396447, -0.396447, 0.125000
+DEAL::-0.396447, -0.396447, 0.625000
+DEAL::-0.396447, -0.396447, 0.375000
+DEAL::-0.396447, -0.396447, 0.875000
+DEAL::-0.396447, 0.396447, 0.125000
+DEAL::-0.396447, 0.396447, 0.375000
+DEAL::-0.396447, 0.396447, 0.625000
+DEAL::-0.396447, 0.396447, 0.875000
+DEAL::-0.301777, -0.713388, 0.00000
+DEAL::-0.301777, -0.713388, 0.250000
+DEAL::-0.301777, -0.713388, 0.500000
+DEAL::-0.301777, -0.713388, 0.750000
+DEAL::-0.301777, -0.713388, 1.00000
+DEAL::-0.301777, 0.713388, 0.00000
+DEAL::-0.301777, 0.713388, 0.500000
+DEAL::-0.301777, 0.713388, 1.00000
+DEAL::-0.301777, 0.713388, 0.250000
+DEAL::-0.301777, 0.713388, 0.750000
+DEAL::-0.250000, -0.573223, 0.125000
+DEAL::-0.250000, -0.573223, 0.375000
+DEAL::-0.250000, -0.573223, 0.625000
+DEAL::-0.250000, -0.573223, 0.875000
+DEAL::-0.250000, 0.573223, 0.125000
+DEAL::-0.250000, 0.573223, 0.375000
+DEAL::-0.250000, 0.573223, 0.625000
+DEAL::-0.250000, 0.573223, 0.875000
+DEAL::-0.198223, -0.433058, 0.00000
+DEAL::-0.198223, -0.433058, 0.500000
+DEAL::-0.198223, -0.433058, 1.00000
+DEAL::-0.198223, 0.433058, 0.00000
+DEAL::-0.198223, 0.433058, 0.500000
+DEAL::-0.198223, 0.433058, 1.00000
+DEAL::-0.198223, -0.433058, 0.250000
+DEAL::-0.198223, -0.433058, 0.750000
+DEAL::-0.198223, 0.433058, 0.250000
+DEAL::-0.198223, 0.433058, 0.750000
+DEAL::-4.16334e-17, -0.823223, 0.125000
+DEAL::-4.16334e-17, -0.823223, 0.625000
+DEAL::-4.16334e-17, -0.823223, 0.375000
+DEAL::-4.16334e-17, -0.823223, 0.875000
+DEAL::-1.38778e-17, -0.469670, 0.125000
+DEAL::-1.38778e-17, -0.469670, 0.375000
+DEAL::-1.38778e-17, -0.469670, 0.625000
+DEAL::-1.38778e-17, -0.469670, 0.875000
+DEAL::1.11022e-16, 0.469670, 0.125000
+DEAL::1.11022e-16, 0.469670, 0.375000
+DEAL::1.11022e-16, 0.469670, 0.625000
+DEAL::1.11022e-16, 0.469670, 0.875000
+DEAL::3.05311e-16, 0.823223, 0.375000
+DEAL::3.05311e-16, 0.823223, 0.875000
+DEAL::3.05311e-16, 0.823223, 0.125000
+DEAL::3.05311e-16, 0.823223, 0.625000
+DEAL::0.198223, -0.433058, 0.250000
+DEAL::0.198223, -0.433058, 0.750000
+DEAL::0.198223, -0.433058, 0.00000
+DEAL::0.198223, -0.433058, 0.500000
+DEAL::0.198223, -0.433058, 1.00000
+DEAL::0.198223, 0.433058, 0.00000
+DEAL::0.198223, 0.433058, 0.500000
+DEAL::0.198223, 0.433058, 1.00000
+DEAL::0.198223, 0.433058, 0.250000
+DEAL::0.198223, 0.433058, 0.750000
+DEAL::0.250000, -0.573223, 0.125000
+DEAL::0.250000, -0.573223, 0.375000
+DEAL::0.250000, -0.573223, 0.625000
+DEAL::0.250000, -0.573223, 0.875000
+DEAL::0.250000, 0.573223, 0.125000
+DEAL::0.250000, 0.573223, 0.375000
+DEAL::0.250000, 0.573223, 0.625000
+DEAL::0.250000, 0.573223, 0.875000
+DEAL::0.301777, -0.713388, 0.250000
+DEAL::0.301777, -0.713388, 0.750000
+DEAL::0.301777, -0.713388, 0.00000
+DEAL::0.301777, -0.713388, 0.500000
+DEAL::0.301777, -0.713388, 1.00000
+DEAL::0.301777, 0.713388, 0.00000
+DEAL::0.301777, 0.713388, 0.500000
+DEAL::0.301777, 0.713388, 1.00000
+DEAL::0.301777, 0.713388, 0.250000
+DEAL::0.301777, 0.713388, 0.750000
+DEAL::0.396447, -0.396447, 0.375000
+DEAL::0.396447, -0.396447, 0.875000
+DEAL::0.396447, -0.396447, 0.125000
+DEAL::0.396447, -0.396447, 0.625000
+DEAL::0.396447, 0.396447, 0.375000
+DEAL::0.396447, 0.396447, 0.875000
+DEAL::0.396447, 0.396447, 0.125000
+DEAL::0.396447, 0.396447, 0.625000
+DEAL::0.433058, -0.198223, 0.00000
+DEAL::0.433058, -0.198223, 0.500000
+DEAL::0.433058, -0.198223, 1.00000
+DEAL::0.433058, -0.198223, 0.250000
+DEAL::0.433058, -0.198223, 0.750000
+DEAL::0.433058, 0.198223, 0.250000
+DEAL::0.433058, 0.198223, 0.750000
+DEAL::0.433058, 0.198223, 0.00000
+DEAL::0.433058, 0.198223, 0.500000
+DEAL::0.433058, 0.198223, 1.00000
+DEAL::0.469670, -6.93889e-18, 0.125000
+DEAL::0.469670, -6.93889e-18, 0.375000
+DEAL::0.469670, -6.93889e-18, 0.625000
+DEAL::0.469670, -6.93889e-18, 0.875000
+DEAL::0.573223, -0.250000, 0.125000
+DEAL::0.573223, -0.250000, 0.375000
+DEAL::0.573223, -0.250000, 0.625000
+DEAL::0.573223, -0.250000, 0.875000
+DEAL::0.573223, 0.250000, 0.125000
+DEAL::0.573223, 0.250000, 0.375000
+DEAL::0.573223, 0.250000, 0.625000
+DEAL::0.573223, 0.250000, 0.875000
+DEAL::0.603553, -0.603553, 0.125000
+DEAL::0.603553, -0.603553, 0.375000
+DEAL::0.603553, -0.603553, 0.625000
+DEAL::0.603553, -0.603553, 0.875000
+DEAL::0.603553, 0.603553, 0.125000
+DEAL::0.603553, 0.603553, 0.375000
+DEAL::0.603553, 0.603553, 0.625000
+DEAL::0.603553, 0.603553, 0.875000
+DEAL::0.713388, -0.301777, 0.00000
+DEAL::0.713388, -0.301777, 0.250000
+DEAL::0.713388, -0.301777, 0.500000
+DEAL::0.713388, -0.301777, 0.750000
+DEAL::0.713388, -0.301777, 1.00000
+DEAL::0.713388, 0.301777, 0.250000
+DEAL::0.713388, 0.301777, 0.750000
+DEAL::0.713388, 0.301777, 0.00000
+DEAL::0.713388, 0.301777, 0.500000
+DEAL::0.713388, 0.301777, 1.00000
+DEAL::0.823223, -3.46945e-17, 0.375000
+DEAL::0.823223, -3.46945e-17, 0.875000
+DEAL::0.823223, -3.46945e-17, 0.125000
+DEAL::0.823223, -3.46945e-17, 0.625000
+DEAL::manifold id: 4294967295
+DEAL::-0.292893, -0.146447, 0.125000
+DEAL::-0.292893, -0.146447, 0.375000
+DEAL::-0.292893, -0.146447, 0.625000
+DEAL::-0.292893, -0.146447, 0.875000
+DEAL::-0.292893, 0.146447, 0.125000
+DEAL::-0.292893, 0.146447, 0.375000
+DEAL::-0.292893, 0.146447, 0.625000
+DEAL::-0.292893, 0.146447, 0.875000
+DEAL::-0.146447, -0.292893, 0.125000
+DEAL::-0.146447, -0.292893, 0.375000
+DEAL::-0.146447, -0.292893, 0.625000
+DEAL::-0.146447, -0.292893, 0.875000
+DEAL::-0.146447, -0.146447, 0.00000
+DEAL::-0.146447, -0.146447, 0.250000
+DEAL::-0.146447, -0.146447, 0.500000
+DEAL::-0.146447, -0.146447, 0.750000
+DEAL::-0.146447, -0.146447, 1.00000
+DEAL::-0.146447, 0.00000, 0.125000
+DEAL::-0.146447, 0.00000, 0.375000
+DEAL::-0.146447, 0.00000, 0.625000
+DEAL::-0.146447, 0.00000, 0.875000
+DEAL::-0.146447, 0.146447, 0.00000
+DEAL::-0.146447, 0.146447, 0.250000
+DEAL::-0.146447, 0.146447, 0.500000
+DEAL::-0.146447, 0.146447, 0.750000
+DEAL::-0.146447, 0.146447, 1.00000
+DEAL::-0.146447, 0.292893, 0.125000
+DEAL::-0.146447, 0.292893, 0.375000
+DEAL::-0.146447, 0.292893, 0.625000
+DEAL::-0.146447, 0.292893, 0.875000
+DEAL::0.00000, -0.146447, 0.125000
+DEAL::0.00000, -0.146447, 0.375000
+DEAL::0.00000, -0.146447, 0.625000
+DEAL::0.00000, -0.146447, 0.875000
+DEAL::0.00000, 0.146447, 0.125000
+DEAL::0.00000, 0.146447, 0.375000
+DEAL::0.00000, 0.146447, 0.625000
+DEAL::0.00000, 0.146447, 0.875000
+DEAL::0.146447, -0.292893, 0.125000
+DEAL::0.146447, -0.292893, 0.375000
+DEAL::0.146447, -0.292893, 0.625000
+DEAL::0.146447, -0.292893, 0.875000
+DEAL::0.146447, -0.146447, 0.00000
+DEAL::0.146447, -0.146447, 0.250000
+DEAL::0.146447, -0.146447, 0.500000
+DEAL::0.146447, -0.146447, 0.750000
+DEAL::0.146447, -0.146447, 1.00000
+DEAL::0.146447, 0.00000, 0.125000
+DEAL::0.146447, 0.00000, 0.375000
+DEAL::0.146447, 0.00000, 0.625000
+DEAL::0.146447, 0.00000, 0.875000
+DEAL::0.146447, 0.146447, 0.00000
+DEAL::0.146447, 0.146447, 0.250000
+DEAL::0.146447, 0.146447, 0.500000
+DEAL::0.146447, 0.146447, 0.750000
+DEAL::0.146447, 0.146447, 1.00000
+DEAL::0.146447, 0.292893, 0.125000
+DEAL::0.146447, 0.292893, 0.375000
+DEAL::0.146447, 0.292893, 0.625000
+DEAL::0.146447, 0.292893, 0.875000
+DEAL::0.292893, -0.146447, 0.125000
+DEAL::0.292893, -0.146447, 0.375000
+DEAL::0.292893, -0.146447, 0.625000
+DEAL::0.292893, -0.146447, 0.875000
+DEAL::0.292893, 0.146447, 0.125000
+DEAL::0.292893, 0.146447, 0.375000
+DEAL::0.292893, 0.146447, 0.625000
+DEAL::0.292893, 0.146447, 0.875000
+DEAL::line centers:
+DEAL::manifold id: 0
+DEAL::-1.00000, 2.22045e-16, 0.125000
+DEAL::-1.00000, 2.22045e-16, 0.375000
+DEAL::-1.00000, 2.22045e-16, 0.625000
+DEAL::-1.00000, 2.22045e-16, 0.875000
+DEAL::-0.853553, -0.353553, 0.00000
+DEAL::-0.853553, -0.353553, 0.500000
+DEAL::-0.853553, -0.353553, 1.00000
+DEAL::-0.853553, 0.353553, 0.00000
+DEAL::-0.853553, 0.353553, 0.500000
+DEAL::-0.853553, 0.353553, 1.00000
+DEAL::-0.853553, -0.353553, 0.250000
+DEAL::-0.853553, -0.353553, 0.750000
+DEAL::-0.853553, 0.353553, 0.250000
+DEAL::-0.853553, 0.353553, 0.750000
+DEAL::-0.707107, -0.707107, 0.125000
+DEAL::-0.707107, -0.707107, 0.375000
+DEAL::-0.707107, -0.707107, 0.625000
+DEAL::-0.707107, -0.707107, 0.875000
+DEAL::-0.707107, 0.707107, 0.125000
+DEAL::-0.707107, 0.707107, 0.375000
+DEAL::-0.707107, 0.707107, 0.625000
+DEAL::-0.707107, 0.707107, 0.875000
+DEAL::-0.353553, -0.853553, 0.00000
+DEAL::-0.353553, -0.853553, 0.500000
+DEAL::-0.353553, -0.853553, 1.00000
+DEAL::-0.353553, -0.853553, 0.250000
+DEAL::-0.353553, -0.853553, 0.750000
+DEAL::-0.353553, 0.853553, 0.00000
+DEAL::-0.353553, 0.853553, 0.500000
+DEAL::-0.353553, 0.853553, 1.00000
+DEAL::-0.353553, 0.853553, 0.250000
+DEAL::-0.353553, 0.853553, 0.750000
+DEAL::-5.55112e-17, -1.00000, 0.125000
+DEAL::-5.55112e-17, -1.00000, 0.375000
+DEAL::-5.55112e-17, -1.00000, 0.625000
+DEAL::-5.55112e-17, -1.00000, 0.875000
+DEAL::3.88578e-16, 1.00000, 0.125000
+DEAL::3.88578e-16, 1.00000, 0.375000
+DEAL::3.88578e-16, 1.00000, 0.625000
+DEAL::3.88578e-16, 1.00000, 0.875000
+DEAL::0.353553, -0.853553, 0.250000
+DEAL::0.353553, -0.853553, 0.750000
+DEAL::0.353553, -0.853553, 0.00000
+DEAL::0.353553, -0.853553, 0.500000
+DEAL::0.353553, -0.853553, 1.00000
+DEAL::0.353553, 0.853553, 0.00000
+DEAL::0.353553, 0.853553, 0.500000
+DEAL::0.353553, 0.853553, 1.00000
+DEAL::0.353553, 0.853553, 0.250000
+DEAL::0.353553, 0.853553, 0.750000
+DEAL::0.707107, -0.707107, 0.125000
+DEAL::0.707107, -0.707107, 0.375000
+DEAL::0.707107, -0.707107, 0.625000
+DEAL::0.707107, -0.707107, 0.875000
+DEAL::0.707107, 0.707107, 0.125000
+DEAL::0.707107, 0.707107, 0.375000
+DEAL::0.707107, 0.707107, 0.625000
+DEAL::0.707107, 0.707107, 0.875000
+DEAL::0.853553, -0.353553, 0.250000
+DEAL::0.853553, -0.353553, 0.750000
+DEAL::0.853553, 0.353553, 0.250000
+DEAL::0.853553, 0.353553, 0.750000
+DEAL::0.853553, -0.353553, 0.00000
+DEAL::0.853553, -0.353553, 0.500000
+DEAL::0.853553, -0.353553, 1.00000
+DEAL::0.853553, 0.353553, 0.00000
+DEAL::0.853553, 0.353553, 0.500000
+DEAL::0.853553, 0.353553, 1.00000
+DEAL::1.00000, -5.55112e-17, 0.125000
+DEAL::1.00000, -5.55112e-17, 0.375000
+DEAL::1.00000, -5.55112e-17, 0.625000
+DEAL::1.00000, -5.55112e-17, 0.875000
+DEAL::manifold id: 1
+DEAL::-0.823223, 1.66533e-16, 0.00000
+DEAL::-0.823223, 1.66533e-16, 0.500000
+DEAL::-0.823223, 1.66533e-16, 1.00000
+DEAL::-0.823223, 1.66533e-16, 0.250000
+DEAL::-0.823223, 1.66533e-16, 0.750000
+DEAL::-0.646447, 1.11022e-16, 0.125000
+DEAL::-0.646447, 1.11022e-16, 0.375000
+DEAL::-0.646447, 1.11022e-16, 0.625000
+DEAL::-0.646447, 1.11022e-16, 0.875000
+DEAL::-0.603553, -0.603553, 0.00000
+DEAL::-0.603553, -0.603553, 0.500000
+DEAL::-0.603553, -0.603553, 1.00000
+DEAL::-0.603553, 0.603553, 0.00000
+DEAL::-0.603553, 0.603553, 0.500000
+DEAL::-0.603553, 0.603553, 1.00000
+DEAL::-0.603553, -0.603553, 0.250000
+DEAL::-0.603553, -0.603553, 0.750000
+DEAL::-0.603553, 0.603553, 0.250000
+DEAL::-0.603553, 0.603553, 0.750000
+DEAL::-0.573223, -0.250000, 0.00000
+DEAL::-0.573223, -0.250000, 0.500000
+DEAL::-0.573223, -0.250000, 1.00000
+DEAL::-0.573223, -0.250000, 0.250000
+DEAL::-0.573223, -0.250000, 0.750000
+DEAL::-0.573223, 0.250000, 0.00000
+DEAL::-0.573223, 0.250000, 0.250000
+DEAL::-0.573223, 0.250000, 0.500000
+DEAL::-0.573223, 0.250000, 0.750000
+DEAL::-0.573223, 0.250000, 1.00000
+DEAL::-0.500000, -0.500000, 0.125000
+DEAL::-0.500000, -0.500000, 0.375000
+DEAL::-0.500000, -0.500000, 0.625000
+DEAL::-0.500000, -0.500000, 0.875000
+DEAL::-0.500000, 0.500000, 0.125000
+DEAL::-0.500000, 0.500000, 0.375000
+DEAL::-0.500000, 0.500000, 0.625000
+DEAL::-0.500000, 0.500000, 0.875000
+DEAL::-0.469670, 5.55112e-17, 0.00000
+DEAL::-0.469670, 5.55112e-17, 0.250000
+DEAL::-0.469670, 5.55112e-17, 0.500000
+DEAL::-0.469670, 5.55112e-17, 0.750000
+DEAL::-0.469670, 5.55112e-17, 1.00000
+DEAL::-0.396447, -0.396447, 0.00000
+DEAL::-0.396447, -0.396447, 0.500000
+DEAL::-0.396447, -0.396447, 1.00000
+DEAL::-0.396447, 0.396447, 0.00000
+DEAL::-0.396447, 0.396447, 0.500000
+DEAL::-0.396447, 0.396447, 1.00000
+DEAL::-0.396447, -0.396447, 0.250000
+DEAL::-0.396447, -0.396447, 0.750000
+DEAL::-0.396447, 0.396447, 0.250000
+DEAL::-0.396447, 0.396447, 0.750000
+DEAL::-0.250000, -0.573223, 0.00000
+DEAL::-0.250000, -0.573223, 0.500000
+DEAL::-0.250000, -0.573223, 1.00000
+DEAL::-0.250000, -0.573223, 0.250000
+DEAL::-0.250000, -0.573223, 0.750000
+DEAL::-0.250000, 0.573223, 0.00000
+DEAL::-0.250000, 0.573223, 0.500000
+DEAL::-0.250000, 0.573223, 1.00000
+DEAL::-0.250000, 0.573223, 0.250000
+DEAL::-0.250000, 0.573223, 0.750000
+DEAL::-4.16334e-17, -0.823223, 0.00000
+DEAL::-4.16334e-17, -0.823223, 0.500000
+DEAL::-4.16334e-17, -0.823223, 1.00000
+DEAL::-4.16334e-17, -0.823223, 0.250000
+DEAL::-4.16334e-17, -0.823223, 0.750000
+DEAL::-2.77556e-17, -0.646447, 0.125000
+DEAL::-2.77556e-17, -0.646447, 0.375000
+DEAL::-2.77556e-17, -0.646447, 0.625000
+DEAL::-2.77556e-17, -0.646447, 0.875000
+DEAL::-1.38778e-17, -0.469670, 0.00000
+DEAL::-1.38778e-17, -0.469670, 0.250000
+DEAL::-1.38778e-17, -0.469670, 0.500000
+DEAL::-1.38778e-17, -0.469670, 0.750000
+DEAL::-1.38778e-17, -0.469670, 1.00000
+DEAL::5.55112e-17, 0.469670, 0.00000
+DEAL::5.55112e-17, 0.469670, 0.500000
+DEAL::5.55112e-17, 0.469670, 1.00000
+DEAL::1.38778e-16, 0.823223, 0.00000
+DEAL::1.38778e-16, 0.823223, 0.500000
+DEAL::1.38778e-16, 0.823223, 1.00000
+DEAL::1.66533e-16, 0.469670, 0.250000
+DEAL::1.66533e-16, 0.469670, 0.750000
+DEAL::2.22045e-16, 0.646447, 0.125000
+DEAL::2.22045e-16, 0.646447, 0.375000
+DEAL::2.22045e-16, 0.646447, 0.625000
+DEAL::2.22045e-16, 0.646447, 0.875000
+DEAL::4.71845e-16, 0.823223, 0.250000
+DEAL::4.71845e-16, 0.823223, 0.750000
+DEAL::0.250000, -0.573223, 0.250000
+DEAL::0.250000, -0.573223, 0.750000
+DEAL::0.250000, -0.573223, 0.00000
+DEAL::0.250000, -0.573223, 0.500000
+DEAL::0.250000, -0.573223, 1.00000
+DEAL::0.250000, 0.573223, 0.00000
+DEAL::0.250000, 0.573223, 0.500000
+DEAL::0.250000, 0.573223, 1.00000
+DEAL::0.250000, 0.573223, 0.250000
+DEAL::0.250000, 0.573223, 0.750000
+DEAL::0.396447, -0.396447, 0.250000
+DEAL::0.396447, -0.396447, 0.750000
+DEAL::0.396447, -0.396447, 0.00000
+DEAL::0.396447, -0.396447, 0.500000
+DEAL::0.396447, -0.396447, 1.00000
+DEAL::0.396447, 0.396447, 0.250000
+DEAL::0.396447, 0.396447, 0.750000
+DEAL::0.396447, 0.396447, 0.00000
+DEAL::0.396447, 0.396447, 0.500000
+DEAL::0.396447, 0.396447, 1.00000
+DEAL::0.469670, -1.38778e-17, 0.00000
+DEAL::0.469670, -1.38778e-17, 0.500000
+DEAL::0.469670, -1.38778e-17, 1.00000
+DEAL::0.469670, 0.00000, 0.250000
+DEAL::0.469670, 0.00000, 0.750000
+DEAL::0.500000, -0.500000, 0.125000
+DEAL::0.500000, -0.500000, 0.375000
+DEAL::0.500000, -0.500000, 0.625000
+DEAL::0.500000, -0.500000, 0.875000
+DEAL::0.500000, 0.500000, 0.125000
+DEAL::0.500000, 0.500000, 0.375000
+DEAL::0.500000, 0.500000, 0.625000
+DEAL::0.500000, 0.500000, 0.875000
+DEAL::0.573223, -0.250000, 0.00000
+DEAL::0.573223, -0.250000, 0.500000
+DEAL::0.573223, -0.250000, 1.00000
+DEAL::0.573223, -0.250000, 0.250000
+DEAL::0.573223, -0.250000, 0.750000
+DEAL::0.573223, 0.250000, 0.250000
+DEAL::0.573223, 0.250000, 0.750000
+DEAL::0.573223, 0.250000, 0.00000
+DEAL::0.573223, 0.250000, 0.500000
+DEAL::0.573223, 0.250000, 1.00000
+DEAL::0.603553, -0.603553, 0.250000
+DEAL::0.603553, -0.603553, 0.750000
+DEAL::0.603553, -0.603553, 0.00000
+DEAL::0.603553, -0.603553, 0.500000
+DEAL::0.603553, -0.603553, 1.00000
+DEAL::0.603553, 0.603553, 0.250000
+DEAL::0.603553, 0.603553, 0.750000
+DEAL::0.603553, 0.603553, 0.00000
+DEAL::0.603553, 0.603553, 0.500000
+DEAL::0.603553, 0.603553, 1.00000
+DEAL::0.646447, -1.38778e-17, 0.125000
+DEAL::0.646447, -1.38778e-17, 0.375000
+DEAL::0.646447, -1.38778e-17, 0.625000
+DEAL::0.646447, -1.38778e-17, 0.875000
+DEAL::0.823223, -2.77556e-17, 0.250000
+DEAL::0.823223, -2.77556e-17, 0.750000
+DEAL::0.823223, -4.16334e-17, 0.00000
+DEAL::0.823223, -4.16334e-17, 0.500000
+DEAL::0.823223, -4.16334e-17, 1.00000
+DEAL::manifold id: 4294967295
+DEAL::-0.292893, -0.292893, 0.125000
+DEAL::-0.292893, -0.292893, 0.375000
+DEAL::-0.292893, -0.292893, 0.625000
+DEAL::-0.292893, -0.292893, 0.875000
+DEAL::-0.292893, -0.146447, 0.00000
+DEAL::-0.292893, -0.146447, 0.250000
+DEAL::-0.292893, -0.146447, 0.500000
+DEAL::-0.292893, -0.146447, 0.750000
+DEAL::-0.292893, -0.146447, 1.00000
+DEAL::-0.292893, 0.00000, 0.125000
+DEAL::-0.292893, 0.00000, 0.375000
+DEAL::-0.292893, 0.00000, 0.625000
+DEAL::-0.292893, 0.00000, 0.875000
+DEAL::-0.292893, 0.146447, 0.00000
+DEAL::-0.292893, 0.146447, 0.250000
+DEAL::-0.292893, 0.146447, 0.500000
+DEAL::-0.292893, 0.146447, 0.750000
+DEAL::-0.292893, 0.146447, 1.00000
+DEAL::-0.292893, 0.292893, 0.125000
+DEAL::-0.292893, 0.292893, 0.375000
+DEAL::-0.292893, 0.292893, 0.625000
+DEAL::-0.292893, 0.292893, 0.875000
+DEAL::-0.146447, -0.292893, 0.00000
+DEAL::-0.146447, -0.292893, 0.250000
+DEAL::-0.146447, -0.292893, 0.500000
+DEAL::-0.146447, -0.292893, 0.750000
+DEAL::-0.146447, -0.292893, 1.00000
+DEAL::-0.146447, 0.00000, 0.00000
+DEAL::-0.146447, 0.00000, 0.250000
+DEAL::-0.146447, 0.00000, 0.500000
+DEAL::-0.146447, 0.00000, 0.750000
+DEAL::-0.146447, 0.00000, 1.00000
+DEAL::-0.146447, 0.292893, 0.00000
+DEAL::-0.146447, 0.292893, 0.250000
+DEAL::-0.146447, 0.292893, 0.500000
+DEAL::-0.146447, 0.292893, 0.750000
+DEAL::-0.146447, 0.292893, 1.00000
+DEAL::0.00000, -0.292893, 0.125000
+DEAL::0.00000, -0.292893, 0.375000
+DEAL::0.00000, -0.292893, 0.625000
+DEAL::0.00000, -0.292893, 0.875000
+DEAL::0.00000, -0.146447, 0.00000
+DEAL::0.00000, -0.146447, 0.250000
+DEAL::0.00000, -0.146447, 0.500000
+DEAL::0.00000, -0.146447, 0.750000
+DEAL::0.00000, -0.146447, 1.00000
+DEAL::0.00000, 0.00000, 0.125000
+DEAL::0.00000, 0.00000, 0.375000
+DEAL::0.00000, 0.00000, 0.625000
+DEAL::0.00000, 0.00000, 0.875000
+DEAL::0.00000, 0.146447, 0.00000
+DEAL::0.00000, 0.146447, 0.250000
+DEAL::0.00000, 0.146447, 0.500000
+DEAL::0.00000, 0.146447, 0.750000
+DEAL::0.00000, 0.146447, 1.00000
+DEAL::0.00000, 0.292893, 0.125000
+DEAL::0.00000, 0.292893, 0.375000
+DEAL::0.00000, 0.292893, 0.625000
+DEAL::0.00000, 0.292893, 0.875000
+DEAL::0.146447, -0.292893, 0.00000
+DEAL::0.146447, -0.292893, 0.250000
+DEAL::0.146447, -0.292893, 0.500000
+DEAL::0.146447, -0.292893, 0.750000
+DEAL::0.146447, -0.292893, 1.00000
+DEAL::0.146447, 0.00000, 0.00000
+DEAL::0.146447, 0.00000, 0.250000
+DEAL::0.146447, 0.00000, 0.500000
+DEAL::0.146447, 0.00000, 0.750000
+DEAL::0.146447, 0.00000, 1.00000
+DEAL::0.146447, 0.292893, 0.00000
+DEAL::0.146447, 0.292893, 0.250000
+DEAL::0.146447, 0.292893, 0.500000
+DEAL::0.146447, 0.292893, 0.750000
+DEAL::0.146447, 0.292893, 1.00000
+DEAL::0.292893, -0.292893, 0.125000
+DEAL::0.292893, -0.292893, 0.375000
+DEAL::0.292893, -0.292893, 0.625000
+DEAL::0.292893, -0.292893, 0.875000
+DEAL::0.292893, -0.146447, 0.00000
+DEAL::0.292893, -0.146447, 0.250000
+DEAL::0.292893, -0.146447, 0.500000
+DEAL::0.292893, -0.146447, 0.750000
+DEAL::0.292893, -0.146447, 1.00000
+DEAL::0.292893, 0.00000, 0.125000
+DEAL::0.292893, 0.00000, 0.375000
+DEAL::0.292893, 0.00000, 0.625000
+DEAL::0.292893, 0.00000, 0.875000
+DEAL::0.292893, 0.146447, 0.00000
+DEAL::0.292893, 0.146447, 0.250000
+DEAL::0.292893, 0.146447, 0.500000
+DEAL::0.292893, 0.146447, 0.750000
+DEAL::0.292893, 0.146447, 1.00000
+DEAL::0.292893, 0.292893, 0.125000
+DEAL::0.292893, 0.292893, 0.375000
+DEAL::0.292893, 0.292893, 0.625000
+DEAL::0.292893, 0.292893, 0.875000

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.