]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implement FiniteElement::shape_function_belongs_to().
authorWolfgang Bangerth <bangerth@colostate.edu>
Thu, 19 Dec 2024 01:29:40 +0000 (18:29 -0700)
committerWolfgang Bangerth <bangerth@colostate.edu>
Thu, 19 Dec 2024 01:29:40 +0000 (18:29 -0700)
include/deal.II/fe/fe.h

index 9ba23410461c77605f5970b9d0e082e216bb0af1..c727e2be94379ed05002e5da6f59fcbd37297c79 100644 (file)
@@ -1825,10 +1825,84 @@ public:
   unsigned int
   component_to_block_index(const unsigned int component) const;
 
+  /**
+   * For vector-valued finite elements, return whether the shape
+   * function indicated by the first argument is part of the scalar
+   * vector component indicated by the second argument. In practice,
+   * what this means is that the set of nonzero vector components of
+   * the indicated shape function includes the one vector component
+   * represented by the second argument; in other words, whether the
+   * shape function is (potentially) nonzero in the given vector
+   * component.
+   *
+   * A typical example would be in a program like step-22, where we
+   * might want to know whether a shape function of the combined,
+   * vector-valued element, corresponds to the pressure part of the
+   * solution.
+   */
+  bool
+  shape_function_belongs_to(const unsigned int                shape_function,
+                            const FEValuesExtractors::Scalar &component) const;
+
+  /**
+   * For vector-valued finite elements, return whether the shape
+   * function indicated by the first argument is part of the `dim`
+   * vector components indicated by the second argument. In practice,
+   * what this means is that the set of nonzero vector components of
+   * the indicated shape function overlaps the vector components
+   * represented by the second argument; in other words, whether the
+   * shape function is (potentially) nonzero any of the given vector
+   * components.
+   *
+   * A typical example would be in a program like step-22, where we
+   * might want to know whether a shape function of the combined,
+   * vector-valued element, corresponds to the velocity part of the
+   * solution.
+   */
+  bool
+  shape_function_belongs_to(const unsigned int                shape_function,
+                            const FEValuesExtractors::Vector &components) const;
+
+
+  /**
+   * For vector-valued finite elements, return whether the shape
+   * function indicated by the first argument is part of the `dim*(dim+1)/2`
+   * vector components indicated by the second argument. In practice,
+   * what this means is that the set of nonzero vector components of
+   * the indicated shape function overlaps the vector components
+   * represented by the second argument; in other words, whether the
+   * shape function is (potentially) nonzero any of the given
+   * components that represent a symmetric tensor.
+   *
+   * An example would be in a program like step-71, where we
+   * might want to know whether a shape function of the combined,
+   * vector-valued element, corresponds to the $H$ or $C$ part of the
+   * solution, where $H$ is a vector and $C$ is a symmetric tensor.
+   */
+  bool
+  shape_function_belongs_to(
+    const unsigned int                            shape_function,
+    const FEValuesExtractors::SymmetricTensor<2> &components) const;
+
+
+  /**
+   * For vector-valued finite elements, return whether the shape
+   * function indicated by the first argument is part of the `dim*dim`
+   * vector components indicated by the second argument. In practice,
+   * what this means is that the set of nonzero vector components of
+   * the indicated shape function overlaps the vector components
+   * represented by the second argument; in other words, whether the
+   * shape function is (potentially) nonzero any of the given
+   * components that represent a tensor.
+   */
+  bool
+  shape_function_belongs_to(
+    const unsigned int                   shape_function,
+    const FEValuesExtractors::Tensor<2> &components) const;
   /** @} */
 
   /**
-   * @name Component and block matrices
+   * @name Component and block masks
    * @{
    */
 
