--- /dev/null
+New: The function GridGenerator::convert_hypercube_to_simplex_mesh allows to
+convert a given triangulation based on quadrilaterals (2D) or hexaedra (3D)
+to a triangulation based on simplices or tetraedra, respectively. Thereby,
+material_IDs and boundary_IDs are inherited from the given triangulation.
+<br>
+(Elias Dejene, Peter Munch, 2020/11/23)
flatten_triangulation(const Triangulation<dim, spacedim1> &in_tria,
Triangulation<dim, spacedim2> & out_tria);
+ /**
+ * Convert a triangulation consisting only of hypercube cells
+ * (quadrilaterals, hexahedrons) to a triangulation only consisting of
+ * simplices (triangles, tetrahedra).
+ *
+ * As an example, the following image shows how a set of three hexahedra
+ * meshing one eighths of a sphere are subdivided into tetrahedra, and how
+ * the curved surface is taken into account. Colors indicate how boundary
+ * indicators are inherited:
+ * @image html "convert_hypercube_to_simplex_mesh_visualization.png"
+ *
+ * Thereby, in 2D a quadrilateral cell is converted into 8 triangle cells,
+ * whereas in 3D a hexahedron cell is converted into 24 tetrahedra cells.
+ *
+ * Material ID and boundary IDs will be inherited after conversion.
+ *
+ * @param in_tria The triangulation containing hex elements.
+ * @param out_tria The converted triangulation containing tet elements.
+ */
+ template <int dim, int spacedim>
+ void
+ convert_hypercube_to_simplex_mesh(const Triangulation<dim, spacedim> &in_tria,
+ Triangulation<dim, spacedim> &out_tria);
+
+ /**
+ * Specialization of the above function for 1D: simply copy triangulation.
+ */
+ template <int spacedim>
+ void
+ convert_hypercube_to_simplex_mesh(const Triangulation<1, spacedim> &in_tria,
+ Triangulation<1, spacedim> & out_tria);
+
/**
* Namespace Airfoil contains classes and functions in order to create a
#include <deal.II/grid/tria_accessor.h>
#include <deal.II/grid/tria_iterator.h>
+#include <array>
#include <cmath>
#include <limits>
+ template <int dim, int spacedim>
+ void
+ convert_hypercube_to_simplex_mesh(const Triangulation<dim, spacedim> &in_tria,
+ Triangulation<dim, spacedim> &out_tria)
+ {
+ Assert(in_tria.n_global_levels() == 1,
+ ExcMessage("Number of global levels has to be 1."));
+
+ Assert(dim > 1, ExcNotImplemented());
+
+ // static tables with the definitions of cells and boundary faces for 2D
+ // and 3D
+ static const std::array<std::array<unsigned int, 3>, 8> table_2D_cell = {
+ {{{0, 6, 4}},
+ {{8, 4, 6}},
+ {{8, 6, 5}},
+ {{1, 5, 6}},
+ {{2, 4, 7}},
+ {{8, 7, 4}},
+ {{8, 5, 7}},
+ {{3, 7, 5}}}};
+
+ static const std::array<std::array<unsigned int, 4>, 24> table_3D_cell = {
+ {{{0, 1, 12, 10}}, {{2, 3, 11, 12}}, {{7, 6, 11, 13}},
+ {{5, 4, 13, 10}}, {{0, 2, 8, 12}}, {{4, 6, 13, 8}},
+ {{5, 13, 7, 9}}, {{1, 9, 3, 12}}, {{0, 8, 4, 10}},
+ {{1, 5, 9, 10}}, {{3, 7, 11, 9}}, {{2, 6, 8, 11}},
+ {{12, 13, 10, 9}}, {{12, 13, 9, 11}}, {{12, 13, 11, 8}},
+ {{12, 13, 8, 10}}, {{13, 8, 10, 4}}, {{13, 10, 9, 5}},
+ {{13, 9, 11, 7}}, {{13, 11, 8, 6}}, {{10, 12, 9, 1}},
+ {{9, 12, 11, 3}}, {{11, 12, 8, 2}}, {{8, 12, 10, 0}}}};
+
+ static const std::array<std::array<std::array<unsigned int, 2>, 2>, 4>
+ table_2D_boundary_faces = {{{{{{0, 4}}, {{4, 2}}}},
+ {{{{1, 5}}, {{5, 3}}}},
+ {{{{0, 6}}, {{6, 1}}}},
+ {{{{2, 7}}, {{7, 3}}}}}};
+
+ static const std::array<std::array<std::array<unsigned int, 3>, 4>, 6>
+ table_3D_boundary_faces = {
+ {{{{{0, 4, 8}}, {{4, 8, 6}}, {{8, 6, 2}}, {{0, 2, 8}}}},
+ {{{{1, 3, 9}}, {{3, 9, 7}}, {{9, 7, 5}}, {{1, 9, 5}}}},
+ {{{{0, 1, 10}}, {{1, 10, 5}}, {{10, 5, 4}}, {{0, 10, 4}}}},
+ {{{{2, 3, 11}}, {{3, 11, 7}}, {{11, 7, 6}}, {{2, 11, 6}}}},
+ {{{{0, 1, 12}}, {{1, 12, 3}}, {{12, 3, 2}}, {{0, 12, 2}}}},
+ {{{{4, 5, 13}}, {{5, 13, 7}}, {{13, 7, 6}}, {{4, 13, 6}}}}}};
+
+ std::vector<Point<spacedim>> vertices;
+ std::vector<CellData<dim>> cells;
+ SubCellData subcell_data;
+
+ // store for each vertex and face the assigned index so that we only
+ // assign them a value once
+ std::vector<unsigned int> old_to_new_vertex_indices(
+ in_tria.n_vertices(), numbers::invalid_unsigned_int);
+ std::vector<unsigned int> face_to_new_vertex_indices(
+ in_tria.n_faces(), numbers::invalid_unsigned_int);
+
+ // We first have to create all of the new vertices. To do this, we loop over
+ // all cells and on each cell
+ // (i) copy the existing vertex locations (and record their new indices in
+ // the 'old_to_new_vertex_indices' vector),
+ // (ii) create new midpoint vertex locations for each face (and record their
+ // new indices in the 'face_to_new_vertex_indices' vector),
+ // (iii) create new midpoint vertex locations for each cell (dim = 2 only)
+ for (const auto &cell : in_tria)
+ {
+ // temporary array storing the global indices of each cell entity in the
+ // sequence: vertices, edges/faces, cell
+ std::array<unsigned int, dim == 2 ? 9 : 14> local_vertex_indices;
+
+ // (i) copy the existing vertex locations
+ for (const auto v : cell.vertex_indices())
+ {
+ const auto v_global = cell.vertex_index(v);
+
+ if (old_to_new_vertex_indices[v_global] ==
+ numbers::invalid_unsigned_int)
+ {
+ old_to_new_vertex_indices[v_global] = vertices.size();
+ vertices.push_back(cell.vertex(v));
+ }
+
+ local_vertex_indices[v] = old_to_new_vertex_indices[v_global];
+ }
+
+ // (ii) create new midpoint vertex locations for each face
+ for (const auto f : cell.face_indices())
+ {
+ const auto f_global = cell.face_index(f);
+
+ if (face_to_new_vertex_indices[f_global] ==
+ numbers::invalid_unsigned_int)
+ {
+ face_to_new_vertex_indices[f_global] = vertices.size();
+ vertices.push_back(
+ cell.face(f)->center(/*respect_manifold*/ true));
+ }
+
+ local_vertex_indices[cell.n_vertices() + f] =
+ face_to_new_vertex_indices[f_global];
+ }
+
+ // (iii) create new midpoint vertex locations for each cell
+ if (dim == 2)
+ {
+ local_vertex_indices[cell.n_vertices() + cell.n_faces()] =
+ vertices.size();
+ vertices.push_back(cell.center(/*respect_manifold*/ true));
+ }
+
+ // helper function for creating cells and subcells
+ const auto add_cell = [&](const unsigned int struct_dim,
+ const auto & index_vertices,
+ const unsigned int material_or_boundary_id) {
+ if (struct_dim == dim) // cells
+ {
+ CellData<dim> cell_data(index_vertices.size());
+ for (unsigned int i = 0; i < index_vertices.size(); ++i)
+ {
+ cell_data.vertices[i] =
+ local_vertex_indices[index_vertices[i]];
+ cell_data.material_id = material_or_boundary_id;
+ }
+ cells.push_back(cell_data);
+ }
+ else if (dim == 2 && struct_dim == 1) // an edge of a triangle
+ {
+ CellData<1> boundary_line(2);
+ boundary_line.boundary_id = material_or_boundary_id;
+ for (unsigned int i = 0; i < index_vertices.size(); ++i)
+ {
+ boundary_line.vertices[i] =
+ local_vertex_indices[index_vertices[i]];
+ }
+ subcell_data.boundary_lines.push_back(boundary_line);
+ }
+ else if (dim == 3 && struct_dim == 2) // a face of a thetrahedron
+ {
+ CellData<2> boundary_quad(3);
+ boundary_quad.boundary_id = material_or_boundary_id;
+ for (unsigned int i = 0; i < index_vertices.size(); ++i)
+ {
+ boundary_quad.vertices[i] =
+ local_vertex_indices[index_vertices[i]];
+ }
+ subcell_data.boundary_quads.push_back(boundary_quad);
+ }
+ else
+ {
+ Assert(false, ExcNotImplemented());
+ }
+ };
+
+ const auto current_material_id = cell.material_id();
+
+ // create cells one by one
+
+ if (dim == 2)
+ for (const auto &cell : table_2D_cell)
+ add_cell(dim, cell, current_material_id);
+ else if (dim == 3)
+ for (const auto &cell : table_3D_cell)
+ add_cell(dim, cell, current_material_id);
+ else
+ Assert(false, ExcNotImplemented());
+
+ // Set up sub-cell data.
+ for (const auto f : cell.face_indices())
+ if (cell.face(f)->boundary_id() !=
+ numbers::internal_face_boundary_id) // face at boundary
+ {
+ const auto current_boundary_id = cell.face(f)->boundary_id();
+ if (dim == 2) // 2D boundary_lines
+ for (const auto &face : table_2D_boundary_faces[f])
+ add_cell(1, face, current_boundary_id);
+ else if (dim == 3) // 3D boundary_quads
+ for (const auto &face : table_3D_boundary_faces[f])
+ add_cell(2, face, current_boundary_id);
+ else
+ Assert(false, ExcNotImplemented());
+ }
+ }
+
+ out_tria.create_triangulation(vertices, cells, subcell_data);
+ }
+
+
+
+ template <int spacedim>
+ void
+ convert_hypercube_to_simplex_mesh(const Triangulation<1, spacedim> &in_tria,
+ Triangulation<1, spacedim> & out_tria)
+ {
+ out_tria.copy_triangulation(in_tria);
+ return;
+ }
+
+
+
template <template <int, int> class MeshType, int dim, int spacedim>
# ifndef _MSC_VER
std::map<typename MeshType<dim - 1, spacedim>::cell_iterator,
#endif
\}
}
+
+for (deal_II_dimension : DIMENSIONS; deal_II_space_dimension : SPACE_DIMENSIONS)
+ {
+ namespace GridGenerator
+ \{
+#if deal_II_dimension <= deal_II_space_dimension
+ template void
+ convert_hypercube_to_simplex_mesh(
+ const Triangulation<deal_II_dimension, deal_II_space_dimension> &,
+ Triangulation<deal_II_dimension, deal_II_space_dimension> &);
+#endif
+ \}
+ }
--- /dev/null
+/* ---------------------------------------------------------------------
+ *
+ * Copyright (C) 2020 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.
+ *
+ * ---------------------------------------------------------------------
+ *
+ * Test function GridGenerator::convert_hypercube_to_simplex_mesh() in 2D
+ * and 3D (dim = spacedim = 2, 3) on a quarter_hyper_ball() triangulation and
+ * for dim = 2 < spacedim = 3 on a subdivided_hyper_cube() triangulation.
+ * Therefore, the output is written in .vtk format for entire triangulations
+ * as well as surface-only triangulations.
+ */
+
+#include <deal.II/grid/grid_generator.h>
+#include <deal.II/grid/grid_out.h>
+#include <deal.II/grid/grid_tools.h>
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/tria_iterator.h>
+
+#include <deal.II/simplex/grid_generator.h>
+
+#include <fstream>
+#include <iostream>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim, int spacedim>
+void
+create_triangulation(Triangulation<dim, spacedim> &triangulation)
+{
+ GridGenerator::subdivided_hyper_cube(triangulation, 4);
+}
+
+template <int dim>
+void
+create_triangulation(Triangulation<dim, dim> &triangulation)
+{
+ GridGenerator::quarter_hyper_ball(triangulation);
+}
+
+template <int dim, int spacedim>
+void
+check_file() // for dim = spaceim
+{
+ Triangulation<dim, spacedim> in_tria, out_tria;
+ create_triangulation(in_tria);
+
+ // make each cell a different material id
+ unsigned int m_id = 0;
+ for (const auto &cell : in_tria)
+ {
+ cell.set_material_id(m_id++);
+ }
+
+ // set different boundary ids and output
+ unsigned int b_id = 0;
+ for (const auto &cell : in_tria)
+ {
+ for (const auto f : cell.face_indices())
+ {
+ if (cell.face(f)->at_boundary())
+ {
+ cell.face(f)->set_boundary_id(b_id);
+ b_id++;
+ }
+ }
+ }
+
+ GridGenerator::convert_hypercube_to_simplex_mesh(in_tria, out_tria);
+
+ // write 2 outputs (total mesh and only surface mesh)
+ const auto grid_out = [](const auto &tria,
+ const bool surface_mesh_only = false) {
+ GridOutFlags::Vtk flags;
+
+ if (surface_mesh_only)
+ {
+ flags.output_cells = false;
+ flags.output_faces = true;
+ flags.output_edges = false;
+ flags.output_only_relevant = false;
+ }
+
+ GridOut grid_out;
+ grid_out.set_flags(flags);
+
+ grid_out.write_vtk(tria, deallog.get_file_stream());
+ };
+
+ grid_out(out_tria); // total mesh
+ grid_out(out_tria, true); // only surface mesh
+
+ deallog << "OK!" << std::endl;
+}
+
+
+int
+main()
+{
+ initlog();
+ // TRIANGULAR ELEMENTS
+ // dim = spacedim = 2
+ deallog.push(
+ "2D: conversion triangulation with quad elements to tri elements: ");
+ check_file<2, 2>();
+ deallog.pop();
+
+ // TETRAHEDRAL ELEMENTS
+ // dim = 2, spacedim = 2
+ deallog.push(
+ "2D: conversion triangulation with quad elements to tri elements: ");
+ check_file<2, 3>();
+ deallog.pop();
+
+ // TETRAHEDRAL ELEMENTS
+ // dim = spacedim = 3
+ deallog.push(
+ "3D: conversion triangulation with tet elements to hex elements: ");
+ check_file<3, 3>();
+ deallog.pop();
+}
--- /dev/null
+
+# vtk DataFile Version 3.0
+Triangulation generated with deal.II
+ASCII
+DATASET UNSTRUCTURED_GRID
+POINTS 19 double
+0.00000 0.00000 0
+0.556470 0.00000 0
+0.00000 0.556470 0
+0.428830 0.428830 0
+0.00000 0.278235 0
+0.492650 0.214415 0
+0.278235 0.00000 0
+0.214415 0.492650 0
+0.246325 0.246325 0
+1.00000 0.00000 0
+0.707107 0.707107 0
+0.778235 0.00000 0
+0.567968 0.567968 0
+0.923880 0.382683 0
+0.681892 0.287625 0
+0.00000 1.00000 0
+0.382683 0.923880 0
+0.00000 0.778235 0
+0.287625 0.681892 0
+
+CELLS 34 126
+3 0 6 4
+3 8 4 6
+3 8 6 5
+3 1 5 6
+3 2 4 7
+3 8 7 4
+3 8 5 7
+3 3 7 5
+3 9 13 11
+3 14 11 13
+3 14 13 12
+3 10 12 13
+3 1 11 5
+3 14 5 11
+3 14 12 5
+3 3 5 12
+3 15 17 16
+3 18 16 17
+3 18 17 7
+3 2 7 17
+3 10 16 12
+3 18 12 16
+3 18 7 12
+3 3 12 7
+2 0 6
+2 1 11
+2 6 1
+2 9 13
+2 10 16
+2 11 9
+2 13 10
+2 15 17
+2 16 15
+2 17 2
+
+CELL_TYPES 34
+5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5
+3 3 3 3 3 3 3 3 3 3
+
+
+CELL_DATA 34
+SCALARS MaterialID int 1
+LOOKUP_TABLE default
+0 0 0 0 0 0 0 0 1 1 1 1 1 1 1 1 2 2 2 2 2 2 2 2
+1 2 1 3 4 2 3 5 4 5
+
+
+SCALARS ManifoldID int 1
+LOOKUP_TABLE default
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+
+# vtk DataFile Version 3.0
+Triangulation generated with deal.II
+ASCII
+DATASET UNSTRUCTURED_GRID
+POINTS 19 double
+0.00000 0.00000 0
+0.556470 0.00000 0
+0.00000 0.556470 0
+0.428830 0.428830 0
+0.00000 0.278235 0
+0.492650 0.214415 0
+0.278235 0.00000 0
+0.214415 0.492650 0
+0.246325 0.246325 0
+1.00000 0.00000 0
+0.707107 0.707107 0
+0.778235 0.00000 0
+0.567968 0.567968 0
+0.923880 0.382683 0
+0.681892 0.287625 0
+0.00000 1.00000 0
+0.382683 0.923880 0
+0.00000 0.778235 0
+0.287625 0.681892 0
+
+CELLS 12 36
+2 0 6
+2 1 11
+2 2 4
+2 4 0
+2 6 1
+2 9 13
+2 10 16
+2 11 9
+2 13 10
+2 15 17
+2 16 15
+2 17 2
+
+CELL_TYPES 12
+3 3 3 3 3 3 3 3 3 3 3 3
+
+
+CELL_DATA 12
+SCALARS MaterialID int 1
+LOOKUP_TABLE default
+1 2 0 0 1 3 4 2 3 5 4 5
+
+
+SCALARS ManifoldID int 1
+LOOKUP_TABLE default
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+DEAL:2D: conversion triangulation with quad elements to tri elements: ::OK!
+# vtk DataFile Version 3.0
+Triangulation generated with deal.II
+ASCII
+DATASET UNSTRUCTURED_GRID
+POINTS 81 double
+0.00000 0.00000 0.00000
+0.250000 0.00000 0.00000
+0.00000 0.250000 0.00000
+0.250000 0.250000 0.00000
+0.00000 0.125000 0.00000
+0.250000 0.125000 0.00000
+0.125000 0.00000 0.00000
+0.125000 0.250000 0.00000
+0.125000 0.125000 0.00000
+0.500000 0.00000 0.00000
+0.500000 0.250000 0.00000
+0.500000 0.125000 0.00000
+0.375000 0.00000 0.00000
+0.375000 0.250000 0.00000
+0.375000 0.125000 0.00000
+0.750000 0.00000 0.00000
+0.750000 0.250000 0.00000
+0.750000 0.125000 0.00000
+0.625000 0.00000 0.00000
+0.625000 0.250000 0.00000
+0.625000 0.125000 0.00000
+1.00000 0.00000 0.00000
+1.00000 0.250000 0.00000
+1.00000 0.125000 0.00000
+0.875000 0.00000 0.00000
+0.875000 0.250000 0.00000
+0.875000 0.125000 0.00000
+0.00000 0.500000 0.00000
+0.250000 0.500000 0.00000
+0.00000 0.375000 0.00000
+0.250000 0.375000 0.00000
+0.125000 0.500000 0.00000
+0.125000 0.375000 0.00000
+0.500000 0.500000 0.00000
+0.500000 0.375000 0.00000
+0.375000 0.500000 0.00000
+0.375000 0.375000 0.00000
+0.750000 0.500000 0.00000
+0.750000 0.375000 0.00000
+0.625000 0.500000 0.00000
+0.625000 0.375000 0.00000
+1.00000 0.500000 0.00000
+1.00000 0.375000 0.00000
+0.875000 0.500000 0.00000
+0.875000 0.375000 0.00000
+0.00000 0.750000 0.00000
+0.250000 0.750000 0.00000
+0.00000 0.625000 0.00000
+0.250000 0.625000 0.00000
+0.125000 0.750000 0.00000
+0.125000 0.625000 0.00000
+0.500000 0.750000 0.00000
+0.500000 0.625000 0.00000
+0.375000 0.750000 0.00000
+0.375000 0.625000 0.00000
+0.750000 0.750000 0.00000
+0.750000 0.625000 0.00000
+0.625000 0.750000 0.00000
+0.625000 0.625000 0.00000
+1.00000 0.750000 0.00000
+1.00000 0.625000 0.00000
+0.875000 0.750000 0.00000
+0.875000 0.625000 0.00000
+0.00000 1.00000 0.00000
+0.250000 1.00000 0.00000
+0.00000 0.875000 0.00000
+0.250000 0.875000 0.00000
+0.125000 1.00000 0.00000
+0.125000 0.875000 0.00000
+0.500000 1.00000 0.00000
+0.500000 0.875000 0.00000
+0.375000 1.00000 0.00000
+0.375000 0.875000 0.00000
+0.750000 1.00000 0.00000
+0.750000 0.875000 0.00000
+0.625000 1.00000 0.00000
+0.625000 0.875000 0.00000
+1.00000 1.00000 0.00000
+1.00000 0.875000 0.00000
+0.875000 1.00000 0.00000
+0.875000 0.875000 0.00000
+
+CELLS 158 602
+3 0 6 4
+3 8 4 6
+3 8 6 5
+3 1 5 6
+3 2 4 7
+3 8 7 4
+3 8 5 7
+3 3 7 5
+3 1 12 5
+3 14 5 12
+3 14 12 11
+3 9 11 12
+3 3 5 13
+3 14 13 5
+3 14 11 13
+3 10 13 11
+3 9 18 11
+3 20 11 18
+3 20 18 17
+3 15 17 18
+3 10 11 19
+3 20 19 11
+3 20 17 19
+3 16 19 17
+3 15 24 17
+3 26 17 24
+3 26 24 23
+3 21 23 24
+3 16 17 25
+3 26 25 17
+3 26 23 25
+3 22 25 23
+3 2 7 29
+3 32 29 7
+3 32 7 30
+3 3 30 7
+3 27 29 31
+3 32 31 29
+3 32 30 31
+3 28 31 30
+3 3 13 30
+3 36 30 13
+3 36 13 34
+3 10 34 13
+3 28 30 35
+3 36 35 30
+3 36 34 35
+3 33 35 34
+3 10 19 34
+3 40 34 19
+3 40 19 38
+3 16 38 19
+3 33 34 39
+3 40 39 34
+3 40 38 39
+3 37 39 38
+3 16 25 38
+3 44 38 25
+3 44 25 42
+3 22 42 25
+3 37 38 43
+3 44 43 38
+3 44 42 43
+3 41 43 42
+3 27 31 47
+3 50 47 31
+3 50 31 48
+3 28 48 31
+3 45 47 49
+3 50 49 47
+3 50 48 49
+3 46 49 48
+3 28 35 48
+3 54 48 35
+3 54 35 52
+3 33 52 35
+3 46 48 53
+3 54 53 48
+3 54 52 53
+3 51 53 52
+3 33 39 52
+3 58 52 39
+3 58 39 56
+3 37 56 39
+3 51 52 57
+3 58 57 52
+3 58 56 57
+3 55 57 56
+3 37 43 56
+3 62 56 43
+3 62 43 60
+3 41 60 43
+3 55 56 61
+3 62 61 56
+3 62 60 61
+3 59 61 60
+3 45 49 65
+3 68 65 49
+3 68 49 66
+3 46 66 49
+3 63 65 67
+3 68 67 65
+3 68 66 67
+3 64 67 66
+3 46 53 66
+3 72 66 53
+3 72 53 70
+3 51 70 53
+3 64 66 71
+3 72 71 66
+3 72 70 71
+3 69 71 70
+3 51 57 70
+3 76 70 57
+3 76 57 74
+3 55 74 57
+3 69 70 75
+3 76 75 70
+3 76 74 75
+3 73 75 74
+3 55 61 74
+3 80 74 61
+3 80 61 78
+3 59 78 61
+3 73 74 79
+3 80 79 74
+3 80 78 79
+3 77 79 78
+2 0 6
+2 1 12
+2 6 1
+2 9 18
+2 12 9
+2 15 24
+2 18 15
+2 21 23
+2 22 42
+2 23 22
+2 24 21
+2 27 29
+2 29 2
+2 41 60
+2 42 41
+2 45 47
+2 47 27
+2 59 78
+2 60 59
+2 63 65
+2 64 67
+2 65 45
+2 67 63
+2 69 71
+2 71 64
+2 73 75
+2 75 69
+2 77 79
+2 78 77
+2 79 73
+
+CELL_TYPES 158
+5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5
+3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3
+
+
+CELL_DATA 158
+SCALARS MaterialID int 1
+LOOKUP_TABLE default
+0 0 0 0 0 0 0 0 1 1 1 1 1 1 1 1 2 2 2 2 2 2 2 2 3 3 3 3 3 3 3 3 4 4 4 4 4 4 4 4 5 5 5 5 5 5 5 5 6 6 6 6 6 6 6 6 7 7 7 7 7 7 7 7 8 8 8 8 8 8 8 8 9 9 9 9 9 9 9 9 10 10 10 10 10 10 10 10 11 11 11 11 11 11 11 11 12 12 12 12 12 12 12 12 13 13 13 13 13 13 13 13 14 14 14 14 14 14 14 14 15 15 15 15 15 15 15 15
+1 2 1 3 2 5 3 4 7 4 5 6 6 9 7 8 8 14 9 10 11 10 11 12 12 13 13 15 14 15
+
+
+SCALARS ManifoldID int 1
+LOOKUP_TABLE default
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+
+# vtk DataFile Version 3.0
+Triangulation generated with deal.II
+ASCII
+DATASET UNSTRUCTURED_GRID
+POINTS 81 double
+0.00000 0.00000 0.00000
+0.250000 0.00000 0.00000
+0.00000 0.250000 0.00000
+0.250000 0.250000 0.00000
+0.00000 0.125000 0.00000
+0.250000 0.125000 0.00000
+0.125000 0.00000 0.00000
+0.125000 0.250000 0.00000
+0.125000 0.125000 0.00000
+0.500000 0.00000 0.00000
+0.500000 0.250000 0.00000
+0.500000 0.125000 0.00000
+0.375000 0.00000 0.00000
+0.375000 0.250000 0.00000
+0.375000 0.125000 0.00000
+0.750000 0.00000 0.00000
+0.750000 0.250000 0.00000
+0.750000 0.125000 0.00000
+0.625000 0.00000 0.00000
+0.625000 0.250000 0.00000
+0.625000 0.125000 0.00000
+1.00000 0.00000 0.00000
+1.00000 0.250000 0.00000
+1.00000 0.125000 0.00000
+0.875000 0.00000 0.00000
+0.875000 0.250000 0.00000
+0.875000 0.125000 0.00000
+0.00000 0.500000 0.00000
+0.250000 0.500000 0.00000
+0.00000 0.375000 0.00000
+0.250000 0.375000 0.00000
+0.125000 0.500000 0.00000
+0.125000 0.375000 0.00000
+0.500000 0.500000 0.00000
+0.500000 0.375000 0.00000
+0.375000 0.500000 0.00000
+0.375000 0.375000 0.00000
+0.750000 0.500000 0.00000
+0.750000 0.375000 0.00000
+0.625000 0.500000 0.00000
+0.625000 0.375000 0.00000
+1.00000 0.500000 0.00000
+1.00000 0.375000 0.00000
+0.875000 0.500000 0.00000
+0.875000 0.375000 0.00000
+0.00000 0.750000 0.00000
+0.250000 0.750000 0.00000
+0.00000 0.625000 0.00000
+0.250000 0.625000 0.00000
+0.125000 0.750000 0.00000
+0.125000 0.625000 0.00000
+0.500000 0.750000 0.00000
+0.500000 0.625000 0.00000
+0.375000 0.750000 0.00000
+0.375000 0.625000 0.00000
+0.750000 0.750000 0.00000
+0.750000 0.625000 0.00000
+0.625000 0.750000 0.00000
+0.625000 0.625000 0.00000
+1.00000 0.750000 0.00000
+1.00000 0.625000 0.00000
+0.875000 0.750000 0.00000
+0.875000 0.625000 0.00000
+0.00000 1.00000 0.00000
+0.250000 1.00000 0.00000
+0.00000 0.875000 0.00000
+0.250000 0.875000 0.00000
+0.125000 1.00000 0.00000
+0.125000 0.875000 0.00000
+0.500000 1.00000 0.00000
+0.500000 0.875000 0.00000
+0.375000 1.00000 0.00000
+0.375000 0.875000 0.00000
+0.750000 1.00000 0.00000
+0.750000 0.875000 0.00000
+0.625000 1.00000 0.00000
+0.625000 0.875000 0.00000
+1.00000 1.00000 0.00000
+1.00000 0.875000 0.00000
+0.875000 1.00000 0.00000
+0.875000 0.875000 0.00000
+
+CELLS 32 96
+2 0 6
+2 1 12
+2 2 4
+2 4 0
+2 6 1
+2 9 18
+2 12 9
+2 15 24
+2 18 15
+2 21 23
+2 22 42
+2 23 22
+2 24 21
+2 27 29
+2 29 2
+2 41 60
+2 42 41
+2 45 47
+2 47 27
+2 59 78
+2 60 59
+2 63 65
+2 64 67
+2 65 45
+2 67 63
+2 69 71
+2 71 64
+2 73 75
+2 75 69
+2 77 79
+2 78 77
+2 79 73
+
+CELL_TYPES 32
+3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3
+
+
+CELL_DATA 32
+SCALARS MaterialID int 1
+LOOKUP_TABLE default
+1 2 0 0 1 3 2 5 3 4 7 4 5 6 6 9 7 8 8 14 9 10 11 10 11 12 12 13 13 15 14 15
+
+
+SCALARS ManifoldID int 1
+LOOKUP_TABLE default
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+DEAL:2D: conversion triangulation with quad elements to tri elements: ::OK!
+# vtk DataFile Version 3.0
+Triangulation generated with deal.II
+ASCII
+DATASET UNSTRUCTURED_GRID
+POINTS 33 double
+0.00000 0.00000 0.00000
+0.528000 0.00000 0.00000
+0.00000 0.528000 0.00000
+0.453300 0.453300 0.00000
+0.00000 0.00000 0.528000
+0.453300 0.00000 0.453300
+0.00000 0.453300 0.453300
+0.375200 0.375200 0.375200
+0.00000 0.245325 0.245325
+0.452450 0.207125 0.207125
+0.245325 0.00000 0.245325
+0.207125 0.452450 0.207125
+0.245325 0.245325 0.00000
+0.207125 0.207125 0.452450
+1.00000 0.00000 0.00000
+0.707107 0.707107 0.00000
+0.707107 0.00000 0.707107
+0.577350 0.577350 0.577350
+0.680892 0.00000 0.293743
+0.532197 0.532197 0.239916
+0.856495 0.364976 0.364976
+0.680892 0.293743 0.00000
+0.532197 0.239916 0.532197
+0.00000 1.00000 0.00000
+0.00000 0.707107 0.707107
+0.364976 0.856495 0.364976
+0.00000 0.680892 0.293743
+0.293743 0.680892 0.00000
+0.239916 0.532197 0.532197
+0.00000 0.00000 1.00000
+0.00000 0.293743 0.680892
+0.293743 0.00000 0.680892
+0.364976 0.364976 0.856495
+
+CELLS 140 656
+4 0 1 12 10
+4 2 3 11 12
+4 7 6 11 13
+4 5 4 13 10
+4 0 2 8 12
+4 4 6 13 8
+4 5 13 7 9
+4 1 9 3 12
+4 0 8 4 10
+4 1 5 9 10
+4 3 7 11 9
+4 2 6 8 11
+4 12 13 10 9
+4 12 13 9 11
+4 12 13 11 8
+4 12 13 8 10
+4 13 8 10 4
+4 13 10 9 5
+4 13 9 11 7
+4 13 11 8 6
+4 10 12 9 1
+4 9 12 11 3
+4 11 12 8 2
+4 8 12 10 0
+4 14 15 21 20
+4 1 3 9 21
+4 7 5 9 22
+4 17 16 22 20
+4 14 1 18 21
+4 16 5 22 18
+4 17 22 7 19
+4 15 19 3 21
+4 14 18 16 20
+4 15 17 19 20
+4 3 7 9 19
+4 1 5 18 9
+4 21 22 20 19
+4 21 22 19 9
+4 21 22 9 18
+4 21 22 18 20
+4 22 18 20 16
+4 22 20 19 17
+4 22 19 9 7
+4 22 9 18 5
+4 20 21 19 15
+4 19 21 9 3
+4 9 21 18 1
+4 18 21 20 14
+4 23 2 27 26
+4 15 3 19 27
+4 7 17 19 28
+4 6 24 28 26
+4 23 15 25 27
+4 24 17 28 25
+4 6 28 7 11
+4 2 11 3 27
+4 23 25 24 26
+4 2 6 11 26
+4 3 7 19 11
+4 15 17 25 19
+4 27 28 26 11
+4 27 28 11 19
+4 27 28 19 25
+4 27 28 25 26
+4 28 25 26 24
+4 28 26 11 6
+4 28 11 19 7
+4 28 19 25 17
+4 26 27 11 2
+4 11 27 19 3
+4 19 27 25 15
+4 25 27 26 23
+4 4 5 13 31
+4 6 7 28 13
+4 17 24 28 32
+4 16 29 32 31
+4 4 6 30 13
+4 29 24 32 30
+4 16 32 17 22
+4 5 22 7 13
+4 4 30 29 31
+4 5 16 22 31
+4 7 17 28 22
+4 6 24 30 28
+4 13 32 31 22
+4 13 32 22 28
+4 13 32 28 30
+4 13 32 30 31
+4 32 30 31 29
+4 32 31 22 16
+4 32 22 28 17
+4 32 28 30 24
+4 31 13 22 5
+4 22 13 28 7
+4 28 13 30 6
+4 30 13 31 4
+3 1 0 10
+3 0 1 12
+3 2 0 12
+3 5 1 10
+3 1 5 18
+3 3 2 12
+3 2 3 27
+3 6 2 26
+3 3 1 21
+3 1 3 12
+3 15 3 21
+3 3 15 27
+3 0 4 10
+3 4 6 30
+3 5 4 31
+3 4 5 10
+3 24 6 26
+3 6 24 30
+3 14 1 18
+3 1 14 21
+3 15 14 20
+3 14 15 21
+3 17 15 20
+3 15 17 25
+3 16 5 31
+3 5 16 18
+3 14 16 20
+3 29 16 31
+3 16 29 32
+3 16 17 20
+3 14 18 16
+3 23 2 27
+3 2 23 26
+3 23 15 25
+3 15 23 27
+3 24 17 32
+3 17 24 25
+3 23 24 26
+3 23 25 24
+3 4 29 31
+3 29 24 32
+3 24 29 30
+3 4 30 29
+3 16 32 17
+
+CELL_TYPES 140
+10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10 10
+5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5
+
+
+CELL_DATA 140
+SCALARS MaterialID int 1
+LOOKUP_TABLE default
+0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 0 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 1 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 2 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3 3
+1 2 2 1 3 2 8 7 5 2 5 8 1 9 10 1 7 9 3 5 4 5 4 6 10 3 4 10 11 4 3 8 7 6 8 11 6 7 6 10 11 9 9 11
+
+
+SCALARS ManifoldID int 1
+LOOKUP_TABLE default
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+
+# vtk DataFile Version 3.0
+Triangulation generated with deal.II
+ASCII
+DATASET UNSTRUCTURED_GRID
+POINTS 33 double
+0.00000 0.00000 0.00000
+0.528000 0.00000 0.00000
+0.00000 0.528000 0.00000
+0.453300 0.453300 0.00000
+0.00000 0.00000 0.528000
+0.453300 0.00000 0.453300
+0.00000 0.453300 0.453300
+0.375200 0.375200 0.375200
+0.00000 0.245325 0.245325
+0.452450 0.207125 0.207125
+0.245325 0.00000 0.245325
+0.207125 0.452450 0.207125
+0.245325 0.245325 0.00000
+0.207125 0.207125 0.452450
+1.00000 0.00000 0.00000
+0.707107 0.707107 0.00000
+0.707107 0.00000 0.707107
+0.577350 0.577350 0.577350
+0.680892 0.00000 0.293743
+0.532197 0.532197 0.239916
+0.856495 0.364976 0.364976
+0.680892 0.293743 0.00000
+0.532197 0.239916 0.532197
+0.00000 1.00000 0.00000
+0.00000 0.707107 0.707107
+0.364976 0.856495 0.364976
+0.00000 0.680892 0.293743
+0.293743 0.680892 0.00000
+0.239916 0.532197 0.532197
+0.00000 0.00000 1.00000
+0.00000 0.293743 0.680892
+0.293743 0.00000 0.680892
+0.364976 0.364976 0.856495
+
+CELLS 48 192
+3 1 0 10
+3 0 1 12
+3 0 2 8
+3 2 0 12
+3 5 1 10
+3 1 5 18
+3 3 2 12
+3 2 3 27
+3 6 2 26
+3 2 6 8
+3 3 1 21
+3 1 3 12
+3 15 3 21
+3 3 15 27
+3 0 4 10
+3 6 4 8
+3 4 6 30
+3 5 4 31
+3 4 5 10
+3 24 6 26
+3 6 24 30
+3 0 8 4
+3 14 1 18
+3 1 14 21
+3 15 14 20
+3 14 15 21
+3 17 15 20
+3 15 17 25
+3 16 5 31
+3 5 16 18
+3 14 16 20
+3 29 16 31
+3 16 29 32
+3 16 17 20
+3 14 18 16
+3 23 2 27
+3 2 23 26
+3 23 15 25
+3 15 23 27
+3 24 17 32
+3 17 24 25
+3 23 24 26
+3 23 25 24
+3 4 29 31
+3 29 24 32
+3 24 29 30
+3 4 30 29
+3 16 32 17
+
+CELL_TYPES 48
+5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5 5
+
+
+CELL_DATA 48
+SCALARS MaterialID int 1
+LOOKUP_TABLE default
+1 2 0 2 1 3 2 8 7 0 5 2 5 8 1 0 9 10 1 7 9 0 3 5 4 5 4 6 10 3 4 10 11 4 3 8 7 6 8 11 6 7 6 10 11 9 9 11
+
+
+SCALARS ManifoldID int 1
+LOOKUP_TABLE default
+-1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1 -1
+DEAL:3D: conversion triangulation with tet elements to hex elements: ::OK!