]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Introduced FECollection::find_face_dominating_fe_in_subset(). 7284/head
authorMarc Fehling <marc.fehling@gmx.net>
Fri, 5 Oct 2018 17:47:52 +0000 (11:47 -0600)
committerMarc Fehling <marc.fehling@gmx.net>
Wed, 10 Oct 2018 00:26:38 +0000 (18:26 -0600)
Introduced get_active_fe_indices() in DoFAccessor namespace.

doc/news/changes/minor/20181005MarcFehling [new file with mode: 0644]
doc/news/changes/minor/20181008MarcFehling [new file with mode: 0644]
include/deal.II/dofs/dof_accessor.h
include/deal.II/dofs/dof_accessor.templates.h
include/deal.II/hp/fe_collection.h
source/dofs/dof_handler_policy.cc
source/fe/fe_enriched.cc
source/hp/fe_collection.cc

diff --git a/doc/news/changes/minor/20181005MarcFehling b/doc/news/changes/minor/20181005MarcFehling
new file mode 100644 (file)
index 0000000..4e24cf2
--- /dev/null
@@ -0,0 +1,5 @@
+New: Member function hp::FECollection::find_face_dominating_fe_in_subset()
+returns the index of the most dominating finite element out of a given
+set of indices.
+<br>
+(Marc Fehling, 2018/10/05)
diff --git a/doc/news/changes/minor/20181008MarcFehling b/doc/news/changes/minor/20181008MarcFehling
new file mode 100644 (file)
index 0000000..879a836
--- /dev/null
@@ -0,0 +1,5 @@
+New: Functions DoFAccessor::get_active_fe_indices() and
+internal::DoFAccessorImplementation::get_active_vertex_fe_indices()
+return all active finite element indices on the corresponding object.
+<br>
+(Marc Fehling, 2018/10/08)
index 9bd462f582ee1c61d5196b738328ebce9e08c849..4bd746724a0ed63dabcdfc34301ba48d3ee8fb8b 100644 (file)
@@ -555,6 +555,15 @@ public:
   unsigned int
   nth_active_fe_index(const unsigned int n) const;
 
