]> https://gitweb.dealii.org/ - dealii.git/commitdiff
FiniteElement: remove 2d orientation assumptions.
authorDavid Wells <drwells@email.unc.edu>
Sat, 13 Apr 2024 11:29:26 +0000 (07:29 -0400)
committerDavid Wells <drwells@email.unc.edu>
Sat, 13 Apr 2024 11:44:58 +0000 (07:44 -0400)
For mixed meshes, we need to support quadrilaterals whose faces may have either
valid line orientation.

include/deal.II/fe/fe.h
source/fe/fe.cc

index 4a5d0f83f99e23cb63698bd0139db9733b4be4ef..f55dbcf96bd10f5fca557e403b969b1116f8a3d4 100644 (file)
@@ -1529,8 +1529,9 @@ public:
    * combined_orientation of that line, return the local dof which accounts for
    * @p combined_orientation.
    *
-   * @note In both 1d and all-quadrilateral meshes in 2d all lines have the
-   * standard orientation.
+   * @note In both 1d and 2d all-quadrilateral meshes all lines have the
+   * standard orientation. However, since 2d meshes may contain both
+   * quadrilaterals and triangles, this assumption cannot be made in this class.
    */
   unsigned int
   adjust_line_dof_index_for_line_orientation(
@@ -2471,13 +2472,11 @@ protected:
   std::vector<Table<2, int>> adjust_quad_dof_index_for_face_orientation_table;
 
   /**
-   * For lines with non-standard line_orientation in 3d, the dofs on lines
+   * For lines with non-standard orientation in 2d or 3d, the dofs on lines
    * have to be permuted in order to be combined with the correct shape
    * functions. Given a local dof @p index on a line, return the shift in the
    * local index, if the line has non-standard line_orientation, i.e.
-   * <code>old_index + shift = new_index</code>. In 2d and 1d there is no need
-   * for permutation so the vector is empty. In 3d it has the size of
-   * #dofs_per_line.
+   * <code>old_index + shift = new_index</code>.
    *
    * The constructor of this class fills this table with zeros, i.e.,
    * no permutation at all. Derived finite element classes have to
index 6e3bf6edf629af2c1563fb9169c73f157e16205b..a114eff2f76ade4e35b5eea2494e2eef80293512 100644 (file)
@@ -569,15 +569,6 @@ FiniteElement<dim, spacedim>::face_to_cell_index(
   AssertIndexRange(face_index, this->n_dofs_per_face(face));
   AssertIndexRange(face, this->reference_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 2d case, we
-  // can't use adjust_line_dof_index_for_line_orientation_table since that array
-  // is not populated for elements with quadrilateral reference cells
-  // (presumably because we thought that there are no flipped edges in 2d, but
-  // these can happen in DoFTools::make_periodicity_constraints(), for example).
-  // so we would need to either fill this field, or rely on derived classes
-  // implementing this function, as we currently do
-
   // see the function's documentation for an explanation of this
   // assertion -- in essence, derived classes have to implement
   // an overloaded version of this function if we are to use any
@@ -591,7 +582,7 @@ FiniteElement<dim, spacedim>::face_to_cell_index(
              "an overloaded version but apparently hasn't done so. See "
              "the documentation of this function for more information."));
 
-  // we need to distinguish between DoFs on vertices, lines and in 3d quads.
+  // we need to distinguish between DoFs on vertices, lines and (in 3d) quads.
   // do so in a sequence of if-else statements
   if (face_index < this->get_first_face_line_index(face))
     // DoF is on a vertex
@@ -685,19 +676,11 @@ FiniteElement<dim, spacedim>::adjust_line_dof_index_for_line_orientation(
   const unsigned int  index,
   const unsigned char combined_orientation) const
 {
-  // We orient quads (and 1D meshes are always oriented) so always skip those
-  // cases
-  //
-  // TODO - we may want to change this in the future: see also the notes in
-  // face_to_cell_index()
   Assert(combined_orientation ==
              ReferenceCell::default_combined_face_orientation() ||
            combined_orientation ==
              ReferenceCell::reversed_combined_line_orientation(),
          ExcInternalError());
-  if (this->reference_cell() == ReferenceCells::Line ||
-      this->reference_cell() == ReferenceCells::Quadrilateral)
-    return index;
 
   AssertIndexRange(index, this->n_dofs_per_line());
   Assert(adjust_line_dof_index_for_line_orientation_table.size() ==

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.