Reference documentation for deal.II version 9.0.0
fe.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 1998 - 2018 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 at
12 // the top level of the deal.II distribution.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_fe_h
17 #define dealii_fe_h
18 
19 #include <deal.II/base/config.h>
20 #include <deal.II/fe/fe_base.h>
21 #include <deal.II/fe/fe_values_extractors.h>
22 #include <deal.II/fe/fe_update_flags.h>
23 #include <deal.II/fe/component_mask.h>
24 #include <deal.II/fe/block_mask.h>
25 #include <deal.II/fe/mapping.h>
26 
27 #include <memory>
28 
29 
30 DEAL_II_NAMESPACE_OPEN
31 
32 template <int dim, int spacedim> class FEValuesBase;
33 template <int dim, int spacedim> class FEValues;
34 template <int dim, int spacedim> class FEFaceValues;
35 template <int dim, int spacedim> class FESubfaceValues;
36 template <int dim, int spacedim> class FESystem;
37 
632 template <int dim, int spacedim=dim>
633 class FiniteElement : public Subscriptor,
634  public FiniteElementData<dim>
635 {
636 public:
640  static const unsigned int space_dimension = spacedim;
641 
669  {
670  private:
674  InternalDataBase (const InternalDataBase &) = delete;
675 
676  public:
681  InternalDataBase ();
682 
686  virtual ~InternalDataBase () = default;
687 
703 
707  virtual std::size_t memory_consumption () const;
708  };
709 
710 public:
753  FiniteElement (const FiniteElementData<dim> &fe_data,
754  const std::vector<bool> &restriction_is_additive_flags,
755  const std::vector<ComponentMask> &nonzero_components);
756 
760  FiniteElement (FiniteElement<dim, spacedim> &&) = default; // NOLINT
761 
765  FiniteElement (const FiniteElement<dim, spacedim> &) = default;
766 
771  virtual ~FiniteElement () = default;
772 
781  std::pair<std::unique_ptr<FiniteElement<dim, spacedim> >, unsigned int>
782  operator^ (const unsigned int multiplicity) const;
783 
795  virtual
796  std::unique_ptr<FiniteElement<dim,spacedim> >
797  clone() const = 0;
798 
809  virtual std::string get_name () const = 0;
810 
832  const FiniteElement<dim,spacedim> &operator[] (const unsigned int fe_index) const;
833 
859  virtual double shape_value (const unsigned int i,
860  const Point<dim> &p) const;
861 
868  virtual double shape_value_component (const unsigned int i,
869  const Point<dim> &p,
870  const unsigned int component) const;
871 
893  virtual Tensor<1,dim> shape_grad (const unsigned int i,
894  const Point<dim> &p) const;
895 
902  virtual Tensor<1,dim> shape_grad_component (const unsigned int i,
903  const Point<dim> &p,
904  const unsigned int component) const;
905 
927  virtual Tensor<2,dim> shape_grad_grad (const unsigned int i,
928  const Point<dim> &p) const;
929 
936  virtual Tensor<2,dim> shape_grad_grad_component (const unsigned int i,
937  const Point<dim> &p,
938  const unsigned int component) const;
939 
961  virtual Tensor<3,dim> shape_3rd_derivative (const unsigned int i,
962  const Point<dim> &p) const;
963 
970  virtual Tensor<3,dim> shape_3rd_derivative_component (const unsigned int i,
971  const Point<dim> &p,
972  const unsigned int component) const;
973 
995  virtual Tensor<4,dim> shape_4th_derivative (const unsigned int i,
996  const Point<dim> &p) const;
997 
1004  virtual Tensor<4,dim> shape_4th_derivative_component (const unsigned int i,
1005  const Point<dim> &p,
1006  const unsigned int component) const;
1017  virtual bool has_support_on_face (const unsigned int shape_index,
1018  const unsigned int face_index) const;
1019 
1021 
1041  virtual const FullMatrix<double> &
1042  get_restriction_matrix (const unsigned int child,
1043  const RefinementCase<dim> &refinement_case=RefinementCase<dim>::isotropic_refinement) const;
1044 
1074  virtual const FullMatrix<double> &
1075  get_prolongation_matrix (const unsigned int child,
1076  const RefinementCase<dim> &refinement_case=RefinementCase<dim>::isotropic_refinement) const;
1077 
1099  bool prolongation_is_implemented () const;
1100 
1117 
1139  bool restriction_is_implemented () const;
1140 
1157 
1158 
1167  bool restriction_is_additive (const unsigned int index) const;
1168 
1180  const FullMatrix<double> &constraints (const ::internal::SubfaceCase<dim> &subface_case=::internal::SubfaceCase<dim>::case_isotropic) const;
1181 
1197  bool constraints_are_implemented (const ::internal::SubfaceCase<dim> &subface_case=::internal::SubfaceCase<dim>::case_isotropic) const;
1198 
1199 
1221  virtual bool hp_constraints_are_implemented () const;
1222 
1223 
1235  virtual void
1237  FullMatrix<double> &matrix) const;
1239 
1257  virtual void
1259  FullMatrix<double> &matrix) const;
1260 
1261 
1273  virtual void
1275  const unsigned int subface,
1276  FullMatrix<double> &matrix) const;
1278 
1279 
1300  virtual
1301  std::vector<std::pair<unsigned int, unsigned int> >
1302  hp_vertex_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
1303 
1308  virtual
1309  std::vector<std::pair<unsigned int, unsigned int> >
1310  hp_line_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
1311 
1316  virtual
1317  std::vector<std::pair<unsigned int, unsigned int> >
1318  hp_quad_dof_identities (const FiniteElement<dim,spacedim> &fe_other) const;
1319 
1329  virtual
1332 
1334 
1365  virtual
1366  bool operator == (const FiniteElement<dim,spacedim> &fe) const;
1367 
1371  bool operator != (const FiniteElement<dim,spacedim> &) const;
1372 
1407  std::pair<unsigned int, unsigned int>
1408  system_to_component_index (const unsigned int index) const;
1409 
1419  unsigned int component_to_system_index(const unsigned int component,
1420  const unsigned int index) const;
1421 
1431  std::pair<unsigned int, unsigned int>
1432  face_system_to_component_index (const unsigned int index) const;
1433 
1442  unsigned int adjust_quad_dof_index_for_face_orientation (const unsigned int index,
1443  const bool face_orientation,
1444  const bool face_flip,
1445  const bool face_rotation) const;
1446 
1501  virtual
1502  unsigned int face_to_cell_index (const unsigned int face_dof_index,
1503  const unsigned int face,
1504  const bool face_orientation = true,
1505  const bool face_flip = false,
1506  const bool face_rotation = false) const;
1507 
1515  unsigned int adjust_line_dof_index_for_line_orientation (const unsigned int index,
1516  const bool line_orientation) const;
1517 
1534  const ComponentMask &
1535  get_nonzero_components (const unsigned int i) const;
1536 
1547  unsigned int
1548  n_nonzero_components (const unsigned int i) const;
1549 
1558  bool is_primitive () const;
1559 
1570  bool
1571  is_primitive (const unsigned int i) const;
1572 
1584  unsigned int n_base_elements () const;
1585 
1590  virtual
1592  base_element (const unsigned int index) const;
1593 
1600  unsigned int
1601  element_multiplicity (const unsigned int index) const;
1602 
1675  get_sub_fe (const ComponentMask &mask) const;
1676 
1684  virtual
1686  get_sub_fe (const unsigned int first_component,
1687  const unsigned int n_selected_components) const;
1688 
1710  std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1711  system_to_base_index (const unsigned int index) const;
1712 
1721  std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1722  face_system_to_base_index (const unsigned int index) const;
1723 
1728  types::global_dof_index first_block_of_base (const unsigned int b) const;
1729 
1742  std::pair<unsigned int, unsigned int>
1743  component_to_base_index (const unsigned int component) const;
1744 
1745 
1750  std::pair<unsigned int,unsigned int>
1751  block_to_base_index (const unsigned int block) const;
1752 
1756  std::pair<unsigned int,types::global_dof_index>
1757  system_to_block_index (const unsigned int component) const;
1758 
1762  unsigned int
1763  component_to_block_index (const unsigned int component) const;
1764 
1766 
1785  component_mask (const FEValuesExtractors::Scalar &scalar) const;
1786 
1800  component_mask (const FEValuesExtractors::Vector &vector) const;
1801 
1816  component_mask (const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
1817 
1833  component_mask (const BlockMask &block_mask) const;
1834 
1855  BlockMask
1856  block_mask (const FEValuesExtractors::Scalar &scalar) const;
1857 
1874  BlockMask
1875  block_mask (const FEValuesExtractors::Vector &vector) const;
1876 
1894  BlockMask
1895  block_mask (const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
1896 
1919  BlockMask
1920  block_mask (const ComponentMask &component_mask) const;
1921 
1937  virtual std::pair<Table<2,bool>,std::vector<unsigned int> >
1938  get_constant_modes () const;
1939 
1941 
1977  const std::vector<Point<dim> > &
1978  get_unit_support_points () const;
1979 
1997  bool has_support_points () const;
1998 
2011  virtual
2012  Point<dim>
2013  unit_support_point (const unsigned int index) const;
2014 
2041  const std::vector<Point<dim-1> > &
2043 
2052  bool has_face_support_points () const;
2053 
2058  virtual
2059  Point<dim-1>
2060  unit_face_support_point (const unsigned int index) const;
2061 
2074  const std::vector<Point<dim> > &
2076 
2086  bool has_generalized_support_points () const;
2087 
2096  DEAL_II_DEPRECATED
2097  const std::vector<Point<dim-1> > &
2099 
2112  DEAL_II_DEPRECATED
2113  bool
2115 
2160  get_associated_geometry_primitive (const unsigned int cell_dof_index) const;
2161 
2162 
2240  virtual
2241  void
2242  convert_generalized_support_point_values_to_dof_values (const std::vector<Vector<double> > &support_point_values,
2243  std::vector<double> &nodal_values) const;
2244 
2246 
2255  virtual std::size_t memory_consumption () const;
2256 
2263  int,
2264  << "The shape function with index " << arg1
2265  << " is not primitive, i.e. it is vector-valued and "
2266  << "has more than one non-zero vector component. This "
2267  << "function cannot be called for these shape functions. "
2268  << "Maybe you want to use the same function with the "
2269  << "_component suffix?");
2282  "You are trying to access the values or derivatives of shape functions "
2283  "on the reference cell of an element that does not define its shape "
2284  "functions through mapping from the reference cell. Consequently, "
2285  "you cannot ask for shape function values or derivatives on the "
2286  "reference cell.");
2287 
2295  "You are trying to access the support points of a finite "
2296  "element that either has no support points at all, or for "
2297  "which the corresponding tables have not been implemented.");
2298 
2306  "You are trying to access the matrices that describe how "
2307  "to embed a finite element function on one cell into the "
2308  "finite element space on one of its children (i.e., the "
2309  "'embedding' or 'prolongation' matrices). However, the "
2310  "current finite element can either not define this sort of "
2311  "operation, or it has not yet been implemented.");
2312 
2321  "You are trying to access the matrices that describe how "
2322  "to restrict a finite element function from the children "
2323  "of one cell to the finite element space defined on their "
2324  "parent (i.e., the 'restriction' or 'projection' matrices). "
2325  "However, the current finite element can either not define "
2326  "this sort of operation, or it has not yet been "
2327  "implemented.");
2328 
2334  int, int,
2335  << "The interface matrix has a size of " << arg1
2336  << "x" << arg2
2337  << ", which is not reasonable for the current element "
2338  "in the present dimension.");
2344 
2345 protected:
2346 
2360  void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false,
2361  const bool isotropic_prolongation_only=false);
2362 
2375  std::vector<std::vector<FullMatrix<double> > > restriction;
2376 
2389  std::vector<std::vector<FullMatrix<double> > > prolongation;
2390 
2402 
2413  std::vector<Point<dim> > unit_support_points;
2414 
2420  std::vector<Point<dim-1> > unit_face_support_points;
2421 
2426  std::vector<Point<dim> > generalized_support_points;
2427 
2433 
2450 
2465 
2469  std::vector<std::pair<unsigned int, unsigned int> > system_to_component_table;
2470 
2480  std::vector<std::pair<unsigned int, unsigned int> > face_system_to_component_table;
2481 
2498  std::vector<std::pair<std::pair<unsigned int,unsigned int>,unsigned int> >
2500 
2504  std::vector<std::pair<std::pair<unsigned int,unsigned int>,unsigned int> >
2506 
2512 
2533  std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int> >
2535 
2541  const std::vector<bool> restriction_is_additive_flags;
2542 
2550  const std::vector<ComponentMask> nonzero_components;
2551 
2559  const std::vector<unsigned int> n_nonzero_components_table;
2560 
2567 
2580  interface_constraints_size () const;
2581 
2587  static
2588  std::vector<unsigned int>
2589  compute_n_nonzero_components (const std::vector<ComponentMask> &nonzero_components);
2590 
2611  virtual
2612  UpdateFlags
2613  requires_update_flags (const UpdateFlags update_flags) const = 0;
2614 
2691  virtual
2692  std::unique_ptr<InternalDataBase>
2693  get_data (const UpdateFlags update_flags,
2694  const Mapping<dim,spacedim> &mapping,
2695  const Quadrature<dim> &quadrature,
2697 
2739  virtual
2740  std::unique_ptr<InternalDataBase>
2741  get_face_data (const UpdateFlags update_flags,
2742  const Mapping<dim,spacedim> &mapping,
2743  const Quadrature<dim-1> &quadrature,
2745 
2787  virtual
2788  std::unique_ptr<InternalDataBase>
2789  get_subface_data (const UpdateFlags update_flags,
2790  const Mapping<dim,spacedim> &mapping,
2791  const Quadrature<dim-1> &quadrature,
2793 
2874  virtual
2875  void
2877  const CellSimilarity::Similarity cell_similarity,
2878  const Quadrature<dim> &quadrature,
2879  const Mapping<dim,spacedim> &mapping,
2880  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
2881  const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
2882  const InternalDataBase &fe_internal,
2884 
2927  virtual
2928  void
2930  const unsigned int face_no,
2931  const Quadrature<dim-1> &quadrature,
2932  const Mapping<dim,spacedim> &mapping,
2933  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
2934  const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
2935  const InternalDataBase &fe_internal,
2937 
2983  virtual
2984  void
2986  const unsigned int face_no,
2987  const unsigned int sub_no,
2988  const Quadrature<dim-1> &quadrature,
2989  const Mapping<dim,spacedim> &mapping,
2990  const typename Mapping<dim,spacedim>::InternalDataBase &mapping_internal,
2991  const ::internal::FEValuesImplementation::MappingRelatedData<dim, spacedim> &mapping_data,
2992  const InternalDataBase &fe_internal,
2994 
2995  friend class InternalDataBase;
2996  friend class FEValuesBase<dim,spacedim>;
2997  friend class FEValues<dim,spacedim>;
2998  friend class FEFaceValues<dim,spacedim>;
2999  friend class FESubfaceValues<dim,spacedim>;
3000  friend class FESystem<dim,spacedim>;
3001 
3002  // explicitly check for sensible template arguments, but not on windows
3003  // because MSVC creates bogus warnings during normal compilation
3004 #ifndef DEAL_II_MSVC
3005  static_assert (dim<=spacedim,
3006  "The dimension <dim> of a FiniteElement must be less than or "
3007  "equal to the space dimension <spacedim> in which it lives.");
3008 #endif
3009 
3010 };
3011 
3012 
3013 //----------------------------------------------------------------------//
3014 
3015 
3016 template <int dim, int spacedim>
3017 inline
3019 FiniteElement<dim,spacedim>::operator[] (const unsigned int fe_index) const
3020 {
3021  (void)fe_index;
3022  Assert (fe_index == 0,
3023  ExcMessage ("A fe_index of zero is the only index allowed here"));
3024  return *this;
3025 }
3026 
3027 
3028 
3029 template <int dim, int spacedim>
3030 inline
3031 std::pair<unsigned int,unsigned int>
3033 {
3034  Assert (index < system_to_component_table.size(),
3035  ExcIndexRange(index, 0, system_to_component_table.size()));
3036  Assert (is_primitive (index),
3038  return system_to_component_table[index];
3039 }
3040 
3041 
3042 
3043 template <int dim, int spacedim>
3044 inline
3045 unsigned int
3047 {
3048  return base_to_block_indices.size();
3049 }
3050 
3051 
3052 
3053 template <int dim, int spacedim>
3054 inline
3055 unsigned int
3057 {
3058  return static_cast<unsigned int>(base_to_block_indices.block_size(index));
3059 }
3060 
3061 
3062 
3063 template <int dim, int spacedim>
3064 inline
3065 unsigned int
3067  const unsigned int index) const
3068 {
3069  AssertIndexRange(component, this->n_components());
3070  const std::vector<std::pair<unsigned int, unsigned int> >::const_iterator
3071  it = std::find(system_to_component_table.begin(), system_to_component_table.end(),
3072  std::pair<unsigned int, unsigned int>(component, index));
3073 
3074  Assert(it != system_to_component_table.end(),
3075  ExcMessage ("You are asking for the number of the shape function "
3076  "within a system element that corresponds to vector "
3077  "component " + Utilities::int_to_string(component) + " and within this to "
3078  "index " + Utilities::int_to_string(index) + ". But no such "
3079  "shape function exists."));
3080  return std::distance(system_to_component_table.begin(), it);
3081 }
3082 
3083 
3084 
3085 template <int dim, int spacedim>
3086 inline
3087 std::pair<unsigned int,unsigned int>
3089 {
3090  Assert(index < face_system_to_component_table.size(),
3091  ExcIndexRange(index, 0, face_system_to_component_table.size()));
3092 
3093  // in debug mode, check whether the
3094  // function is primitive, since
3095  // otherwise the result may have no
3096  // meaning
3097  //
3098  // since the primitivity tables are
3099  // all geared towards cell dof
3100  // indices, rather than face dof
3101  // indices, we have to work a
3102  // little bit...
3103  //
3104  // in 1d, the face index is equal
3105  // to the cell index
3106  Assert (is_primitive(this->face_to_cell_index(index, 0)),
3108 
3109  return face_system_to_component_table[index];
3110 }
3111 
3112 
3113 
3114 
3115 template <int dim, int spacedim>
3116 inline
3117 std::pair<std::pair<unsigned int,unsigned int>,unsigned int>
3119 {
3120  Assert (index < system_to_base_table.size(),
3121  ExcIndexRange(index, 0, system_to_base_table.size()));
3122  return system_to_base_table[index];
3123 }
3124 
3125 
3126 
3127 
3128 template <int dim, int spacedim>
3129 inline
3130 std::pair<std::pair<unsigned int,unsigned int>,unsigned int>
3132 {
3133  Assert(index < face_system_to_base_table.size(),
3134  ExcIndexRange(index, 0, face_system_to_base_table.size()));
3135  return face_system_to_base_table[index];
3136 }
3137 
3138 
3139 
3140 template <int dim, int spacedim>
3141 inline
3144 {
3145  return base_to_block_indices.block_start(index);
3146 }
3147 
3148 
3149 
3150 template <int dim, int spacedim>
3151 inline
3152 std::pair<unsigned int,unsigned int>
3154 {
3155  Assert(index < component_to_base_table.size(),
3156  ExcIndexRange(index, 0, component_to_base_table.size()));
3157 
3158  return component_to_base_table[index].first;
3159 }
3160 
3161 
3162 
3163 template <int dim, int spacedim>
3164 inline
3165 std::pair<unsigned int,unsigned int>
3167 {
3168  return base_to_block_indices.global_to_local(index);
3169 }
3170 
3171 
3172 
3173 template <int dim, int spacedim>
3174 inline
3175 std::pair<unsigned int,types::global_dof_index>
3177 {
3178  Assert (index < this->dofs_per_cell,
3179  ExcIndexRange(index, 0, this->dofs_per_cell));
3180  // The block is computed simply as
3181  // first block of this base plus
3182  // the index within the base blocks
3183  return std::pair<unsigned int, types::global_dof_index>(
3184  first_block_of_base(system_to_base_table[index].first.first)
3185  + system_to_base_table[index].first.second,
3186  system_to_base_table[index].second);
3187 }
3188 
3189 
3190 
3191 template <int dim, int spacedim>
3192 inline
3193 bool
3195 {
3196  Assert(index < this->dofs_per_cell,
3197  ExcIndexRange(index, 0, this->dofs_per_cell));
3198  return restriction_is_additive_flags[index];
3199 }
3200 
3201 
3202 
3203 template <int dim, int spacedim>
3204 inline
3205 const ComponentMask &
3207 {
3208  Assert (i < this->dofs_per_cell, ExcIndexRange (i, 0, this->dofs_per_cell));
3209  return nonzero_components[i];
3210 }
3211 
3212 
3213 
3214 template <int dim, int spacedim>
3215 inline
3216 unsigned int
3218 {
3219  Assert (i < this->dofs_per_cell, ExcIndexRange (i, 0, this->dofs_per_cell));
3220  return n_nonzero_components_table[i];
3221 }
3222 
3223 
3224 
3225 template <int dim, int spacedim>
3226 inline
3227 bool
3229 {
3230  return cached_primitivity;
3231 }
3232 
3233 
3234 
3235 template <int dim, int spacedim>
3236 inline
3237 bool
3238 FiniteElement<dim,spacedim>::is_primitive (const unsigned int i) const
3239 {
3240  Assert (i < this->dofs_per_cell, ExcIndexRange (i, 0, this->dofs_per_cell));
3241 
3242  // return primitivity of a shape
3243  // function by checking whether it
3244  // has more than one non-zero
3245  // component or not. we could cache
3246  // this value in an array of bools,
3247  // but accessing a bit-vector (as
3248  // std::vector<bool> is) is
3249  // probably more expensive than
3250  // just comparing against 1
3251  //
3252  // for good measure, short circuit the test
3253  // if the entire FE is primitive
3254  return (is_primitive() ||
3255  (n_nonzero_components_table[i] == 1));
3256 }
3257 
3258 
3259 
3260 template <int dim, int spacedim>
3261 inline
3264 {
3265  Assert (cell_dof_index < this->dofs_per_cell,
3266  ExcIndexRange (cell_dof_index, 0, this->dofs_per_cell));
3267 
3268  // just go through the usual cases, taking into account how DoFs
3269  // are enumerated on the reference cell
3270  if (cell_dof_index < this->first_line_index)
3272  else if (cell_dof_index < this->first_quad_index)
3273  return GeometryPrimitive::line;
3274  else if (cell_dof_index < this->first_hex_index)
3275  return GeometryPrimitive::quad;
3276  else
3277  return GeometryPrimitive::hex;
3278 }
3279 
3280 
3281 
3282 DEAL_II_NAMESPACE_CLOSE
3283 
3284 #endif
static ::ExceptionBase & ExcFEHasNoSupportPoints()
const unsigned int first_hex_index
Definition: fe_base.h:271
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe.cc:1214
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe.cc:439
bool restriction_is_additive(const unsigned int index) const
Definition: fe.h:3194
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:217
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2375
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:358
virtual std::size_t memory_consumption() const
Definition: fe.cc:47
const std::vector< Point< dim-1 > > & get_unit_face_support_points() const
Definition: fe.cc:1026
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:863
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
Definition: fe.cc:939
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2426
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:229
FullMatrix< double > interface_constraints
Definition: fe.h:2401
bool operator!=(const FiniteElement< dim, spacedim > &) const
Definition: fe.cc:957
unsigned int n_nonzero_components(const unsigned int i) const
Definition: fe.h:3217
const std::vector< Point< dim-1 > > & get_generalized_face_support_points() const
Definition: fe.cc:1051
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
Definition: fe.cc:56
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
Definition: fe.cc:1100
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:171
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
Definition: fe.cc:847
static ::ExceptionBase & ExcUnitShapeValuesDoNotExist()
virtual Point< dim > unit_support_point(const unsigned int index) const
Definition: fe.cc:1013
bool restriction_is_implemented() const
Definition: fe.cc:695
const std::vector< Point< dim > > & get_generalized_support_points() const
Definition: fe.cc:991
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:183
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:240
#define AssertIndexRange(index, range)
Definition: exceptions.h:1284
virtual ~InternalDataBase()=default
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
Definition: fe.h:2449
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:918
bool isotropic_prolongation_is_implemented() const
Definition: fe.cc:721
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition: fe.h:2534
bool has_generalized_support_points() const
Definition: fe.cc:1004
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
Definition: fe.cc:357
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
Definition: fe.cc:1199
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
Definition: fe.h:2505
size_type block_size(const unsigned int i) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:322
bool is_primitive() const
Definition: fe.h:3228
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
static ::ExceptionBase & ExcInterpolationNotImplemented()
const std::vector< unsigned int > n_nonzero_components_table
Definition: fe.h:2559
static ::ExceptionBase & ExcIndexRange(int arg1, int arg2, int arg3)
bool has_face_support_points() const
Definition: fe.cc:1042
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2413
bool has_support_points() const
Definition: fe.cc:982
unsigned int component_to_block_index(const unsigned int component) const
Definition: fe.cc:344
static ::ExceptionBase & ExcFENotPrimitive()
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const =0
static ::ExceptionBase & ExcMessage(std::string arg1)
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
Definition: fe.cc:794
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:194
const unsigned int first_quad_index
Definition: fe_base.h:266
virtual std::size_t memory_consumption() const
Definition: fe.cc:1179
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:346
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2389
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
Definition: fe.h:3153
unsigned int global_dof_index
Definition: types.h:88
#define Assert(cond, exc)
Definition: exceptions.h:1142
unsigned int element_multiplicity(const unsigned int index) const
Definition: fe.h:3056
UpdateFlags
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
Definition: fe.cc:302
Abstract base class for mapping classes.
Definition: dof_tools.h:46
GeometryPrimitive get_associated_geometry_primitive(const unsigned int cell_dof_index) const
Definition: fe.h:3263
std::vector< Point< dim-1 > > unit_face_support_points
Definition: fe.h:2420
Definition: fe.h:36
const ComponentMask & get_nonzero_components(const unsigned int i) const
Definition: fe.h:3206
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:335
#define DeclException0(Exception0)
Definition: exceptions.h:323
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
Definition: fe.cc:879
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
Definition: fe.cc:773
static ::ExceptionBase & ExcWrongInterfaceMatrixSize(int arg1, int arg2)
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index) const
Definition: fe.h:3088
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
Definition: fe.cc:1242
virtual std::string get_name() const =0
static ::ExceptionBase & ExcProjectionVoid()
bool has_generalized_face_support_points() const
Definition: fe.cc:1066
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:896
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
Definition: fe.cc:525
virtual Point< dim-1 > unit_face_support_point(const unsigned int index) const
Definition: fe.cc:1075
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:98
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
Definition: fe.h:3166
const unsigned int dofs_per_cell
Definition: fe_base.h:295
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
std::vector< Point< dim-1 > > generalized_face_support_points
Definition: fe.h:2432
types::global_dof_index first_block_of_base(const unsigned int b) const
Definition: fe.h:3143
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const bool line_orientation) const
Definition: fe.cc:646
const std::vector< Point< dim > > & get_unit_support_points() const
Definition: fe.cc:966
unsigned int component_to_system_index(const unsigned int component, const unsigned int index) const
Definition: fe.h:3066
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:907
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
Definition: fe.h:3032
const FiniteElement< dim, spacedim > & operator[](const unsigned int fe_index) const
Definition: fe.h:3019
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index) const
Definition: fe.h:3131
static ::ExceptionBase & ExcEmbeddingVoid()
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:160
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index(const unsigned int index) const
Definition: fe.h:3118
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
Definition: fe.h:3176
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
Definition: fe.h:2469
virtual std::unique_ptr< InternalDataBase > get_subface_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim-1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
Definition: fe.cc:1228
virtual bool hp_constraints_are_implemented() const
Definition: fe.cc:785
size_type block_start(const unsigned int i) const
Definition: fe.h:33
bool prolongation_is_implemented() const
Definition: fe.cc:669
const unsigned int first_line_index
Definition: fe_base.h:261
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
Definition: fe.cc:263
const bool cached_primitivity
Definition: fe.h:2566
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:206
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2550
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
Definition: fe.cc:274
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
Definition: fe.h:2480
bool isotropic_restriction_is_implemented() const
Definition: fe.cc:747
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const bool face_orientation=true, const bool face_flip=false, const bool face_rotation=false) const
Definition: fe_system.cc:824
BlockIndices base_to_block_indices
Definition: fe.h:2511
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
Definition: fe.cc:252
unsigned int n_base_elements() const
Definition: fe.h:3046
unsigned int size() const
static const unsigned int space_dimension
Definition: fe.h:640
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
Definition: fe.cc:614
virtual void fill_fe_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellSimilarity::Similarity cell_similarity, const Quadrature< dim > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
std::vector< int > adjust_line_dof_index_for_line_orientation_table
Definition: fe.h:2464
TableIndices< 2 > interface_constraints_size() const
Definition: fe.cc:820
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
virtual ~FiniteElement()=default
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
Definition: fe.cc:1150
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2541
UpdateFlags update_each
Definition: fe.h:702
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
Definition: fe.cc:151
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
Definition: fe.cc:1163
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const
Definition: fe.cc:929
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
Definition: fe.cc:1088
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
Definition: fe.h:2499
virtual std::unique_ptr< InternalDataBase > get_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0