Reference documentation for deal.II version 9.2.0
\(\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 - 2020 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 
25 #include <memory>
26 
28 
29 namespace hp
30 {
53  template <int dim, int spacedim = dim>
54  class FECollection : public Subscriptor
55  {
56  public:
67  {
74  static unsigned int
75  next_index(const typename hp::FECollection<dim, spacedim> &fe_collection,
76  const unsigned int fe_index)
77  {
78  return ((fe_index + 1) < fe_collection.size()) ? fe_index + 1 :
79  fe_index;
80  }
81 
88  static unsigned int
90  const typename hp::FECollection<dim, spacedim> &fe_collection,
91  const unsigned int fe_index)
92  {
93  (void)fe_collection;
94  return (fe_index > 0) ? fe_index - 1 : fe_index;
95  }
96  };
97 
103  FECollection();
104 
111  explicit FECollection(const FiniteElement<dim, spacedim> &fe);
112 
119  template <class... FETypes>
120  explicit FECollection(const FETypes &... fes);
121 
129  FECollection(const std::vector<const FiniteElement<dim, spacedim> *> &fes);
130 
134  FECollection(const FECollection<dim, spacedim> &) = default;
135 
146  std::is_nothrow_move_constructible<
147  std::vector<std::shared_ptr<const FiniteElement<dim, spacedim>>>>::value
148  &&std::is_nothrow_move_constructible<std::function<
149  unsigned int(const typename hp::FECollection<dim, spacedim> &,
150  const unsigned int)>>::value) = default;
151 
155  FECollection<dim, spacedim> &
156  operator=(FECollection<dim, spacedim> &&) = default; // NOLINT
157 
162  bool
163  operator==(const FECollection<dim, spacedim> &fe_collection) const;
164 
169  bool
170  operator!=(const FECollection<dim, spacedim> &fe_collection) const;
171 
181  void
182  push_back(const FiniteElement<dim, spacedim> &new_fe);
183 
190  const FiniteElement<dim, spacedim> &
191  operator[](const unsigned int index) const;
192 
196  unsigned int
197  size() const;
198 
208  unsigned int
209  n_components() const;
210 
227  unsigned int
228  n_blocks() const;
229 
234  unsigned int
235  max_degree() const;
236 
241  unsigned int
242  max_dofs_per_vertex() const;
243 
248  unsigned int
249  max_dofs_per_line() const;
250 
255  unsigned int
256  max_dofs_per_quad() const;
257 
262  unsigned int
263  max_dofs_per_hex() const;
264 
269  unsigned int
270  max_dofs_per_face() const;
271 
276  unsigned int
277  max_dofs_per_cell() const;
278 
282  std::size_t
283  memory_consumption() const;
284 
285 
303  bool
305 
353  DEAL_II_DEPRECATED unsigned int
354  find_least_face_dominating_fe(const std::set<unsigned int> &fes) const;
355 
377  std::set<unsigned int>
378  find_common_fes(const std::set<unsigned int> &fes,
379  const unsigned int codim = 0) const;
380 
402  std::set<unsigned int>
403  find_enclosing_fes(const std::set<unsigned int> &fes,
404  const unsigned int codim = 0) const;
405 
438  unsigned int
439  find_dominating_fe(const std::set<unsigned int> &fes,
440  const unsigned int codim = 0) const;
441 
474  unsigned int
475  find_dominated_fe(const std::set<unsigned int> &fes,
476  const unsigned int codim = 0) const;
477 
501  unsigned int
502  find_dominating_fe_extended(const std::set<unsigned int> &fes,
503  const unsigned int codim = 0) const;
504 
527  unsigned int
528  find_dominated_fe_extended(const std::set<unsigned int> &fes,
529  const unsigned int codim = 0) const;
530 
543  void
544  set_hierarchy(const std::function<unsigned int(
545  const typename hp::FECollection<dim, spacedim> &,
546  const unsigned int)> &next,
547  const std::function<unsigned int(
548  const typename hp::FECollection<dim, spacedim> &,
549  const unsigned int)> &prev);
550 
558  void
560 
570  unsigned int
571  next_in_hierarchy(const unsigned int fe_index) const;
572 
582  unsigned int
583  previous_in_hierarchy(const unsigned int fe_index) const;
584 
602  component_mask(const FEValuesExtractors::Scalar &scalar) const;
603 
621  component_mask(const FEValuesExtractors::Vector &vector) const;
622 
642  const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
643 
665  component_mask(const BlockMask &block_mask) const;
666 
693  BlockMask
694  block_mask(const FEValuesExtractors::Scalar &scalar) const;
695 
718  BlockMask
719  block_mask(const FEValuesExtractors::Vector &vector) const;
720 
744  BlockMask
745  block_mask(const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
746 
775  BlockMask
777 
789 
797  "No FiniteElement has been found in your FECollection that is "
798  "dominated by all children of a cell you are trying to coarsen!");
799 
804  private:
808  std::vector<std::shared_ptr<const FiniteElement<dim, spacedim>>>
810 
815  std::function<unsigned int(const typename hp::FECollection<dim, spacedim> &,
816  const unsigned int)>
818 
823  std::function<unsigned int(const typename hp::FECollection<dim, spacedim> &,
824  const unsigned int)>
826  };
827 
828 
829 
830  /* --------------- inline functions ------------------- */
831 
832  template <int dim, int spacedim>
833  template <class... FETypes>
834  FECollection<dim, spacedim>::FECollection(const FETypes &... fes)
835  {
836  static_assert(
838  "Not all of the input arguments of this function "
839  "are derived from FiniteElement<dim,spacedim>!");
840 
841  // loop over all of the given arguments and add the finite elements to
842  // this collection. Inlining the definition of fe_pointers causes internal
843  // compiler errors on GCC 7.1.1 so we define it separately:
844  const auto fe_pointers = {
845  (static_cast<const FiniteElement<dim, spacedim> *>(&fes))...};
846  for (const auto p : fe_pointers)
847  push_back(*p);
848  }
849 
850 
851  template <int dim, int spacedim>
852  inline unsigned int
854  {
855  return finite_elements.size();
856  }
857 
858 
859  template <int dim, int spacedim>
860  inline unsigned int
862  {
864 
865  // note that there is no need
866  // here to enforce that indeed
867  // all elements have the same
868  // number of components since we
869  // have already done this when
870  // adding a new element to the
871  // collection.
872 
873  return finite_elements[0]->n_components();
874  }
875 
876 
877 
878  template <int dim, int spacedim>
879  inline bool
881  operator==(const FECollection<dim, spacedim> &fe_collection) const
882  {
883  const unsigned int n_elements = size();
884  if (n_elements != fe_collection.size())
885  return false;
886 
887  for (unsigned int i = 0; i < n_elements; ++i)
888  if (!(*finite_elements[i] == fe_collection[i]))
889  return false;
890 
891  return true;
892  }
893 
894 
895 
896  template <int dim, int spacedim>
897  inline bool
899  operator!=(const FECollection<dim, spacedim> &fe_collection) const
900  {
901  return !(*this == fe_collection);
902  }
903 
904 
905 
906  template <int dim, int spacedim>
908  operator[](const unsigned int index) const
909  {
910  AssertIndexRange(index, finite_elements.size());
911  return *finite_elements[index];
912  }
913 
914 
915 
916  template <int dim, int spacedim>
917  unsigned int
919  {
921 
922  unsigned int max = 0;
923  for (unsigned int i = 0; i < finite_elements.size(); ++i)
924  max = std::max(max, finite_elements[i]->degree);
925 
926  return max;
927  }
928 
929 
930 
931  template <int dim, int spacedim>
932  unsigned int
934  {
936 
937  unsigned int max = 0;
938  for (unsigned int i = 0; i < finite_elements.size(); ++i)
939  if (finite_elements[i]->dofs_per_vertex > max)
940  max = finite_elements[i]->dofs_per_vertex;
941 
942  return max;
943  }
944 
945 
946 
947  template <int dim, int spacedim>
948  unsigned int
950  {
952 
953  unsigned int max = 0;
954  for (unsigned int i = 0; i < finite_elements.size(); ++i)
955  if (finite_elements[i]->dofs_per_line > max)
956  max = finite_elements[i]->dofs_per_line;
957 
958  return max;
959  }
960 
961 
962 
963  template <int dim, int spacedim>
964  unsigned int
966  {
968 
969  unsigned int max = 0;
970  for (unsigned int i = 0; i < finite_elements.size(); ++i)
971  if (finite_elements[i]->dofs_per_quad > max)
972  max = finite_elements[i]->dofs_per_quad;
973 
974  return max;
975  }
976 
977 
978 
979  template <int dim, int spacedim>
980  unsigned int
982  {
984 
985  unsigned int max = 0;
986  for (unsigned int i = 0; i < finite_elements.size(); ++i)
987  if (finite_elements[i]->dofs_per_hex > max)
988  max = finite_elements[i]->dofs_per_hex;
989 
990  return max;
991  }
992 
993 
994 
995  template <int dim, int spacedim>
996  unsigned int
998  {
1000 
1001  unsigned int max = 0;
1002  for (unsigned int i = 0; i < finite_elements.size(); ++i)
1003  if (finite_elements[i]->dofs_per_face > max)
1004  max = finite_elements[i]->dofs_per_face;
1005 
1006  return max;
1007  }
1008 
1009 
1010 
1011  template <int dim, int spacedim>
1012  unsigned int
1014  {
1016 
1017  unsigned int max = 0;
1018  for (unsigned int i = 0; i < finite_elements.size(); ++i)
1019  if (finite_elements[i]->dofs_per_cell > max)
1020  max = finite_elements[i]->dofs_per_cell;
1021 
1022  return max;
1023  }
1024 
1025 
1026  template <int dim, int spacedim>
1027  bool
1029  {
1031  return std::all_of(
1032  finite_elements.cbegin(),
1033  finite_elements.cend(),
1034  [](const std::shared_ptr<const FiniteElement<dim, spacedim>> &fe) {
1035  return fe->hp_constraints_are_implemented();
1036  });
1037  }
1038 
1039 
1040 } // namespace hp
1041 
1043 
1044 #endif
DeclExceptionMsg
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:496
BlockMask
Definition: block_mask.h:74
FEValuesExtractors
Definition: fe_values_extractors.h:82
LinearAlgebraDealII::Vector
Vector< double > Vector
Definition: generic_linear_algebra.h:43
hp::FECollection::operator[]
const FiniteElement< dim, spacedim > & operator[](const unsigned int index) const
Definition: fe_collection.h:908
hp::FECollection::find_dominated_fe
unsigned int find_dominated_fe(const std::set< unsigned int > &fes, const unsigned int codim=0) const
Definition: fe_collection.cc:164
SymmetricTensor
Definition: symmetric_tensor.h:611
is_base_of_all
Definition: template_constraints.h:60
hp::FECollection::next_in_hierarchy
unsigned int next_in_hierarchy(const unsigned int fe_index) const
Definition: fe_collection.cc:330
hp::FECollection::size
unsigned int size() const
Definition: fe_collection.h:853
hp::FECollection::block_mask
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe_collection.cc:441
hp::FECollection::memory_consumption
std::size_t memory_consumption() const
Definition: fe_collection.cc:548
hp::FECollection::n_components
unsigned int n_components() const
Definition: fe_collection.h:861
hp::FECollection::previous_in_hierarchy
unsigned int previous_in_hierarchy(const unsigned int fe_index) const
Definition: fe_collection.cc:345
hp::FECollection::hierarchy_prev
std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> hierarchy_prev
Definition: fe_collection.h:825
hp::FECollection::ExcNoDominatedFiniteElementAmongstChildren
static ::ExceptionBase & ExcNoDominatedFiniteElementAmongstChildren()
AssertIndexRange
#define AssertIndexRange(index, range)
Definition: exceptions.h:1649
hp::FECollection::operator!=
bool operator!=(const FECollection< dim, spacedim > &fe_collection) const
Definition: fe_collection.h:899
ComponentMask
Definition: component_mask.h:83
LinearAlgebra::CUDAWrappers::kernel::set
__global__ void set(Number *val, const Number s, const size_type N)
GridRefinement::coarsen
void coarsen(Triangulation< dim, spacedim > &tria, const Vector< Number > &criteria, const double threshold)
Definition: grid_refinement.cc:88
Subscriptor
Definition: subscriptor.h:62
fe_values_extractors.h
hp::FECollection::set_hierarchy
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)
Definition: fe_collection.cc:302
hp::FECollection::find_dominating_fe_extended
unsigned int find_dominating_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
Definition: fe_collection.cc:210
hp::FECollection::FECollection
FECollection()
Definition: fe_collection.cc:249
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
hp::FECollection::push_back
void push_back(const FiniteElement< dim, spacedim > &new_fe)
Definition: fe_collection.cc:282
fe.h
hp::FECollection::n_blocks
unsigned int n_blocks() const
Definition: fe_collection.cc:531
DEAL_II_DEPRECATED
#define DEAL_II_DEPRECATED
Definition: config.h:98
hp::FECollection::ExcNoFiniteElements
static ::ExceptionBase & ExcNoFiniteElements()
hp::FECollection::hierarchy_next
std::function< unsigned int(const typename hp::FECollection< dim, spacedim > &, const unsigned int)> hierarchy_next
Definition: fe_collection.h:817
component_mask.h
hp::FECollection::max_dofs_per_vertex
unsigned int max_dofs_per_vertex() const
Definition: fe_collection.h:933
hp::FECollection::max_dofs_per_hex
unsigned int max_dofs_per_hex() const
Definition: fe_collection.h:981
FiniteElement
Definition: fe.h:648
hp::FECollection::hp_constraints_are_implemented
bool hp_constraints_are_implemented() const
Definition: fe_collection.h:1028
hp::FECollection::DefaultHierarchy::next_index
static unsigned int next_index(const typename hp::FECollection< dim, spacedim > &fe_collection, const unsigned int fe_index)
Definition: fe_collection.h:75
value
static const bool value
Definition: dof_tools_constraints.cc:433
hp::FECollection::max_dofs_per_cell
unsigned int max_dofs_per_cell() const
Definition: fe_collection.h:1013
hp::FECollection::set_default_hierarchy
void set_default_hierarchy()
Definition: fe_collection.cc:319
hp::FECollection::max_degree
unsigned int max_degree() const
Definition: fe_collection.h:918
Assert
#define Assert(cond, exc)
Definition: exceptions.h:1419
hp::FECollection::finite_elements
std::vector< std::shared_ptr< const FiniteElement< dim, spacedim > > > finite_elements
Definition: fe_collection.h:809
hp::FECollection::DefaultHierarchy::previous_index
static unsigned int previous_index(const typename hp::FECollection< dim, spacedim > &fe_collection, const unsigned int fe_index)
Definition: fe_collection.h:89
hp::FECollection::find_dominated_fe_extended
unsigned int find_dominated_fe_extended(const std::set< unsigned int > &fes, const unsigned int codim=0) const
Definition: fe_collection.cc:230
DeclException0
#define DeclException0(Exception0)
Definition: exceptions.h:473
hp::FECollection::find_enclosing_fes
std::set< unsigned int > find_enclosing_fes(const std::set< unsigned int > &fes, const unsigned int codim=0) const
Definition: fe_collection.cc:78
hp
Definition: hp.h:117
hp::FECollection::max_dofs_per_face
unsigned int max_dofs_per_face() const
Definition: fe_collection.h:997
hp::FECollection::max_dofs_per_quad
unsigned int max_dofs_per_quad() const
Definition: fe_collection.h:965
hp::FECollection::max_dofs_per_line
unsigned int max_dofs_per_line() const
Definition: fe_collection.h:949
hp::FECollection::component_mask
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe_collection.cc:360
config.h
hp::FECollection::find_dominating_fe
unsigned int find_dominating_fe(const std::set< unsigned int > &fes, const unsigned int codim=0) const
Definition: fe_collection.cc:118
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
hp::FECollection::operator==
bool operator==(const FECollection< dim, spacedim > &fe_collection) const
Definition: fe_collection.h:881
hp::FECollection::find_common_fes
std::set< unsigned int > find_common_fes(const std::set< unsigned int > &fes, const unsigned int codim=0) const
Definition: fe_collection.cc:38
hp::FECollection::find_least_face_dominating_fe
unsigned int find_least_face_dominating_fe(const std::set< unsigned int > &fes) const
Definition: fe_collection.cc:27
hp::FECollection::DefaultHierarchy
Definition: fe_collection.h:66
Utilities::MPI::max
T max(const T &t, const MPI_Comm &mpi_communicator)
hp::FECollection
Definition: fe_collection.h:54