+  /**
+   * Returns all active fe indices on this object.
+   *
+   * The size of the returned set equals the number of finite elements that
+   * are active on this object.
+   */
+  std::set<unsigned int>
+  get_active_fe_indices() const;
+
   /**
    * Return true if the finite element with given index is active on the
    * present object. For non-hp DoF accessors, this is of course the case only
index 55442f112186b0e4e3327c4b416f68b0726f8b9e..a7c898c9537f134d5962a7a2c2b79a231721689d 100644 (file)
@@ -1323,6 +1323,29 @@ namespace internal
 
 
 
+      /**
+       * Returns all active fe indices on a given vertex.
+       *
+       * The size of the returned set equals the number of finite elements that
+       * are active on this vertex.
+       */
+      template <int dim, int spacedim>
+      static std::set<unsigned int>
+      get_active_vertex_fe_indices(
+        const dealii::hp::DoFHandler<dim, spacedim> &dof_handler,
+        const unsigned int                           vertex_index)
+      {
+        std::set<unsigned int> active_fe_indices;
+        for (unsigned int i = 0;
+             i < n_active_vertex_fe_indices(dof_handler, vertex_index);
+             ++i)
+          active_fe_indices.insert(
+            nth_active_vertex_fe_index(dof_handler, vertex_index, i));
+        return active_fe_indices;
+      }
+
+
+
       /**
        * Return whether a particular finite element index is active on the
        * specified vertex.
@@ -1620,6 +1643,19 @@ DoFAccessor<dim, DoFHandlerType, level_dof_access>::nth_active_fe_index(
 
 
 
+template <int dim, typename DoFHandlerType, bool level_dof_access>
+inline std::set<unsigned int>
+DoFAccessor<dim, DoFHandlerType, level_dof_access>::get_active_fe_indices()
+  const
+{
+  std::set<unsigned int> active_fe_indices;
+  for (unsigned int i = 0; i < n_active_fe_indices(); ++i)
+    active_fe_indices.insert(nth_active_fe_index(i));
+  return active_fe_indices;
+}
+
+
+
 template <int dim, typename DoFHandlerType, bool level_dof_access>
 inline bool
 DoFAccessor<dim, DoFHandlerType, level_dof_access>::fe_index_is_active(
index 7284bdd119e85b5db9d81f08c33b3e840a204f11..10e09287e7326bfb458651488cad9b1faf2bb30f 100644 (file)
@@ -301,6 +301,33 @@ namespace hp
     DEAL_II_DEPRECATED unsigned int
     find_least_face_dominating_fe(const std::set<unsigned int> &fes) const;
 
+    /**
+     * Try to find a most face dominating finite element inside the subset of
+     * fe_indices @p fes as part of this FECollection. For example, if an
+     * FECollection consists of `{FE_Q(1),FE_Q(2),FE_Q(3),FE_Q(4)}` elements
+     * and we are looking for the most face dominating finite element among the
+     * last two elements of this collection (i.e., @p fes is `{2,3}`), then the
+     * answer is FE_Q(3) and therefore this function will return its index in
+     * the FECollection, namely `2`.
+     *
+     * This function differs from find_least_face_dominating_fe() in such a way
+     * that it looks for the most dominating finite element within the given
+     * subset @p fes, instead of finding a finite element in the whole
+     * FECollection that dominates all elements of the subset @p fes.
+     *
+     * For the purpose of this function by domination we consider either
+     * FiniteElementDomination::Domination::this_element_dominates or
+     * FiniteElementDomination::Domination::either_element_can_dominate;
+     * therefore the element can dominate itself. Thus, if an FECollection
+     * contains `{FE_Q(1),FE_Q(2),FE_Q(3),FE_Q(4)}` and @p fes only has
+     * a single element `{3}`, then the function returns 3.
+     *
+     * If the function is not able to find a finite element, the function
+     * returns numbers::invalid_unsigned_int.
+     */
+    unsigned int
+    find_face_dominating_fe_in_subset(const std::set<unsigned int> &fes) const;
+
     /**
      * Return a component mask with as many elements as this object has vector
      * components and of which exactly the one component is true that
index b8c4ed6a79e949b07fa16ffe878cd60476c49bbb..7297b431634460d384b933a6c10febf4fccc9bad 100644 (file)
@@ -221,65 +221,6 @@ namespace internal
                 }
             }
         }
-
-
-
-        /**
-         * For an object, such as a line or a quad iterator, determine
-         * the fe_index of the most dominating finite element that
-         * lives on this object.
-         *
-         * Return numbers::invalid_unsigned_int if we couldn't find one.
-         */
-        template <int dim, int spacedim, typename iterator>
-        unsigned int
-        get_most_dominating_fe_index(const iterator &object)
-        {
-          unsigned int dominating_fe_index = 0;
-          for (; dominating_fe_index < object->n_active_fe_indices();
-               ++dominating_fe_index)
-            {
-              const FiniteElement<dim, spacedim> &this_fe = object->get_fe(
-                object->nth_active_fe_index(dominating_fe_index));
-
-              FiniteElementDomination::Domination domination =
-                FiniteElementDomination::either_element_can_dominate;
-              for (unsigned int other_fe_index = 0;
-                   other_fe_index < object->n_active_fe_indices();
-                   ++other_fe_index)
-                if (other_fe_index != dominating_fe_index)
-                  {
-                    const FiniteElement<dim, spacedim> &that_fe =
-                      object->get_fe(
-                        object->nth_active_fe_index(other_fe_index));
-
-                    domination =
-                      domination & this_fe.compare_for_face_domination(that_fe);
-                  }
-
-              // see if this element is able to dominate all the other
-              // ones, and if so take it
-              if ((domination ==
-                   FiniteElementDomination::this_element_dominates) ||
-                  (domination ==
-                   FiniteElementDomination::either_element_can_dominate) ||
-                  (domination == FiniteElementDomination::no_requirements))
-                break;
-            }
-
-          // check that we have found one such fe
-          if (dominating_fe_index != object->n_active_fe_indices())
-            {
-              // return the finite element index used on it. note that
-              // only a single fe can be active on such subfaces
-              return object->nth_active_fe_index(dominating_fe_index);
-            }
-          else
-            {
-              // if we couldn't find the most dominating object
-              return numbers::invalid_unsigned_int;
-            }
-        }
       } // namespace
 
 
