]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Generalize FiniteElement::face_to_cell_index 10730/head
authorPeter Munch <peterrmuench@gmail.com>
Mon, 20 Jul 2020 16:36:37 +0000 (18:36 +0200)
committerPeter Munch <peterrmuench@gmail.com>
Mon, 20 Jul 2020 16:36:37 +0000 (18:36 +0200)
include/deal.II/base/geometry_info.h
include/deal.II/grid/reference_cell.h
source/fe/fe.cc

index 61afc04ccbcbc63e0330b4b31484805581fda857..41de5c5338347235f596c4a40de02824d76d3a52 100644 (file)
@@ -1307,6 +1307,57 @@ struct GeometryInfo<0>
   static std::array<unsigned int, vertices_per_cell>
   vertex_indices();
 
+  /**
+   * Map face vertex number to cell vertex number, i.e. give the cell vertex
+   * number of the <tt>vertex</tt>th vertex of face <tt>face</tt>, e.g.
+   * <tt>GeometryInfo<2>::face_to_cell_vertices(3,0)=2</tt>, see the image
+   * under point N4 in the 2d section of this class's documentation.
+   *
+   * Through the <tt>face_orientation</tt>, <tt>face_flip</tt> and
+   * <tt>face_rotation</tt> arguments this function handles faces oriented in
+   * the standard and non-standard orientation. <tt>face_orientation</tt>
+   * defaults to <tt>true</tt>, <tt>face_flip</tt> and <tt>face_rotation</tt>
+   * default to <tt>false</tt> (standard orientation). In 2d only
+   * <tt>face_flip</tt> is considered. See this
+   * @ref GlossFaceOrientation "glossary"
+   * article for more information.
+   *
+   * As the children of a cell are ordered according to the vertices of the
+   * cell, this call is passed down to the child_cell_on_face() function.
+   * Hence this function is simply a wrapper of child_cell_on_face() giving it
+   * a suggestive name.
+   *
+   * Of course, since this class is for the case `dim==0`, this function
+   * is not implemented.
+   */
+  static unsigned int
+  face_to_cell_vertices(const unsigned int face,
+                        const unsigned int vertex,
+                        const bool         face_orientation = true,
+                        const bool         face_flip        = false,
+                        const bool         face_rotation    = false);
+
+  /**
+   * Map face line number to cell line number, i.e. give the cell line number
+   * of the <tt>line</tt>th line of face <tt>face</tt>, e.g.
+   * <tt>GeometryInfo<3>::face_to_cell_lines(5,0)=4</tt>.
+   *
+   * Through the <tt>face_orientation</tt>, <tt>face_flip</tt> and
+   * <tt>face_rotation</tt> arguments this function handles faces oriented in
+   * the standard and non-standard orientation. <tt>face_orientation</tt>
+   * defaults to <tt>true</tt>, <tt>face_flip</tt> and <tt>face_rotation</tt>
+   * default to <tt>false</tt> (standard orientation) and has no effect in 2d.
+   *
+   * Of course, since this class is for the case `dim==0`, this function
+   * is not implemented.
+   */
+  static unsigned int
+  face_to_cell_lines(const unsigned int face,
+                     const unsigned int line,
+                     const bool         face_orientation = true,
+                     const bool         face_flip        = false,
+                     const bool         face_rotation    = false);
+
   /**
    * Number of vertices each face has. Since this is not useful in one
    * dimension, we provide a useless number (in the hope that a compiler may
@@ -4592,6 +4643,19 @@ GeometryInfo<3>::face_to_cell_lines(const unsigned int face,
 
 
 
+inline unsigned int
+GeometryInfo<0>::face_to_cell_lines(const unsigned int,
+                                    const unsigned int,
+                                    const bool,
+                                    const bool,
+                                    const bool)
+{
+  Assert(false, ExcNotImplemented());
+  return numbers::invalid_unsigned_int;
+}
+
+
+
 template <int dim>
 inline unsigned int
 GeometryInfo<dim>::face_to_cell_lines(const unsigned int,
@@ -4624,6 +4688,19 @@ GeometryInfo<dim>::face_to_cell_vertices(const unsigned int face,
 
 
 
+inline unsigned int
+GeometryInfo<0>::face_to_cell_vertices(const unsigned int,
+                                       const unsigned int,
+                                       const bool,
+                                       const bool,
+                                       const bool)
+{
+  Assert(false, ExcNotImplemented());
+  return numbers::invalid_unsigned_int;
+}
+
+
+
 template <int dim>
 inline Point<dim>
 GeometryInfo<dim>::project_to_unit_cell(const Point<dim> &q)
index 8832704e63af3542c8453213ed2afef34baea433..fbeca6ec1bde2ed45837f4bdce3323ef1b74b742 100644 (file)
@@ -241,6 +241,38 @@ namespace ReferenceCell
 
           return ReferenceCell::Type::Invalid;
         }
+
+        /**
+         * Map face line number to cell line number.
+         */
+        virtual unsigned int
+        face_to_cell_lines(const unsigned int  face,
+                           const unsigned int  line,
+                           const unsigned char face_orientation) const
+        {
+          Assert(false, ExcNotImplemented());
+          (void)face;
+          (void)line;
+          (void)face_orientation;
+
+          return 0;
+        }
+
+        /**
+         * Map face vertex number to cell vertex number.
+         */
+        virtual unsigned int
+        face_to_cell_vertices(const unsigned int  face,
+                              const unsigned int  vertex,
+                              const unsigned char face_orientation) const
+        {
+          Assert(false, ExcNotImplemented());
+          (void)face;
+          (void)vertex;
+          (void)face_orientation;
+
+          return 0;
+        }
       };
 
 
