*/
bool hp_constraints_are_implemented () const;
+ /**
+ * Try to find a least dominant finite element inside this FECollection
+ * which dominates other finite elements provided as fe_indices in @p fes .
+ * For example, if FECollection consists of {Q1,Q2,Q3,Q4} and we are looking
+ * for the least dominant FE for Q3 and Q4 (@p fes is {2,3}), then the
+ * answer is Q3 and therefore this function will return its index in
+ * FECollection, namely 2.
+ *
+ * If we were not able to find a finite element, the function returns
+ * numbers::invalid_unsigned_int .
+ *
+ * Note that for the cases like when FECollection consists of
+ * {FE_Nothing x FE_Nothing, Q1xQ2, Q2xQ1} with @p fes = {1}, the function
+ * will not find the most dominating element as the default behavior of
+ * FE_Nothing is to return FiniteElementDomination::no_requirements when
+ * comparing for face domination. This, therefore, can't be considered as a
+ * dominating element in the sense defined in FiniteElementDomination .
+ *
+ * Finally, for the purpose of this function we consider that an element
+ * dominates itself. Thus if FECollection contains {Q1,Q2,Q4,Q3} and @p fes
+ * = {3}, the function returns 3.
+ */
+ unsigned int
+ find_least_face_dominating_fe (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
namespace hp
{
+ template <int dim, int spacedim>
+ unsigned int
+ FECollection<dim,spacedim>::find_least_face_dominating_fe (const std::set<unsigned int> &fes) const
+ {
+ // If the set of elements to be dominated contains only a single element X,
+ // then by definition the dominating set contains this single element X
+ // (because each element can dominate itself). There may also be others,
+ // say Y1...YN. Next you have to find one or more elements in the dominating
+ // set {X,Y1...YN} that is the weakest. Well, you can't find one that is
+ // weaker than X because if it were, it would not dominate X. In other words,
+ // X is guaranteed to be in the subset of {X,Y1...YN} of weakest dominating
+ // elements. Since we only guarantee that the function returns one of them,
+ // we may as well return X right away.
+ if (fes.size()==1)
+ return *fes.begin();
+
+ const hp::FECollection<dim,spacedim> &fe_collection = *this;
+ std::set<unsigned int> candidate_fes;
+
+ // first loop over all FEs and check which can dominate those given in FEs:
+ for (unsigned int cur_fe = 0; cur_fe < fe_collection.size(); cur_fe++)
+ {
+ FiniteElementDomination::Domination domination = FiniteElementDomination::no_requirements;
+ // check if cur_fe can dominate all FEs in @p fes:
+ for (std::set<unsigned int>::const_iterator it = fes.begin();
+ it!=fes.end(); ++it)
+ {
+ Assert (*it < fe_collection.size(),
+ ExcIndexRangeType<unsigned int> (*it, 0, fe_collection.size()));
+ domination = domination &
+ fe_collection[cur_fe].compare_for_face_domination
+ (fe_collection[*it]);
+ }
+
+ // if we found dominating element, keep them in a set.
+ if (domination == FiniteElementDomination::this_element_dominates)
+ candidate_fes.insert(cur_fe);
+ }
+
+ // among the ones we found, pick one that is dominated by all others and
+ // thus should represent the largest FE space.
+ if (candidate_fes.size() == 1)
+ {
+ return *candidate_fes.begin();
+ }
+ else
+ for (std::set<unsigned int>::const_iterator it = candidate_fes.begin(); it!=candidate_fes.end(); ++it)
+ {
+ FiniteElementDomination::Domination domination = FiniteElementDomination::no_requirements;
+ for (std::set<unsigned int>::const_iterator ito = candidate_fes.begin(); ito!=candidate_fes.end(); ++ito)
+ if (it != ito)
+ {
+ domination = domination &
+ fe_collection[*it].compare_for_face_domination(fe_collection[*ito]);
+ }
+
+ if (domination == FiniteElementDomination::other_element_dominates ||
+ domination == FiniteElementDomination::either_element_can_dominate /*covers cases like candidate_fes={Q1,Q1}*/)
+ return *it;
+ }
+ // We couldn't find the FE, return invalid_unsigned_int :
+ return numbers::invalid_unsigned_int;
+ }
+
template <int dim, int spacedim>
FECollection<dim,spacedim>::FECollection ()
{}
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2005 - 2015 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+
+// test the results of FECollection::find_least_face_dominating_fe(), namely for:
+// {Q1,Q2,Q3,Q4} with {2,3} => Q3 2
+// {Q1xQ1, Q2xQ2, Q3xQ4, Q4xQ3} with {2,3} => Q2xQ2 1
+// {Q1xQ1, Q3xQ4, Q4xQ3} with {1,2} => Q1xQ1 0
+// {0x0, 0x0, Q1x0, 0xQ1} with {2,3} => none invalid_unsigned_int
+// {0x0, 0x0, Q1x0, 0xQ1} with {2,3} => 0x0 0 (with dominating FE_Nothing)
+// {Q1xQ1,Q1xQ1,Q2xQ1,Q1,Q2} with {2,3} => Q1 0
+// {Q4xQ4, Q5xQ5, Q3xQ4, Q4xQ3} with {2,3} => none invalid_unsigned_int
+// {Q1,Q2,Q4,Q3} with {3} => Q3 3
+
+
+#include "../tests.h"
+#include <deal.II/base/logstream.h>
+#include <deal.II/hp/fe_collection.h>
+#include <deal.II/fe/fe_system.h>
+#include <deal.II/fe/fe_q.h>
+#include <deal.II/fe/fe_nothing.h>
+
+#include <fstream>
+
+
+template <int dim>
+void test ()
+{
+ std::set<unsigned int> fes;
+ fes.insert(2);
+ fes.insert(3);
+
+ // {Q1,Q2,Q3,Q4}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FE_Q<dim>(1));
+ fe_collection.push_back (FE_Q<dim>(2));
+ fe_collection.push_back (FE_Q<dim>(3));
+ fe_collection.push_back (FE_Q<dim>(4));
+ deallog << fe_collection.find_least_face_dominating_fe(fes) << std::endl;
+ }
+
+ // {Q1xQ1, Q2xQ2, Q3xQ4, Q4xQ3}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Q<dim>(1),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(2),1,
+ FE_Q<dim>(2),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(3),1,
+ FE_Q<dim>(4),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(4),1,
+ FE_Q<dim>(3),1));
+ deallog << fe_collection.find_least_face_dominating_fe(fes) << std::endl;
+ }
+
+ // {Q1xQ1, Q3xQ4, Q4xQ3}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Q<dim>(1),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(3),1,
+ FE_Q<dim>(4),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(4),1,
+ FE_Q<dim>(3),1));
+ std::set<unsigned int> fes;
+ fes.insert(1);
+ fes.insert(2);
+ deallog << fe_collection.find_least_face_dominating_fe(fes) << std::endl;
+ }
+
+ // {0x0, 0x0, Q1x0, 0xQ1}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FESystem<dim>(FE_Nothing<dim>(),1,
+ FE_Nothing<dim>(),1));
+ fe_collection.push_back (FESystem<dim>(FE_Nothing<dim>(),1,
+ FE_Nothing<dim>(),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Nothing<dim>(),1));
+ fe_collection.push_back (FESystem<dim>(FE_Nothing<dim>(),1,
+ FE_Q<dim>(1),1));
+ const unsigned int ind = fe_collection.find_least_face_dominating_fe(fes);
+ if (ind == numbers::invalid_unsigned_int)
+ deallog << "numbers::invalid_unsigned_int" << std::endl;
+ else
+ deallog << ind << std::endl;
+ }
+
+ // dominating FE_Nothing
+ // {0x0, 0x0, Q1x0, 0xQ1}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FESystem<dim>(FE_Nothing<dim>(1,true),1,
+ FE_Nothing<dim>(1,true),1));
+ fe_collection.push_back (FESystem<dim>(FE_Nothing<dim>(1,true),1,
+ FE_Nothing<dim>(1,true),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Nothing<dim>(1,true),1));
+ fe_collection.push_back (FESystem<dim>(FE_Nothing<dim>(1,true),1,
+ FE_Q<dim>(1),1));
+ deallog << fe_collection.find_least_face_dominating_fe(fes) << std::endl;
+ }
+
+
+ // {Q1xQ1,Q1xQ1,Q2xQ1,Q1,Q2}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Q<dim>(1),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Q<dim>(1),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(2),1,
+ FE_Q<dim>(1),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(1),1,
+ FE_Q<dim>(2),1));
+ deallog << fe_collection.find_least_face_dominating_fe(fes) << std::endl;
+ }
+
+ // {Q4xQ4, Q5xQ5, Q3xQ4, Q4xQ3
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(4),1,
+ FE_Q<dim>(4),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(5),1,
+ FE_Q<dim>(5),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(3),1,
+ FE_Q<dim>(4),1));
+ fe_collection.push_back (FESystem<dim>(FE_Q<dim>(4),1,
+ FE_Q<dim>(3),1));
+ const unsigned int ind = fe_collection.find_least_face_dominating_fe(fes);
+ if (ind == numbers::invalid_unsigned_int)
+ deallog << "numbers::invalid_unsigned_int" << std::endl;
+ else
+ deallog << ind << std::endl;
+ }
+
+ // {Q1,Q2,Q4,Q3}
+ {
+ hp::FECollection<dim> fe_collection;
+ fe_collection.push_back (FE_Q<dim>(1));
+ fe_collection.push_back (FE_Q<dim>(2));
+ fe_collection.push_back (FE_Q<dim>(4));
+ fe_collection.push_back (FE_Q<dim>(3));
+ std::set<unsigned int> fes;
+ fes.insert(3);
+ deallog << fe_collection.find_least_face_dominating_fe(fes) << std::endl;
+ }
+
+
+}
+
+int main ()
+{
+ std::ofstream logfile("output");
+ logfile.precision(2);
+
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ test<2> ();
+ test<3> ();
+
+ deallog << "OK" << std::endl;
+}