@@ -1020,10 +961,14 @@ namespace internal
                   // is not what we describe in the paper!.
                   if ((unique_sets_of_dofs == 2) && (dim == 2))
                     {
+                      const std::set<unsigned int> fe_indices =
+                        line->get_active_fe_indices();
+
                       // find out which is the most dominating finite element of
                       // the ones that are used on this line
                       const unsigned int most_dominating_fe_index =
-                        get_most_dominating_fe_index<dim, spacedim>(line);
+                        dof_handler.get_fe_collection()
+                          .find_face_dominating_fe_in_subset(fe_indices);
 
                       // if we found the most dominating element, then use this
                       // to eliminate some of the degrees of freedom by
@@ -1033,19 +978,12 @@ namespace internal
                       if (most_dominating_fe_index !=
                           numbers::invalid_unsigned_int)
                         {
-                          const unsigned int n_active_fe_indices =
-                            line->n_active_fe_indices();
-
                           // loop over the indices of all the finite elements
                           // that are not dominating, and identify their dofs to
                           // the most dominating one
-                          for (unsigned int f = 0; f < n_active_fe_indices; ++f)
-                            if (line->nth_active_fe_index(f) !=
-                                most_dominating_fe_index)
+                          for (const auto &other_fe_index : fe_indices)
+                            if (other_fe_index != most_dominating_fe_index)
                               {
-                                const unsigned int other_fe_index =
-                                  line->nth_active_fe_index(f);
-
                                 ensure_existence_of_dof_identities<1>(
                                   dof_handler.get_fe(most_dominating_fe_index),
                                   dof_handler.get_fe(other_fe_index),
@@ -1182,10 +1120,14 @@ namespace internal
                     quad = cell->quad(q);
                   quad->set_user_flag();
 
+                  const std::set<unsigned int> fe_indices =
+                    quad->get_active_fe_indices();
+
                   // find out which is the most dominating finite
                   // element of the ones that are used on this quad
                   const unsigned int most_dominating_fe_index =
-                    get_most_dominating_fe_index<dim, spacedim>(quad);
+                    dof_handler.get_fe_collection()
+                      .find_face_dominating_fe_in_subset(fe_indices);
 
                   // if we found the most dominating element, then use
                   // this to eliminate some of the degrees of freedom
@@ -1195,20 +1137,13 @@ namespace internal
                   // along this face/edge
                   if (most_dominating_fe_index != numbers::invalid_unsigned_int)
                     {
-                      const unsigned int n_active_fe_indices =
-                        quad->n_active_fe_indices();
-
                       // loop over the indices of all the finite
                       // elements that are not dominating, and
                       // identify their dofs to the most dominating
                       // one
-                      for (unsigned int f = 0; f < n_active_fe_indices; ++f)
-                        if (quad->nth_active_fe_index(f) !=
-                            most_dominating_fe_index)
+                      for (const auto &other_fe_index : fe_indices)
+                        if (other_fe_index != most_dominating_fe_index)
                           {
-                            const unsigned int other_fe_index =
-                              quad->nth_active_fe_index(f);
-
                             ensure_existence_of_dof_identities<2>(
                               dof_handler.get_fe(most_dominating_fe_index),
                               dof_handler.get_fe(other_fe_index),
@@ -1739,10 +1674,14 @@ namespace internal
                   // is not what we describe in the paper!.
                   if ((unique_sets_of_dofs == 2) && (dim == 2))
                     {
+                      const std::set<unsigned int> fe_indices =
+                        line->get_active_fe_indices();
+
                       // find out which is the most dominating finite element of
                       // the ones that are used on this line
                       const unsigned int most_dominating_fe_index =
-                        get_most_dominating_fe_index<dim, spacedim>(line);
+                        dof_handler.get_fe_collection()
+                          .find_face_dominating_fe_in_subset(fe_indices);
 
                       // if we found the most dominating element, then use this
                       // to eliminate some of the degrees of freedom by
@@ -1752,19 +1691,12 @@ namespace internal
                       if (most_dominating_fe_index !=
                           numbers::invalid_unsigned_int)
                         {
-                          const unsigned int n_active_fe_indices =
-                            line->n_active_fe_indices();
-
                           // loop over the indices of all the finite elements
                           // that are not dominating, and identify their dofs to
                           // the most dominating one
-                          for (unsigned int f = 0; f < n_active_fe_indices; ++f)
-                            if (line->nth_active_fe_index(f) !=
-                                most_dominating_fe_index)
+                          for (const auto &other_fe_index : fe_indices)
+                            if (other_fe_index != most_dominating_fe_index)
                               {
-                                const unsigned int other_fe_index =
-                                  line->nth_active_fe_index(f);
-
                                 ensure_existence_of_dof_identities<1>(
                                   dof_handler.get_fe(most_dominating_fe_index),
                                   dof_handler.get_fe(other_fe_index),
@@ -1884,10 +1816,14 @@ namespace internal
                     quad = cell->quad(q);
                   quad->clear_user_flag();
 
+                  const std::set<unsigned int> fe_indices =
+                    quad->get_active_fe_indices();
+
                   // find out which is the most dominating finite
                   // element of the ones that are used on this quad
                   const unsigned int most_dominating_fe_index =
-                    get_most_dominating_fe_index<dim, spacedim>(quad);
+                    dof_handler.get_fe_collection()
+                      .find_face_dominating_fe_in_subset(fe_indices);
 
                   // if we found the most dominating element, then use
                   // this to eliminate some of the degrees of freedom
@@ -1897,20 +1833,13 @@ namespace internal
                   // along this face/edge
                   if (most_dominating_fe_index != numbers::invalid_unsigned_int)
                     {
-                      const unsigned int n_active_fe_indices =
-                        quad->n_active_fe_indices();
-
                       // loop over the indices of all the finite
                       // elements that are not dominating, and
                       // identify their dofs to the most dominating
                       // one
-                      for (unsigned int f = 0; f < n_active_fe_indices; ++f)
-                        if (quad->nth_active_fe_index(f) !=
-                            most_dominating_fe_index)
+                      for (const auto &other_fe_index : fe_indices)
+                        if (other_fe_index != most_dominating_fe_index)
                           {
-                            const unsigned int other_fe_index =
-                              quad->nth_active_fe_index(f);
-
                             ensure_existence_of_dof_identities<2>(
                               dof_handler.get_fe(most_dominating_fe_index),
                               dof_handler.get_fe(other_fe_index),
index 6eb0b2069c87965e2276febfd242943c9fd4c3b3..e6aaa777056f10a93b31854e40972437c378a58f 100644 (file)
@@ -1258,14 +1258,14 @@ namespace ColorEnriched
        * Each time we build constraints at the
        * interface between two different FE_Enriched, we look for the least
        * dominating FE via
-       * hp::FECollection< dim, spacedim
-       * >::find_least_face_dominating_fe_in_collection(). If we don't take
-       * further actions, we may find a dominating FE that is too restrictive,
-       * i.e. enriched FE consisting of only FE_Nothing. New elements needs to
-       * be added to FECollection object to help find the correct enriched FE
-       * underlying the spaces in the adjacent cells. This is done by creating
-       * an appropriate set in fe_sets and a call to the function
-       * make_fe_collection_from_colored_enrichments at a later stage.
+       * hp::FECollection<dim,
+       * spacedim>::find_least_face_dominating_fe_in_collection(). If we don't
+       * take further actions, we may find a dominating FE that is too
+       * restrictive, i.e. enriched FE consisting of only FE_Nothing. New
+       * elements needs to be added to FECollection object to help find the
+       * correct enriched FE underlying the spaces in the adjacent cells. This
+       * is done by creating an appropriate set in fe_sets and a call to the
+       * function make_fe_collection_from_colored_enrichments at a later stage.
        *
        * Consider a domain with three predicates and hence with three different
        * enrichment functions. Let the enriched finite element of a cell with
index e01b11a101d96b4194bddc7457f201b06eb35ba0..4677a30efbee43d75eaaee2cb25c4b14d9230d08 100644 (file)
@@ -104,6 +104,52 @@ namespace hp
 
 
 
+  template <int dim, int spacedim>
+  unsigned int
+  FECollection<dim, spacedim>::find_face_dominating_fe_in_subset(
+    const std::set<unsigned int> &fes) const
+  {
+    for (auto it = fes.cbegin(); it != fes.cend(); ++it)
+      AssertIndexRange(*it, finite_elements.size());
+
+    // If the set of elements to be dominated contains only a single element X,
+    // then by definition the dominating set contains this single element
+    // (because each element can dominate itself).
+    // There may also be others, in which case we'll check if any of these
+    // elements is able to dominate all others. If one was found, we stop
+    // looking further and return the dominating element.
+    if (fes.size() == 1)
+      return *fes.begin();
+
+    // loop over all finite elements given in the subset
+    // and check which one dominates the whole subset
+    for (const auto &current_fe : fes)
+      {
+        FiniteElementDomination::Domination domination =
+          FiniteElementDomination::either_element_can_dominate;
+
+        for (const auto &other_fe : fes)
+          if (other_fe != current_fe)
+            domination =
+              domination &
+              finite_elements[current_fe]->compare_for_face_domination(
+                *finite_elements[other_fe]);
+
+        // see if this element is able to dominate all the other
+        // ones, and if so take it
+        if ((domination == FiniteElementDomination::this_element_dominates) ||
+            (domination ==
+             FiniteElementDomination::either_element_can_dominate) ||
+            (domination == FiniteElementDomination::no_requirements))
+          return current_fe;
+      }
+
+    // if we couldn't find the most dominating object
+    return numbers::invalid_unsigned_int;
+  }
+
+
+
   template <int dim, int spacedim>
   FECollection<dim, spacedim>::FECollection(
     const FiniteElement<dim, spacedim> &fe)

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.