31 template <
int dim,
int spacedim>
34 set_default_hierarchy();
39 template <
int dim,
int spacedim>
49 template <
int dim,
int spacedim>
55 ExcMessage(
"Need to pass at least one finite element."));
57 for (
unsigned int i = 0; i < fes.size(); ++i)
63 template <
int dim,
int spacedim>
71 Assert(this->size() == 0 ||
72 new_fe.
n_components() == this->operator[](0).n_components(),
73 ExcMessage(
"All elements inside a collection need to have the "
74 "same number of vector components!"));
81 template <
int dim,
int spacedim>
85 Assert(this->size() > 0, ExcNoFiniteElements());
96 if (!reference_cell_default_linear_mapping ||
97 reference_cell_default_linear_mapping->size() != this->size())
102 std::make_shared<MappingCollection<dim, spacedim>>();
104 for (
const auto &fe : *
this)
105 this_nc.reference_cell_default_linear_mapping->push_back(
107 .template get_default_linear_mapping<dim, spacedim>());
110 return *reference_cell_default_linear_mapping;
115 template <
int dim,
int spacedim>
116 std::set<unsigned int>
118 const std::set<unsigned int> &fes,
119 const unsigned int codim)
const
125 for (
const auto &fe : fes)
132 std::set<unsigned int> dominating_fes;
133 for (
unsigned int current_fe = 0; current_fe < this->size(); ++current_fe)
138 for (
const auto &other_fe : fes)
141 this->operator[](current_fe)
142 .compare_for_domination(this->operator[](other_fe), codim);
148 dominating_fes.insert(current_fe);
150 return dominating_fes;
155 template <
int dim,
int spacedim>
156 std::set<unsigned int>
158 const std::set<unsigned int> &fes,
159 const unsigned int codim)
const
165 for (
const auto &fe : fes)
172 std::set<unsigned int> dominated_fes;
173 for (
unsigned int current_fe = 0; current_fe < this->size(); ++current_fe)
178 for (
const auto &other_fe : fes)
181 this->
operator[](current_fe)
182 .compare_for_domination(this->
operator[](other_fe), codim);
188 dominated_fes.insert(current_fe);
190 return dominated_fes;
195 template <
int dim,
int spacedim>
198 const std::set<unsigned int> &fes,
199 const unsigned int codim)
const
210 for (
const auto &fe : fes)
217 for (
const auto ¤t_fe : fes)
222 for (
const auto &other_fe : fes)
223 if (current_fe != other_fe)
226 this->operator[](current_fe)
227 .compare_for_domination(this->
operator[](other_fe), codim);
242 template <
int dim,
int spacedim>
245 const std::set<unsigned int> &fes,
246 const unsigned int codim)
const
257 for (
const auto &fe : fes)
264 for (
const auto ¤t_fe : fes)
269 for (
const auto &other_fe : fes)
270 if (current_fe != other_fe)
273 this->operator[](current_fe)
274 .compare_for_domination(this->
operator[](other_fe), codim);
289 template <
int dim,
int spacedim>
292 const std::set<unsigned int> &fes,
293 const unsigned int codim)
const
295 unsigned int fe_index = find_dominating_fe(fes, codim);
299 const std::set<unsigned int> dominating_fes =
300 find_common_fes(fes, codim);
301 fe_index = find_dominated_fe(dominating_fes, codim);
309 template <
int dim,
int spacedim>
312 const std::set<unsigned int> &fes,
313 const unsigned int codim)
const
315 unsigned int fe_index = find_dominated_fe(fes, codim);
319 const std::set<unsigned int> dominated_fes =
320 find_enclosing_fes(fes, codim);
321 fe_index = find_dominating_fe(dominated_fes, codim);
335 std::vector<std::map<unsigned int, unsigned int>>
336 compute_hp_dof_identities(
337 const std::set<unsigned int> &fes,
338 const std::function<std::vector<std::pair<unsigned int, unsigned int>>(
340 const unsigned int)> & query_identities)
352 const unsigned int fe_index_1 = *fes.begin();
353 const unsigned int fe_index_2 = *(++fes.begin());
354 const auto reduced_identities =
355 query_identities(fe_index_1, fe_index_2);
357 std::vector<std::map<unsigned int, unsigned int>> complete_identities;
359 for (
const auto &reduced_identity : reduced_identities)
364 std::map<unsigned int, unsigned int> complete_identity = {
365 {fe_index_1, reduced_identity.first},
366 {fe_index_2, reduced_identity.second}};
367 complete_identities.emplace_back(std::move(complete_identity));
370 return complete_identities;
381 using Node = std::pair<unsigned int, unsigned int>;
382 using Edge = std::pair<Node, Node>;
383 using Graph = std::set<Edge>;
385 Graph identities_graph;
386 for (
const unsigned int fe_index_1 : fes)
387 for (const unsigned
int fe_index_2 : fes)
388 if (fe_index_1 != fe_index_2)
390 query_identities(fe_index_1, fe_index_2))
400 for (
const auto &edge : identities_graph)
401 Assert(identities_graph.find({edge.second, edge.first}) !=
402 identities_graph.end(),
454 std::vector<std::map<unsigned int, unsigned int>> identities;
455 while (identities_graph.size() > 0)
458 std::set<Node> sub_graph_nodes;
460 sub_graph.emplace(*identities_graph.begin());
461 sub_graph_nodes.emplace(identities_graph.begin()->first);
462 sub_graph_nodes.emplace(identities_graph.begin()->second);
464 for (
const Edge &e : identities_graph)
465 if ((sub_graph_nodes.find(e.first) != sub_graph_nodes.end()) ||
466 (sub_graph_nodes.find(e.second) != sub_graph_nodes.end()))
469 sub_graph_nodes.insert(e.first);
470 sub_graph_nodes.insert(e.second);
475 for (
const Edge &e : sub_graph)
476 identities_graph.erase(e);
482 for (
const auto &edge : sub_graph)
483 Assert(sub_graph.find({edge.second, edge.first}) != sub_graph.end(),
492 for (
const Node &n : sub_graph_nodes)
493 for (const Edge &e : identities_graph)
500 Assert(sub_graph.size() ==
501 sub_graph_nodes.size() * (sub_graph_nodes.size() - 1),
514 identities.emplace_back(sub_graph_nodes.begin(),
515 sub_graph_nodes.end());
516 Assert(identities.back().size() == sub_graph_nodes.size(),
526 template <
int dim,
int spacedim>
527 std::vector<std::map<unsigned int, unsigned int>>
529 const std::set<unsigned int> &fes)
const
531 auto query_vertex_dof_identities = [
this](
const unsigned int fe_index_1,
532 const unsigned int fe_index_2) {
533 return (*
this)[fe_index_1].hp_vertex_dof_identities((*
this)[fe_index_2]);
535 return compute_hp_dof_identities(fes, query_vertex_dof_identities);
540 template <
int dim,
int spacedim>
541 std::vector<std::map<unsigned int, unsigned int>>
543 const std::set<unsigned int> &fes)
const
545 auto query_line_dof_identities = [
this](
const unsigned int fe_index_1,
546 const unsigned int fe_index_2) {
547 return (*
this)[fe_index_1].hp_line_dof_identities((*
this)[fe_index_2]);
549 return compute_hp_dof_identities(fes, query_line_dof_identities);
554 template <
int dim,
int spacedim>
555 std::vector<std::map<unsigned int, unsigned int>>
557 const std::set<unsigned int> &fes,
558 const unsigned int face_no)
const
560 auto query_quad_dof_identities = [
this,
561 face_no](
const unsigned int fe_index_1,
562 const unsigned int fe_index_2) {
563 return (*
this)[fe_index_1].hp_quad_dof_identities((*
this)[fe_index_2],
566 return compute_hp_dof_identities(fes, query_quad_dof_identities);
571 template <
int dim,
int spacedim>
576 const unsigned int)> &next,
579 const unsigned int)> &prev)
582 hierarchy_next = next;
583 hierarchy_prev = prev;
588 template <
int dim,
int spacedim>
593 set_hierarchy(&DefaultHierarchy::next_index,
594 &DefaultHierarchy::previous_index);
599 template <
int dim,
int spacedim>
600 std::vector<unsigned int>
602 const unsigned int fe_index)
const
606 std::deque<unsigned int> sequence = {fe_index};
610 unsigned int front = sequence.front();
611 unsigned int previous;
612 while ((previous = previous_in_hierarchy(front)) != front)
614 sequence.push_front(previous);
617 Assert(sequence.size() <= this->size(),
619 "The registered hierarchy is not terminated: "
620 "previous_in_hierarchy() does not stop at a final index."));
626 unsigned int back = sequence.back();
628 while ((next = next_in_hierarchy(back)) != back)
630 sequence.push_back(next);
633 Assert(sequence.size() <= this->size(),
635 "The registered hierarchy is not terminated: "
636 "next_in_hierarchy() does not stop at a final index."));
640 return {sequence.begin(), sequence.end()};
645 template <
int dim,
int spacedim>
648 const unsigned int fe_index)
const
652 const unsigned int new_fe_index = hierarchy_next(*
this, fe_index);
660 template <
int dim,
int spacedim>
663 const unsigned int fe_index)
const
667 const unsigned int new_fe_index = hierarchy_prev(*
this, fe_index);
675 template <
int dim,
int spacedim>
681 ExcMessage(
"This collection contains no finite element."));
684 const ComponentMask mask = (*this)[0].component_mask(scalar);
688 for (
unsigned int c = 1; c < this->size(); ++c)
695 template <
int dim,
int spacedim>
701 ExcMessage(
"This collection contains no finite element."));
704 const ComponentMask mask = (*this)[0].component_mask(vector);
708 for (
unsigned int c = 1; c < this->size(); ++c)
715 template <
int dim,
int spacedim>
721 ExcMessage(
"This collection contains no finite element."));
724 const ComponentMask mask = (*this)[0].component_mask(sym_tensor);
728 for (
unsigned int c = 1; c < this->size(); ++c)
735 template <
int dim,
int spacedim>
740 ExcMessage(
"This collection contains no finite element."));
743 const ComponentMask mask = (*this)[0].component_mask(block_mask);
747 for (
unsigned int c = 1; c < this->size(); ++c)
748 Assert(mask == (*
this)[c].component_mask(block_mask),
749 ExcMessage(
"Not all elements of this collection agree on what "
750 "the appropriate mask should be."));
756 template <
int dim,
int spacedim>
762 ExcMessage(
"This collection contains no finite element."));
765 const BlockMask mask = (*this)[0].block_mask(scalar);
769 for (
unsigned int c = 1; c < this->size(); ++c)
770 Assert(mask == (*
this)[c].block_mask(scalar),
771 ExcMessage(
"Not all elements of this collection agree on what "
772 "the appropriate mask should be."));
778 template <
int dim,
int spacedim>
784 ExcMessage(
"This collection contains no finite element."));
787 const BlockMask mask = (*this)[0].block_mask(vector);
791 for (
unsigned int c = 1; c < this->size(); ++c)
792 Assert(mask == (*
this)[c].block_mask(vector),
793 ExcMessage(
"Not all elements of this collection agree on what "
794 "the appropriate mask should be."));
800 template <
int dim,
int spacedim>
806 ExcMessage(
"This collection contains no finite element."));
809 const BlockMask mask = (*this)[0].block_mask(sym_tensor);
813 for (
unsigned int c = 1; c < this->size(); ++c)
814 Assert(mask == (*
this)[c].block_mask(sym_tensor),
815 ExcMessage(
"Not all elements of this collection agree on what "
816 "the appropriate mask should be."));
823 template <
int dim,
int spacedim>
829 ExcMessage(
"This collection contains no finite element."));
832 const BlockMask mask = (*this)[0].block_mask(component_mask);
836 for (
unsigned int c = 1; c < this->size(); ++c)
837 Assert(mask == (*
this)[c].block_mask(component_mask),
838 ExcMessage(
"Not all elements of this collection agree on what "
839 "the appropriate mask should be."));
846 template <
int dim,
int spacedim>
850 Assert(this->size() > 0, ExcNoFiniteElements());
852 const unsigned int nb = this->operator[](0).n_blocks();
853 for (
unsigned int i = 1; i < this->size(); ++i)
854 Assert(this->
operator[](i).n_blocks() == nb,
855 ExcMessage(
"Not all finite elements in this collection have "
856 "the same number of components."));
865#include "fe_collection.inst"
unsigned int n_components() const
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const =0
std::vector< std::map< unsigned int, unsigned int > > hp_vertex_dof_identities(const std::set< unsigned int > &fes) const
unsigned int previous_in_hierarchy(const unsigned int fe_index) const
std::vector< unsigned int > get_hierarchy_sequence(const unsigned int fe_index=0) const
unsigned int find_dominating_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
const MappingCollection< dim, spacedim > & get_reference_cell_default_linear_mapping() const
std::set< unsigned int > find_common_fes(const std::set< unsigned int > &fes, const unsigned int codim=0) const
void push_back(const FiniteElement< dim, spacedim > &new_fe)
unsigned int next_in_hierarchy(const unsigned int fe_index) const
void set_default_hierarchy()
std::shared_ptr< MappingCollection< dim, spacedim > > reference_cell_default_linear_mapping
unsigned int find_dominating_fe(const std::set< unsigned int > &fes, const unsigned int codim=0) const
std::set< unsigned int > find_enclosing_fes(const std::set< unsigned int > &fes, const unsigned int codim=0) const
unsigned int find_dominated_fe(const std::set< unsigned int > &fes, const unsigned int codim=0) const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
std::vector< std::map< unsigned int, unsigned int > > hp_line_dof_identities(const std::set< unsigned int > &fes) const
void set_hierarchy(const std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> &next, const std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> &prev)
unsigned int find_dominated_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
unsigned int n_blocks() const
std::vector< std::map< unsigned int, unsigned int > > hp_quad_dof_identities(const std::set< unsigned int > &fes, const unsigned int face_no=0) const
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
static ::ExceptionBase & ExcEmptyObject()
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
@ either_element_can_dominate
@ other_element_dominates
const types::fe_index invalid_fe_index