]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Add two convenience functions to FE_PolyTensor and improve based on suggestions
authorgrahambenharper <grahambenharper@gmail.com>
Tue, 9 Jul 2019 10:14:47 +0000 (04:14 -0600)
committergrahambenharper <grahambenharper@gmail.com>
Tue, 9 Jul 2019 10:21:01 +0000 (04:21 -0600)
include/deal.II/fe/fe_poly_tensor.h
source/fe/fe_poly_tensor.cc

index 89034a6f45a05fb5111d42437977d4f9c2dd9987..4bb34d4dc8a2f340b532000f3bc1e4e6361b03d0 100644 (file)
@@ -216,6 +216,18 @@ protected:
    */
   std::vector<MappingType> mapping_type;
 
+  /**
+   * Returns a boolean that is true when the finite element uses a single
+   * mapping and false when the finite element uses multiple mappings.
+   */
+  bool
+  single_mapping() const;
+
+  /**
+   * Returns MappingType @p i for the finite element.
+   */
+  MappingType
+  get_mapping_type(unsigned int i) const;
 
   /* NOTE: The following function has its definition inlined into the class
      declaration because we otherwise run into a compiler error with MS Visual
@@ -254,19 +266,17 @@ protected:
     // allocate memory
 
     const bool update_transformed_shape_values =
-      (std::find_if(this->mapping_type.begin(),
-                    this->mapping_type.end(),
-                    [](const MappingType t) { return t != mapping_none; }) !=
-       this->mapping_type.end());
+      std::any_of(this->mapping_type.begin(),
+                  this->mapping_type.end(),
+                  [](const MappingType t) { return t != mapping_none; });
 
     const bool update_transformed_shape_grads =
-      (std::find_if(this->mapping_type.begin(),
-                    this->mapping_type.end(),
-                    [](const MappingType t) {
-                      return (t == mapping_raviart_thomas ||
-                              t == mapping_piola || t == mapping_nedelec ||
-                              t == mapping_contravariant);
-                    }) != this->mapping_type.end());
+      std::any_of(this->mapping_type.begin(),
+                  this->mapping_type.end(),
+                  [](const MappingType t) {
+                    return (t == mapping_raviart_thomas || t == mapping_piola ||
+                            t == mapping_nedelec || t == mapping_contravariant);
+                  });
 
     const bool update_transformed_shape_hessian_tensors =
       update_transformed_shape_values;
index bc1d0844a30af4c20a9a72c49da169d352f6ea71..28bc8e5a52d60460757c1622049e4bedd8e81825 100644 (file)
@@ -158,7 +158,8 @@ namespace internal
         // TODO: This is probably only going to work for those elements for
         // which all dofs are face dofs
         for (unsigned int l = 0; l < GeometryInfo<dim>::lines_per_cell; ++l)
-          if (!(cell->line_orientation(l)) & mapping_type[0] == mapping_nedelec)
+          if (!(cell->line_orientation(l)) &&
+              mapping_type[0] == mapping_nedelec)
             face_sign[l] = -1.0;
       }
     } // namespace
@@ -193,6 +194,28 @@ FE_PolyTensor<PolynomialType, dim, spacedim>::FE_PolyTensor(
 
 
 
+template <class PolynomialType, int dim, int spacedim>
+bool
+FE_PolyTensor<PolynomialType, dim, spacedim>::single_mapping() const
+{
+  return mapping_type.size() == 1;
+}
+
+
+
+template <class PolynomialType, int dim, int spacedim>
+MappingType
+FE_PolyTensor<PolynomialType, dim, spacedim>::get_mapping_type(
+  unsigned int i) const
+{
+  if (single_mapping())
+    return mapping_type[0];
+  Assert(i < mapping_type.size(), ExcIndexRange(i, 0, mapping_type.size()));
+  return mapping_type[i];
+}
+
+
+
 template <class PolynomialType, int dim, int spacedim>
 double
 FE_PolyTensor<PolynomialType, dim, spacedim>::shape_value(
@@ -392,9 +415,7 @@ FE_PolyTensor<PolynomialType, dim, spacedim>::fill_fe_values(
 
   for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
     {
-      MappingType mapping_type =
-        (this->mapping_type.size() > 1 ? this->mapping_type[i] :
-                                         this->mapping_type[0]);
+      const MappingType mapping_type = get_mapping_type(i);
 
       const unsigned int first =
         output_data.shape_function_to_row_table[i * this->n_components() +
@@ -972,9 +993,7 @@ FE_PolyTensor<PolynomialType, dim, spacedim>::fill_fe_face_values(
 
   for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
     {
-      const MappingType mapping_type =
-        (this->mapping_type.size() > 1 ? this->mapping_type[i] :
-                                         this->mapping_type[0]);
+      const MappingType mapping_type = get_mapping_type(i);
 
       const unsigned int first =
         output_data.shape_function_to_row_table[i * this->n_components() +
@@ -1604,9 +1623,7 @@ FE_PolyTensor<PolynomialType, dim, spacedim>::fill_fe_subface_values(
 
   for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
     {
-      MappingType mapping_type =
-        (this->mapping_type.size() > 1 ? this->mapping_type[i] :
-                                         this->mapping_type[0]);
+      const MappingType mapping_type = get_mapping_type(i);
 
       const unsigned int first =
         output_data.shape_function_to_row_table[i * this->n_components() +
@@ -2170,9 +2187,7 @@ FE_PolyTensor<PolynomialType, dim, spacedim>::requires_update_flags(
 
   for (unsigned int i = 0; i < this->dofs_per_cell; ++i)
     {
-      MappingType mapping_type =
-        (this->mapping_type.size() > 1 ? this->mapping_type[i] :
-                                         this->mapping_type[0]);
+      const MappingType mapping_type = get_mapping_type(i);
 
       switch (mapping_type)
         {

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.