]> https://gitweb.dealii.org/ - dealii.git/commitdiff
FE_PolyTensor: make a function use the combined orientation. 18533/head
authorDavid Wells <drwells@email.unc.edu>
Sun, 1 Jun 2025 20:13:09 +0000 (16:13 -0400)
committerDavid Wells <drwells@email.unc.edu>
Sun, 1 Jun 2025 20:23:37 +0000 (16:23 -0400)
include/deal.II/fe/fe_poly_tensor.h
source/fe/fe_poly_tensor.cc

index 1a66c682c0a393449aadd7eb1883430a498e8348..8ad405734f0ff0a3268bb37ff907e9611b141790 100644 (file)
@@ -225,23 +225,19 @@ protected:
 
   /**
    * For faces with non-standard face_orientation in 3d, the dofs on faces
-   * (quads) have to be permuted in order to be combined with the correct
-   * shape functions and additionally can change the sign. Given a local
-   * dof @p index on a quad, return the
-   * sign of the permuted shape function, if the face has non-standard
-   * face_orientation, face_flip or face_rotation. In 2d and 1d there is no need
-   * for permutation and consequently it does nothing in this case.
+   * (quads) have to be permuted in order to be combined with the correct shape
+   * functions and additionally can change the sign. Given a local dof @p index
+   * on a quad, return the sign of the permuted shape function.
    *
    * The permutation itself is returned by
    * adjust_quad_dof_index_for_face_orientation implemented in the interface
    * class FiniteElement<dim>.
    */
   bool
-  adjust_quad_dof_sign_for_face_orientation(const unsigned int index,
-                                            const unsigned int face_no,
-                                            const bool         face_orientation,
-                                            const bool         face_flip,
-                                            const bool face_rotation) const;
+  adjust_quad_dof_sign_for_face_orientation(
+    const unsigned int                 index,
+    const unsigned int                 face_no,
+    const types::geometric_orientation combined_orientation) const;
 
   /**
    * For faces with non-standard face_orientation in 3d, the dofs on faces
index c3426e46c776030fc5dfc5d12f8d6f6b0b1aa9a0..bd2bc5d866513c9e55078a2efb471907d688f54d 100644 (file)
@@ -341,11 +341,9 @@ FE_PolyTensor<dim, spacedim>::single_mapping_kind() const
 template <int dim, int spacedim>
 bool
 FE_PolyTensor<dim, spacedim>::adjust_quad_dof_sign_for_face_orientation(
-  const unsigned int index,
-  const unsigned int face,
-  const bool         face_orientation,
-  const bool         face_flip,
-  const bool         face_rotation) const
+  const unsigned int                 index,
+  const unsigned int                 face,
+  const types::geometric_orientation combined_orientation) const
 {
   // do nothing in 1d and 2d
   if (dim < 3)
@@ -366,11 +364,8 @@ FE_PolyTensor<dim, spacedim>::adjust_quad_dof_sign_for_face_orientation(
          ExcInternalError());
 
   return adjust_quad_dof_sign_for_face_orientation_table
-    [this->n_unique_2d_subobjects() == 1 ? 0 : face](
-      index,
-      internal::combined_face_orientation(face_orientation,
-                                          face_rotation,
-                                          face_flip));
+    [this->n_unique_2d_subobjects() == 1 ? 0 : face](index,
+                                                     combined_orientation);
 }
 
 
@@ -630,9 +625,7 @@ FE_PolyTensor<dim, spacedim>::fill_fe_values(
           if (adjust_quad_dof_sign_for_face_orientation(
                 local_quad_dof_index,
                 face_index_from_dof_index,
-                cell->face_orientation(face_index_from_dof_index),
-                cell->face_flip(face_index_from_dof_index),
-                cell->face_rotation(face_index_from_dof_index)))
+                cell->combined_face_orientation(face_index_from_dof_index)))
             dof_sign = -1.0;
         }
 
@@ -1267,9 +1260,7 @@ FE_PolyTensor<dim, spacedim>::fill_fe_face_values(
           if (adjust_quad_dof_sign_for_face_orientation(
                 local_quad_dof_index,
                 face_index_from_dof_index,
-                cell->face_orientation(face_index_from_dof_index),
-                cell->face_flip(face_index_from_dof_index),
-                cell->face_rotation(face_index_from_dof_index)))
+                cell->combined_face_orientation(face_index_from_dof_index)))
             dof_sign = -1.0;
         }
 
@@ -1951,9 +1942,7 @@ FE_PolyTensor<dim, spacedim>::fill_fe_subface_values(
           if (adjust_quad_dof_sign_for_face_orientation(
                 local_quad_dof_index,
                 face_index_from_dof_index,
-                cell->face_orientation(face_index_from_dof_index),
-                cell->face_flip(face_index_from_dof_index),
-                cell->face_rotation(face_index_from_dof_index)))
+                cell->combined_face_orientation(face_index_from_dof_index)))
             dof_sign = -1.0;
         }
 

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.