Reference documentation for deal.II version 9.3.3
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
fe_collection.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2003 - 2021 by the deal.II authors
4//
5// This file is part of the deal.II library.
6//
7// The deal.II library is free software; you can use it, redistribute
8// it, and/or modify it under the terms of the GNU Lesser General
9// Public License as published by the Free Software Foundation; either
10// version 2.1 of the License, or (at your option) any later version.
11// The full text of the license can be found in the file LICENSE.md at
12// the top level directory of deal.II.
13//
14// ---------------------------------------------------------------------
15
16#ifndef dealii_fe_collection_h
17#define dealii_fe_collection_h
18
19#include <deal.II/base/config.h>
20
22#include <deal.II/fe/fe.h>
24
26
27#include <memory>
28
30
31namespace hp
32{
52 template <int dim, int spacedim = dim>
53 class FECollection : public Collection<FiniteElement<dim, spacedim>>
54 {
55 public:
66 {
73 static unsigned int
74 next_index(const typename hp::FECollection<dim, spacedim> &fe_collection,
75 const unsigned int fe_index)
76 {
77 return ((fe_index + 1) < fe_collection.size()) ? fe_index + 1 :
78 fe_index;
79 }
80
87 static unsigned int
89 const typename hp::FECollection<dim, spacedim> &fe_collection,
90 const unsigned int fe_index)
91 {
92 (void)fe_collection;
93 return (fe_index > 0) ? fe_index - 1 : fe_index;
94 }
95 };
96
102 FECollection();
103
110 explicit FECollection(const FiniteElement<dim, spacedim> &fe);
111
118 template <class... FETypes>
119 explicit FECollection(const FETypes &... fes);
120
128 FECollection(const std::vector<const FiniteElement<dim, spacedim> *> &fes);
129
134
145 std::is_nothrow_move_constructible<
146 std::vector<std::shared_ptr<const FiniteElement<dim, spacedim>>>>::value
147 &&std::is_nothrow_move_constructible<std::function<
148 unsigned int(const typename hp::FECollection<dim, spacedim> &,
149 const unsigned int)>>::value) = default;
150
154 FECollection<dim, spacedim> &
155 operator=(FECollection<dim, spacedim> &&) = default; // NOLINT
156
161 bool
162 operator==(const FECollection<dim, spacedim> &fe_collection) const;
163
168 bool
169 operator!=(const FECollection<dim, spacedim> &fe_collection) const;
170
180 void
181 push_back(const FiniteElement<dim, spacedim> &new_fe);
182
192 unsigned int
194
211 unsigned int
212 n_blocks() const;
213
218 unsigned int
219 max_degree() const;
220
225 unsigned int
227
232 unsigned int
234
239 unsigned int
241
246 unsigned int
248
253 unsigned int
255
260 unsigned int
262
263
281 bool
283
305 std::set<unsigned int>
306 find_common_fes(const std::set<unsigned int> &fes,
307 const unsigned int codim = 0) const;
308
330 std::set<unsigned int>
331 find_enclosing_fes(const std::set<unsigned int> &fes,
332 const unsigned int codim = 0) const;
333
366 unsigned int
367 find_dominating_fe(const std::set<unsigned int> &fes,
368 const unsigned int codim = 0) const;
369
402 unsigned int
403 find_dominated_fe(const std::set<unsigned int> &fes,
404 const unsigned int codim = 0) const;
405
429 unsigned int
430 find_dominating_fe_extended(const std::set<unsigned int> &fes,
431 const unsigned int codim = 0) const;
432
455 unsigned int
456 find_dominated_fe_extended(const std::set<unsigned int> &fes,
457 const unsigned int codim = 0) const;
458
471 void
472 set_hierarchy(const std::function<unsigned int(
473 const typename hp::FECollection<dim, spacedim> &,
474 const unsigned int)> &next,
475 const std::function<unsigned int(
476 const typename hp::FECollection<dim, spacedim> &,
477 const unsigned int)> &prev);
478
486 void
488
508 std::vector<unsigned int>
509 get_hierarchy_sequence(const unsigned int fe_index = 0) const;
510
520 unsigned int
521 next_in_hierarchy(const unsigned int fe_index) const;
522
532 unsigned int
533 previous_in_hierarchy(const unsigned int fe_index) const;
534
552 component_mask(const FEValuesExtractors::Scalar &scalar) const;
553
571 component_mask(const FEValuesExtractors::Vector &vector) const;
572
592 const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
593
615 component_mask(const BlockMask &block_mask) const;
616
644 block_mask(const FEValuesExtractors::Scalar &scalar) const;
645
669 block_mask(const FEValuesExtractors::Vector &vector) const;
670
695 block_mask(const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
696
727
739
744 private:
749 std::function<unsigned int(const typename hp::FECollection<dim, spacedim> &,
750 const unsigned int)>
752
757 std::function<unsigned int(const typename hp::FECollection<dim, spacedim> &,
758 const unsigned int)>
760 };
761
762
763
764 /* --------------- inline functions ------------------- */
765
766 template <int dim, int spacedim>
767 template <class... FETypes>
768 FECollection<dim, spacedim>::FECollection(const FETypes &... fes)
769 {
770 static_assert(
772 "Not all of the input arguments of this function "
773 "are derived from FiniteElement<dim,spacedim>!");
774
775 // loop over all of the given arguments and add the finite elements to
776 // this collection. Inlining the definition of fe_pointers causes internal
777 // compiler errors on GCC 7.1.1 so we define it separately:
778 const auto fe_pointers = {
779 (static_cast<const FiniteElement<dim, spacedim> *>(&fes))...};
780 for (const auto p : fe_pointers)
781 push_back(*p);
782 }
783
784
785
786 template <int dim, int spacedim>
787 inline unsigned int
789 {
790 Assert(this->size() > 0, ExcNoFiniteElements());
791
792 // note that there is no need
793 // here to enforce that indeed
794 // all elements have the same
795 // number of components since we
796 // have already done this when
797 // adding a new element to the
798 // collection.
799
800 return this->operator[](0).n_components();
801 }
802
803
804
805 template <int dim, int spacedim>
806 inline bool
808 operator==(const FECollection<dim, spacedim> &fe_collection) const
809 {
810 const unsigned int n_elements = this->size();
811 if (n_elements != fe_collection.size())
812 return false;
813
814 for (unsigned int i = 0; i < n_elements; ++i)
815 if (!(this->operator[](i) == fe_collection[i]))
816 return false;
817
818 return true;
819 }
820
821
822
823 template <int dim, int spacedim>
824 inline bool
826 operator!=(const FECollection<dim, spacedim> &fe_collection) const
827 {
828 return !(*this == fe_collection);
829 }
830
831
832
833 template <int dim, int spacedim>
834 unsigned int
836 {
837 Assert(this->size() > 0, ExcNoFiniteElements());
838
839 unsigned int max = 0;
840 for (unsigned int i = 0; i < this->size(); ++i)
841 max = std::max(max, this->operator[](i).degree);
842
843 return max;
844 }
845
846
847
848 template <int dim, int spacedim>
849 unsigned int
851 {
852 Assert(this->size() > 0, ExcNoFiniteElements());
853
854 unsigned int max = 0;
855 for (unsigned int i = 0; i < this->size(); ++i)
856 max = std::max(max, this->operator[](i).n_dofs_per_vertex());
857
858 return max;
859 }
860
861
862
863 template <int dim, int spacedim>
864 unsigned int
866 {
867 Assert(this->size() > 0, ExcNoFiniteElements());
868
869 unsigned int max = 0;
870 for (unsigned int i = 0; i < this->size(); ++i)
871 max = std::max(max, this->operator[](i).n_dofs_per_line());
872
873 return max;
874 }
875
876
877
878 template <int dim, int spacedim>
879 unsigned int
881 {
882 Assert(this->size() > 0, ExcNoFiniteElements());
883
884 unsigned int max = 0;
885 for (unsigned int i = 0; i < this->size(); ++i)
886 max = std::max(max, this->operator[](i).max_dofs_per_quad());
887
888 return max;
889 }
890
891
892
893 template <int dim, int spacedim>
894 unsigned int
896 {
897 Assert(this->size() > 0, ExcNoFiniteElements());
898
899 unsigned int max = 0;
900 for (unsigned int i = 0; i < this->size(); ++i)
901 max = std::max(max, this->operator[](i).n_dofs_per_hex());
902
903 return max;
904 }
905
906
907
908 template <int dim, int spacedim>
909 unsigned int
911 {
912 Assert(this->size() > 0, ExcNoFiniteElements());
913
914 unsigned int max = 0;
915 for (unsigned int i = 0; i < this->size(); ++i)
916 max = std::max(max, this->operator[](i).max_dofs_per_face());
917
918 return max;
919 }
920
921
922
923 template <int dim, int spacedim>
924 unsigned int
926 {
927 Assert(this->size() > 0, ExcNoFiniteElements());
928
929 unsigned int max = 0;
930 for (unsigned int i = 0; i < this->size(); ++i)
931 max = std::max(max, this->operator[](i).n_dofs_per_cell());
932
933 return max;
934 }
935
936
937 template <int dim, int spacedim>
938 bool
940 {
941 Assert(this->size() > 0, ExcNoFiniteElements());
942
943 for (unsigned int i = 0; i < this->size(); ++i)
944 if (this->operator[](i).hp_constraints_are_implemented() == false)
945 return false;
946
947 return true;
948 }
949
950
951} // namespace hp
952
954
955#endif
unsigned int n_components() const
Definition: vector.h:110
unsigned int size() const
Definition: collection.h:109
const FiniteElement< dim, dim > & operator[](const unsigned int index) const
Definition: collection.h:117
unsigned int previous_in_hierarchy(const unsigned int fe_index) const
unsigned int max_dofs_per_line() const
unsigned int max_dofs_per_hex() const
std::vector< unsigned int > get_hierarchy_sequence(const unsigned int fe_index=0) const
unsigned int max_dofs_per_vertex() const
unsigned int find_dominating_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
bool hp_constraints_are_implemented() const
bool operator==(const FECollection< dim, spacedim > &fe_collection) const
std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> hierarchy_next
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()
FECollection(const FETypes &... fes)
unsigned int max_degree() const
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
unsigned int max_dofs_per_face() const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
bool operator!=(const FECollection< dim, spacedim > &fe_collection) 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)
FECollection(FECollection< dim, spacedim > &&) noexcept(std::is_nothrow_move_constructible< std::vector< std::shared_ptr< const FiniteElement< dim, spacedim > > > >::value &&std::is_nothrow_move_constructible< std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> >::value)=default
unsigned int max_dofs_per_quad() const
unsigned int find_dominated_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> hierarchy_prev
unsigned int n_blocks() const
FECollection(const FECollection< dim, spacedim > &)=default
unsigned int n_components() const
unsigned int max_dofs_per_cell() const
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
__global__ void set(Number *val, const Number s, const size_type N)
#define DeclException0(Exception0)
Definition: exceptions.h:470
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcNoFiniteElements()
Definition: hp.h:118
STL namespace.
::VectorizedArray< Number, width > max(const ::VectorizedArray< Number, width > &, const ::VectorizedArray< Number, width > &)
static unsigned int next_index(const typename hp::FECollection< dim, spacedim > &fe_collection, const unsigned int fe_index)
Definition: fe_collection.h:74
static unsigned int previous_index(const typename hp::FECollection< dim, spacedim > &fe_collection, const unsigned int fe_index)
Definition: fe_collection.h:88