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
}
}
}
-
-
-
- /**
- * 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
// 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
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),
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
// 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),
// 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
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),
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
// 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),
+ 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 ¤t_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)