#include <deal.II/base/config.h>
+#include <deal.II/base/array_view.h>
#include <deal.II/base/geometry_info.h>
void
serialize(Archive &archive, const unsigned int /*version*/);
+ /**
+ * Return a vector of faces a @p vertex belongs to.
+ */
+ ArrayView<const unsigned int>
+ faces_for_given_vertex(const unsigned int vertex) const;
+
private:
/**
* The variable that stores what this object actually corresponds to.
+ inline ArrayView<const unsigned int>
+ Type::faces_for_given_vertex(const unsigned int vertex) const
+ {
+ if (kind == Type::Line)
+ {
+ AssertIndexRange(vertex, GeometryInfo<1>::vertices_per_cell);
+ return {&GeometryInfo<2>::vertex_to_face[vertex][0], 1};
+ }
+ else if (kind == Type::Quad)
+ {
+ AssertIndexRange(vertex, GeometryInfo<2>::vertices_per_cell);
+ return {&GeometryInfo<2>::vertex_to_face[vertex][0], 2};
+ }
+ else if (kind == Type::Hex)
+ {
+ AssertIndexRange(vertex, GeometryInfo<3>::vertices_per_cell);
+ return {&GeometryInfo<3>::vertex_to_face[vertex][0], 3};
+ }
+ else if (kind == Type::Tri)
+ {
+ AssertIndexRange(vertex, 3);
+ static const std::array<std::array<unsigned int, 2>, 3> table = {
+ {{{0, 2}}, {{0, 1}}, {{1, 2}}}};
+
+ return table[vertex];
+ }
+ else if (kind == Type::Tet)
+ {
+ AssertIndexRange(vertex, 4);
+ static const std::array<std::array<unsigned int, 3>, 4> table = {
+ {{{0, 1, 2}}, {{0, 1, 3}}, {{0, 2, 3}}, {{1, 2, 3}}}};
+
+ return table[vertex];
+ }
+ else if (kind == Type::Wedge)
+ {
+ AssertIndexRange(vertex, 6);
+ static const std::array<std::array<unsigned int, 3>, 6> table = {
+ {{{0, 2, 4}},
+ {{0, 2, 3}},
+ {{0, 3, 4}},
+ {{1, 2, 4}},
+ {{1, 2, 3}},
+ {{1, 3, 4}}}};
+
+ return table[vertex];
+ }
+ else if (kind == Type::Pyramid)
+ {
+ AssertIndexRange(vertex, 5);
+ static const unsigned int X = numbers::invalid_unsigned_int;
+ static const std::array<std::array<unsigned int, 4>, 5> table = {
+ {{{0, 1, 3, X}},
+ {{0, 2, 3, X}},
+ {{0, 1, 4, X}},
+ {{0, 2, 4, X}},
+ {{1, 2, 3, 4}}}};
+
+ return {&table[vertex][0], vertex == 4 ? 4u : 3u};
+ }
+
+ Assert(false, ExcNotImplemented());
+
+ return {};
+ }
+
+
+
/**
* Return the dimension of the given reference-cell type @p type.
*/
// (possibly) coarser neighbor. if this is the case,
// then we need to also add this neighbor
if (dim >= 2)
- for (unsigned int vface = 0; vface < dim; vface++)
- {
- const unsigned int face =
- GeometryInfo<dim>::vertex_to_face[v][vface]; // TODO
-
- if (!cell->at_boundary(face) &&
- cell->neighbor(face)->is_active())
- {
- // there is a (possibly) coarser cell behind a
- // face to which the vertex belongs. the
- // vertex we are looking at is then either a
- // vertex of that coarser neighbor, or it is a
- // hanging node on one of the faces of that
- // cell. in either case, it is adjacent to the
- // vertex, so add it to the list as well (if
- // the cell was already in the list then the
- // std::set makes sure that we get it only
- // once)
- adjacent_cells.insert(cell->neighbor(face));
- }
- }
+ for (const auto face :
+ cell->reference_cell_type().faces_for_given_vertex(v))
+ if (!cell->at_boundary(face) &&
+ cell->neighbor(face)->is_active())
+ {
+ // there is a (possibly) coarser cell behind a
+ // face to which the vertex belongs. the
+ // vertex we are looking at is then either a
+ // vertex of that coarser neighbor, or it is a
+ // hanging node on one of the faces of that
+ // cell. in either case, it is adjacent to the
+ // vertex, so add it to the list as well (if
+ // the cell was already in the list then the
+ // std::set makes sure that we get it only
+ // once)
+ adjacent_cells.insert(cell->neighbor(face));
+ }
// in any case, we have found a cell, so go to the next cell
goto next_cell;
// (possibly) coarser neighbor. if this is the case,
// then we need to also add this neighbor
if (dim >= 2)
- for (unsigned int vface = 0; vface < dim; vface++)
- {
- const unsigned int face =
- GeometryInfo<dim>::vertex_to_face[v][vface];
-
- if (!cell->at_boundary(face) &&
- cell->neighbor(face)->is_active())
- {
- // there is a (possibly) coarser cell behind a
- // face to which the vertex belongs. the
- // vertex we are looking at is then either a
- // vertex of that coarser neighbor, or it is a
- // hanging node on one of the faces of that
- // cell. in either case, it is adjacent to the
- // vertex, so add it to the list as well (if
- // the cell was already in the list then the
- // std::set makes sure that we get it only
- // once)
- adjacent_cells.insert(cell->neighbor(face));
- }
- }
+ for (const auto face :
+ cell->reference_cell_type().faces_for_given_vertex(v))
+ if (!cell->at_boundary(face) &&
+ cell->neighbor(face)->is_active())
+ {
+ // there is a (possibly) coarser cell behind a
+ // face to which the vertex belongs. the
+ // vertex we are looking at is then either a
+ // vertex of that coarser neighbor, or it is a
+ // hanging node on one of the faces of that
+ // cell. in either case, it is adjacent to the
+ // vertex, so add it to the list as well (if
+ // the cell was already in the list then the
+ // std::set makes sure that we get it only
+ // once)
+ adjacent_cells.insert(cell->neighbor(face));
+ }
// in any case, we have found a cell, so go to the next cell
goto next_cell;
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2021 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 ReferenceCell::Kind::faces_for_given_vertex().
+
+
+#include <deal.II/grid/reference_cell.h>
+
+#include <deal.II/simplex/quadrature_lib.h>
+
+#include "../tests.h"
+
+using namespace dealii;
+
+template <int dim>
+void
+test(const ReferenceCell::Type &reference_cell)
+{
+ const auto kind = ReferenceCell::Type(reference_cell);
+ const auto &info = ReferenceCell::internal::Info::get_cell(reference_cell);
+
+ for (const auto v : info.vertex_indices())
+ {
+ deallog << v << ": ";
+ for (const auto i : kind.faces_for_given_vertex(v))
+ deallog << i << " ";
+ deallog << std::endl;
+ }
+ deallog << std::endl;
+}
+
+
+
+int
+main()
+{
+ initlog();
+
+ test<2>(ReferenceCell::Type::Line);
+ test<2>(ReferenceCell::Type::Tri);
+ test<2>(ReferenceCell::Type::Quad);
+ test<3>(ReferenceCell::Type::Tet);
+ test<3>(ReferenceCell::Type::Pyramid);
+ test<3>(ReferenceCell::Type::Wedge);
+ test<3>(ReferenceCell::Type::Hex);
+}
--- /dev/null
+
+DEAL::0: 0
+DEAL::1: 1
+DEAL::
+DEAL::0: 0 2
+DEAL::1: 0 1
+DEAL::2: 1 2
+DEAL::
+DEAL::0: 0 2
+DEAL::1: 1 2
+DEAL::2: 0 3
+DEAL::3: 1 3
+DEAL::
+DEAL::0: 0 1 2
+DEAL::1: 0 1 3
+DEAL::2: 0 2 3
+DEAL::3: 1 2 3
+DEAL::
+DEAL::0: 0 1 3
+DEAL::1: 0 2 3
+DEAL::2: 0 1 4
+DEAL::3: 0 2 4
+DEAL::4: 1 2 3 4
+DEAL::
+DEAL::0: 0 2 4
+DEAL::1: 0 2 3
+DEAL::2: 0 3 4
+DEAL::3: 1 2 4
+DEAL::4: 1 2 3
+DEAL::5: 1 3 4
+DEAL::
+DEAL::0: 0 2 4
+DEAL::1: 1 2 4
+DEAL::2: 0 3 4
+DEAL::3: 1 3 4
+DEAL::4: 0 2 5
+DEAL::5: 1 2 5
+DEAL::6: 0 3 5
+DEAL::7: 1 3 5
+DEAL::