@@ -3309,6 +3383,124 @@ FiniteElement<dim, spacedim>::system_to_block_index(
 
 
 
+template <int dim, int spacedim>
+inline bool
+FiniteElement<dim, spacedim>::shape_function_belongs_to(
+  const unsigned int                shape_function,
+  const FEValuesExtractors::Scalar &component) const
+{
+  AssertIndexRange(shape_function, this->n_dofs_per_cell());
+  AssertIndexRange(component.component, this->n_components());
+
+  if (is_primitive(shape_function))
+    return (system_to_component_table[shape_function].first ==
+            component.component);
+  else
+    return nonzero_components[shape_function][component.component];
+}
+
+
+
+template <int dim, int spacedim>
+inline bool
+FiniteElement<dim, spacedim>::shape_function_belongs_to(
+  const unsigned int                shape_function,
+  const FEValuesExtractors::Vector &components) const
+{
+  AssertIndexRange(shape_function, this->n_dofs_per_cell());
+  AssertIndexRange(components.first_vector_component, this->n_components());
+  AssertIndexRange(components.first_vector_component + dim,
+                   this->n_components() + 1);
+
+  if (is_primitive(shape_function))
+    return ((system_to_component_table[shape_function].first >=
+             components.first_vector_component) &&
+            (system_to_component_table[shape_function].first <
+             components.first_vector_component + dim));
+  else
+    // Return whether there is overlap between the nonzero components
+    // of the current shape function and the selected components of
+    // the extractor:
+    {
+      for (unsigned int i = components.first_vector_component;
+           i < components.first_vector_component + dim;
+           ++i)
+        if (nonzero_components[shape_function][i])
+          return true;
+      return false;
+    }
+}
+
+
+
+template <int dim, int spacedim>
+inline bool
+FiniteElement<dim, spacedim>::shape_function_belongs_to(
+  const unsigned int                            shape_function,
+  const FEValuesExtractors::SymmetricTensor<2> &components) const
+{
+  AssertIndexRange(shape_function, this->n_dofs_per_cell());
+  AssertIndexRange(components.first_tensor_component, this->n_components());
+  AssertIndexRange(components.first_tensor_component + dim,
+                   this->n_components() + 1);
+
+  if (is_primitive(shape_function))
+    return ((system_to_component_table[shape_function].first >=
+             components.first_tensor_component) &&
+            (system_to_component_table[shape_function].first <
+             components.first_tensor_component +
+               SymmetricTensor<2, dim>::n_independent_components));
+  else
+    // Return whether there is overlap between the nonzero components
+    // of the current shape function and the selected components of
+    // the extractor:
+    {
+      for (unsigned int i = components.first_tensor_component;
+           i < components.first_tensor_component +
+                 SymmetricTensor<2, dim>::n_independent_components;
+           ++i)
+        if (nonzero_components[shape_function][i])
+          return true;
+      return false;
+    }
+}
+
+
+
+template <int dim, int spacedim>
+inline bool
+FiniteElement<dim, spacedim>::shape_function_belongs_to(
+  const unsigned int                   shape_function,
+  const FEValuesExtractors::Tensor<2> &components) const
+{
+  AssertIndexRange(shape_function, this->n_dofs_per_cell());
+  AssertIndexRange(components.first_tensor_component, this->n_components());
+  AssertIndexRange(components.first_tensor_component + dim,
+                   this->n_components() + 1);
+
+  if (is_primitive(shape_function))
+    return ((system_to_component_table[shape_function].first >=
+             components.first_tensor_component) &&
+            (system_to_component_table[shape_function].first <
+             components.first_tensor_component +
+               Tensor<2, dim>::n_independent_components));
+  else
+    // Return whether there is overlap between the nonzero components
+    // of the current shape function and the selected components of
+    // the extractor:
+    {
+      for (unsigned int i = components.first_tensor_component;
+           i < components.first_tensor_component +
+                 Tensor<2, dim>::n_independent_components;
+           ++i)
+        if (nonzero_components[shape_function][i])
+          return true;
+      return false;
+    }
+}
+
+
+
 template <int dim, int spacedim>
 inline bool
 FiniteElement<dim, spacedim>::restriction_is_additive(

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.