@@ -267,6 +299,33 @@ namespace ReferenceCell
         {
           return GeometryInfo<dim>::faces_per_cell;
         }
+
+        unsigned int
+        face_to_cell_lines(const unsigned int  face,
+                           const unsigned int  line,
+                           const unsigned char face_orientation) const override
+        {
+          return GeometryInfo<dim>::face_to_cell_lines(
+            face,
+            line,
+            get_bit(face_orientation, 0),
+            get_bit(face_orientation, 2),
+            get_bit(face_orientation, 1));
+        }
+
+        unsigned int
+        face_to_cell_vertices(
+          const unsigned int  face,
+          const unsigned int  vertex,
+          const unsigned char face_orientation) const override
+        {
+          return GeometryInfo<dim>::face_to_cell_vertices(
+            face,
+            vertex,
+            get_bit(face_orientation, 0),
+            get_bit(face_orientation, 2),
+            get_bit(face_orientation, 1));
+        }
       };
 
 
@@ -359,6 +418,32 @@ namespace ReferenceCell
 
           return ReferenceCell::Type::Line;
         }
+
+        unsigned int
+        face_to_cell_lines(const unsigned int  face,
+                           const unsigned int  line,
+                           const unsigned char face_orientation) const override
+        {
+          AssertIndexRange(face, n_faces());
+          AssertDimension(line, 0);
+
+          (void)line;
+          (void)face_orientation;
+
+          return face;
+        }
+
+        unsigned int
+        face_to_cell_vertices(
+          const unsigned int  face,
+          const unsigned int  vertex,
+          const unsigned char face_orientation) const override
+        {
+          static const std::array<std::array<unsigned int, 2>, 3> table = {
+            {{0, 1}, {1, 2}, {2, 0}}};
+
+          return table[face][face_orientation ? vertex : (1 - vertex)];
+        }
       };
 
 
@@ -495,6 +580,33 @@ namespace ReferenceCell
 
           return ReferenceCell::Type::Tri;
         }
+
+        unsigned int
+        face_to_cell_lines(const unsigned int  face,
+                           const unsigned int  line,
+                           const unsigned char face_orientation) const override
+        {
+          AssertIndexRange(face, n_faces());
+
+          const static std::array<std::array<unsigned int, 3>, 4> table = {
+            {{0, 1, 2}, {0, 3, 4}, {2, 5, 3}, {1, 4, 5}}};
+
+          return table[face][standard_to_real_face_line(
+            line, face, face_orientation)];
+        }
+
+        unsigned int
+        face_to_cell_vertices(
+          const unsigned int  face,
+          const unsigned int  vertex,
+          const unsigned char face_orientation) const override
+        {
+          static const std::array<std::array<unsigned int, 3>, 4> table = {
+            {{0, 1, 2}, {1, 0, 3}, {0, 2, 3}, {2, 1, 3}}};
+
+          return table[face][standard_to_real_face_vertex(
+            vertex, face, face_orientation)];
+        }
       };
 
 
index f95944c77af5fc01984f6091bdb480adf315b17f..713f0637b902104563cac8e83134adbae9676eae 100644 (file)
@@ -544,8 +544,11 @@ FiniteElement<dim, spacedim>::face_to_cell_index(const unsigned int face_index,
                                                  const bool face_flip,
                                                  const bool face_rotation) const
 {
+  const auto &refence_cell =
+    ReferenceCell::internal::Info::get_cell(this->reference_cell_type());
+
   AssertIndexRange(face_index, this->n_dofs_per_face());
-  AssertIndexRange(face, GeometryInfo<dim>::faces_per_cell);
+  AssertIndexRange(face, refence_cell.n_faces());
 
   // TODO: we could presumably solve the 3d case below using the
   // adjust_quad_dof_index_for_face_orientation_table field. for the
@@ -582,8 +585,11 @@ FiniteElement<dim, spacedim>::face_to_cell_index(const unsigned int face_index,
 
       // then get the number of this vertex on the cell and translate
       // this to a DoF number on the cell
-      return (GeometryInfo<dim>::face_to_cell_vertices(
-                face, face_vertex, face_orientation, face_flip, face_rotation) *
+      return (refence_cell.face_to_cell_vertices(face,
+                                                 face_vertex,
+                                                 face_orientation +
+                                                   2 * face_rotation +
+                                                   4 * face_flip) *
                 this->n_dofs_per_vertex() +
               dof_index_on_vertex);
     }
@@ -598,8 +604,11 @@ FiniteElement<dim, spacedim>::face_to_cell_index(const unsigned int face_index,
       const unsigned int dof_index_on_line = index % this->n_dofs_per_line();
 
       return (this->get_first_line_index() +
-              GeometryInfo<dim>::face_to_cell_lines(
-                face, face_line, face_orientation, face_flip, face_rotation) *
+              refence_cell.face_to_cell_lines(face,
+                                              face_line,
+                                              face_orientation +
+                                                2 * face_rotation +
+                                                4 * face_flip) *
                 this->n_dofs_per_line() +
               dof_index_on_line);
     }

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.