Loading [MathJax]/extensions/TeX/AMSsymbols.js
 Reference documentation for deal.II version 9.4.1
\(\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\}}\)
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages
fe.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 2022 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_h
17#define dealii_fe_h
18
19#include <deal.II/base/config.h>
20
23#include <deal.II/fe/fe_data.h>
26#include <deal.II/fe/mapping.h>
27
29#include <deal.II/lac/vector.h>
30
31#include <memory>
32
33
35
36// Forward declarations:
37#ifndef DOXYGEN
38template <int dim, int spacedim>
39class FEValuesBase;
40template <int dim, int spacedim>
41class FEValues;
42template <int dim, int spacedim>
43class FEFaceValues;
44template <int dim, int spacedim>
45class FESubfaceValues;
46namespace NonMatching
47{
48 template <int dim>
49 class FEImmersedSurfaceValues;
50}
51template <int dim, int spacedim>
52class FESystem;
53#endif
54
653template <int dim, int spacedim = dim>
654class FiniteElement : public Subscriptor, public FiniteElementData<dim>
655{
656public:
660 static constexpr unsigned int space_dimension = spacedim;
661
687 {
688 public:
694
698 virtual ~InternalDataBase() = default;
699
704
720
724 virtual std::size_t
726 };
727
728public:
772 const std::vector<bool> & restriction_is_additive_flags,
773 const std::vector<ComponentMask> &nonzero_components);
774
779
784
789 virtual ~FiniteElement() override = default;
790
799 std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>
800 operator^(const unsigned int multiplicity) const;
801
813 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
814 clone() const = 0;
815
826 virtual std::string
827 get_name() const = 0;
828
854 operator[](const unsigned int fe_index) const;
855
881 virtual double
882 shape_value(const unsigned int i, const Point<dim> &p) const;
883
890 virtual double
891 shape_value_component(const unsigned int i,
892 const Point<dim> & p,
893 const unsigned int component) const;
894
916 virtual Tensor<1, dim>
917 shape_grad(const unsigned int i, const Point<dim> &p) const;
918
925 virtual Tensor<1, dim>
926 shape_grad_component(const unsigned int i,
927 const Point<dim> & p,
928 const unsigned int component) const;
929
951 virtual Tensor<2, dim>
952 shape_grad_grad(const unsigned int i, const Point<dim> &p) const;
953
960 virtual Tensor<2, dim>
961 shape_grad_grad_component(const unsigned int i,
962 const Point<dim> & p,
963 const unsigned int component) const;
964
986 virtual Tensor<3, dim>
987 shape_3rd_derivative(const unsigned int i, const Point<dim> &p) const;
988
995 virtual Tensor<3, dim>
996 shape_3rd_derivative_component(const unsigned int i,
997 const Point<dim> & p,
998 const unsigned int component) const;
999
1021 virtual Tensor<4, dim>
1022 shape_4th_derivative(const unsigned int i, const Point<dim> &p) const;
1023
1030 virtual Tensor<4, dim>
1032 const Point<dim> & p,
1033 const unsigned int component) const;
1044 virtual bool
1045 has_support_on_face(const unsigned int shape_index,
1046 const unsigned int face_index) const;
1047
1049
1069 virtual const FullMatrix<double> &
1070 get_restriction_matrix(const unsigned int child,
1071 const RefinementCase<dim> &refinement_case =
1073
1103 virtual const FullMatrix<double> &
1104 get_prolongation_matrix(const unsigned int child,
1105 const RefinementCase<dim> &refinement_case =
1107
1129 bool
1131
1147 bool
1149
1171 bool
1173
1189 bool
1191
1192
1201 bool
1202 restriction_is_additive(const unsigned int index) const;
1203
1215 const FullMatrix<double> &
1216 constraints(const ::internal::SubfaceCase<dim> &subface_case =
1218
1234 bool
1236 const ::internal::SubfaceCase<dim> &subface_case =
1238
1239
1261 virtual bool
1263
1264
1276 virtual void
1278 FullMatrix<double> & matrix) const;
1280
1298 virtual void
1300 FullMatrix<double> & matrix,
1301 const unsigned int face_no = 0) const;
1302
1303
1315 virtual void
1317 const unsigned int subface,
1318 FullMatrix<double> & matrix,
1319 const unsigned int face_no = 0) const;
1321
1322
1343 virtual std::vector<std::pair<unsigned int, unsigned int>>
1345
1350 virtual std::vector<std::pair<unsigned int, unsigned int>>
1352
1357 virtual std::vector<std::pair<unsigned int, unsigned int>>
1359 const unsigned int face_no = 0) const;
1360
1378 const unsigned int codim = 0) const;
1379
1381
1412 virtual bool
1414
1419 bool
1421
1457 std::pair<unsigned int, unsigned int>
1458 system_to_component_index(const unsigned int index) const;
1459
1469 unsigned int
1470 component_to_system_index(const unsigned int component,
1471 const unsigned int index) const;
1472
1482 std::pair<unsigned int, unsigned int>
1483 face_system_to_component_index(const unsigned int index,
1484 const unsigned int face_no = 0) const;
1485
1494 unsigned int
1496 const unsigned int face_no,
1497 const bool face_orientation,
1498 const bool face_flip,
1499 const bool face_rotation) const;
1500
1555 virtual unsigned int
1556 face_to_cell_index(const unsigned int face_dof_index,
1557 const unsigned int face,
1558 const bool face_orientation = true,
1559 const bool face_flip = false,
1560 const bool face_rotation = false) const;
1561
1569 unsigned int
1571 const bool line_orientation) const;
1572
1589 const ComponentMask &
1590 get_nonzero_components(const unsigned int i) const;
1591
1602 unsigned int
1603 n_nonzero_components(const unsigned int i) const;
1604
1613 bool
1615
1626 bool
1627 is_primitive(const unsigned int i) const;
1628
1640 unsigned int
1642
1647 virtual const FiniteElement<dim, spacedim> &
1648 base_element(const unsigned int index) const;
1649
1656 unsigned int
1657 element_multiplicity(const unsigned int index) const;
1658
1732 get_sub_fe(const ComponentMask &mask) const;
1733
1742 virtual const FiniteElement<dim, spacedim> &
1743 get_sub_fe(const unsigned int first_component,
1744 const unsigned int n_selected_components) const;
1745
1768 std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1769 system_to_base_index(const unsigned int index) const;
1770
1779 std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1780 face_system_to_base_index(const unsigned int index,
1781 const unsigned int face_no = 0) const;
1782
1788 first_block_of_base(const unsigned int b) const;
1789
1802 std::pair<unsigned int, unsigned int>
1803 component_to_base_index(const unsigned int component) const;
1804
1805
1810 std::pair<unsigned int, unsigned int>
1811 block_to_base_index(const unsigned int block) const;
1812
1816 std::pair<unsigned int, types::global_dof_index>
1817 system_to_block_index(const unsigned int component) const;
1818
1822 unsigned int
1823 component_to_block_index(const unsigned int component) const;
1824
1826
1846
1861
1877 const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
1878
1895
1916 BlockMask
1918
1935 BlockMask
1937
1955 BlockMask
1957
1980 BlockMask
1982
1998 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
2000
2002
2038 const std::vector<Point<dim>> &
2040
2058 bool
2060
2073 virtual Point<dim>
2074 unit_support_point(const unsigned int index) const;
2075
2102 const std::vector<Point<dim - 1>> &
2103 get_unit_face_support_points(const unsigned int face_no = 0) const;
2104
2113 bool
2114 has_face_support_points(const unsigned int face_no = 0) const;
2115
2120 virtual Point<dim - 1>
2121 unit_face_support_point(const unsigned int index,
2122 const unsigned int face_no = 0) const;
2123
2136 const std::vector<Point<dim>> &
2138
2148 bool
2150
2195 get_associated_geometry_primitive(const unsigned int cell_dof_index) const;
2196
2197
2275 virtual void
2277 const std::vector<Vector<double>> &support_point_values,
2278 std::vector<double> & nodal_values) const;
2279
2281
2290 virtual std::size_t
2292
2299 int,
2300 << "The shape function with index " << arg1
2301 << " is not primitive, i.e. it is vector-valued and "
2302 << "has more than one non-zero vector component. This "
2303 << "function cannot be called for these shape functions. "
2304 << "Maybe you want to use the same function with the "
2305 << "_component suffix?");
2319 "You are trying to access the values or derivatives of shape functions "
2320 "on the reference cell of an element that does not define its shape "
2321 "functions through mapping from the reference cell. Consequently, "
2322 "you cannot ask for shape function values or derivatives on the "
2323 "reference cell.");
2324
2332 "You are trying to access the support points of a finite "
2333 "element that either has no support points at all, or for "
2334 "which the corresponding tables have not been implemented.");
2335
2343 "You are trying to access the matrices that describe how "
2344 "to embed a finite element function on one cell into the "
2345 "finite element space on one of its children (i.e., the "
2346 "'embedding' or 'prolongation' matrices). However, the "
2347 "current finite element can either not define this sort of "
2348 "operation, or it has not yet been implemented.");
2349
2358 "You are trying to access the matrices that describe how "
2359 "to restrict a finite element function from the children "
2360 "of one cell to the finite element space defined on their "
2361 "parent (i.e., the 'restriction' or 'projection' matrices). "
2362 "However, the current finite element can either not define "
2363 "this sort of operation, or it has not yet been "
2364 "implemented.");
2365
2371 int,
2372 int,
2373 << "The interface matrix has a size of " << arg1 << 'x' << arg2
2374 << ", which is not reasonable for the current element "
2375 "in the present dimension.");
2381
2382protected:
2396 void
2398 const bool isotropic_restriction_only = false,
2399 const bool isotropic_prolongation_only = false);
2400
2413 std::vector<std::vector<FullMatrix<double>>> restriction;
2414
2427 std::vector<std::vector<FullMatrix<double>>> prolongation;
2428
2440
2451 std::vector<Point<dim>> unit_support_points;
2452
2458 std::vector<std::vector<Point<dim - 1>>> unit_face_support_points;
2459
2464 std::vector<Point<dim>> generalized_support_points;
2465
2470 std::vector<std::vector<Point<dim - 1>>> generalized_face_support_points;
2471
2488
2503
2507 std::vector<std::pair<unsigned int, unsigned int>> system_to_component_table;
2508
2518 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
2520
2537 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>
2539
2543 std::vector<
2544 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>>
2546
2552
2573 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>
2575
2581 const std::vector<bool> restriction_is_additive_flags;
2582
2590 const std::vector<ComponentMask> nonzero_components;
2591
2599 const std::vector<unsigned int> n_nonzero_components_table;
2600
2607
2621
2627 static std::vector<unsigned int>
2629 const std::vector<ComponentMask> &nonzero_components);
2630
2651 virtual UpdateFlags
2652 requires_update_flags(const UpdateFlags update_flags) const = 0;
2653
2730 virtual std::unique_ptr<InternalDataBase>
2731 get_data(const UpdateFlags update_flags,
2732 const Mapping<dim, spacedim> &mapping,
2733 const Quadrature<dim> & quadrature,
2734 ::internal::FEValuesImplementation::
2735 FiniteElementRelatedData<dim, spacedim> &output_data) const = 0;
2736
2778 virtual std::unique_ptr<InternalDataBase>
2779 get_face_data(const UpdateFlags update_flags,
2780 const Mapping<dim, spacedim> & mapping,
2781 const hp::QCollection<dim - 1> &quadrature,
2782 ::internal::FEValuesImplementation::
2783 FiniteElementRelatedData<dim, spacedim> &output_data) const;
2784
2788 virtual std::unique_ptr<InternalDataBase>
2790 const UpdateFlags update_flags,
2791 const Mapping<dim, spacedim> &mapping,
2792 const Quadrature<dim - 1> & quadrature,
2794 &output_data) const;
2795
2837 virtual std::unique_ptr<InternalDataBase>
2839 const UpdateFlags update_flags,
2840 const Mapping<dim, spacedim> &mapping,
2841 const Quadrature<dim - 1> & quadrature,
2843 spacedim>
2844 &output_data) const;
2845
2926 virtual void
2929 const CellSimilarity::Similarity cell_similarity,
2930 const Quadrature<dim> & quadrature,
2931 const Mapping<dim, spacedim> & mapping,
2932 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
2933 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
2934 spacedim>
2935 & mapping_data,
2936 const InternalDataBase &fe_internal,
2938 spacedim>
2939 &output_data) const = 0;
2940
2983 virtual void
2986 const unsigned int face_no,
2987 const hp::QCollection<dim - 1> & quadrature,
2988 const Mapping<dim, spacedim> & mapping,
2989 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
2990 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
2991 spacedim>
2992 & mapping_data,
2993 const InternalDataBase &fe_internal,
2995 spacedim>
2996 &output_data) const;
2997
3001 virtual void
3004 const unsigned int face_no,
3005 const Quadrature<dim - 1> & quadrature,
3006 const Mapping<dim, spacedim> & mapping,
3007 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3009 & mapping_data,
3010 const InternalDataBase &fe_internal,
3012 &output_data) const;
3013
3059 virtual void
3062 const unsigned int face_no,
3063 const unsigned int sub_no,
3064 const Quadrature<dim - 1> & quadrature,
3065 const Mapping<dim, spacedim> & mapping,
3066 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3067 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
3068 spacedim>
3069 & mapping_data,
3070 const InternalDataBase &fe_internal,
3072 spacedim>
3073 &output_data) const = 0;
3074
3075 friend class InternalDataBase;
3076 friend class FEValuesBase<dim, spacedim>;
3077 friend class FEValues<dim, spacedim>;
3078 friend class FEFaceValues<dim, spacedim>;
3079 friend class FESubfaceValues<dim, spacedim>;
3080 friend class NonMatching::FEImmersedSurfaceValues<dim>;
3081 friend class FESystem<dim, spacedim>;
3082
3083 // explicitly check for sensible template arguments, but not on windows
3084 // because MSVC creates bogus warnings during normal compilation
3085#ifndef DEAL_II_MSVC
3086 static_assert(dim <= spacedim,
3087 "The dimension <dim> of a FiniteElement must be less than or "
3088 "equal to the space dimension <spacedim> in which it lives.");
3089#endif
3090};
3091
3092
3093//----------------------------------------------------------------------//
3094
3095
3096template <int dim, int spacedim>
3097inline const FiniteElement<dim, spacedim> &
3098FiniteElement<dim, spacedim>::operator[](const unsigned int fe_index) const
3099{
3100 (void)fe_index;
3101 Assert(fe_index == 0,
3102 ExcMessage("A fe_index of zero is the only index allowed here"));
3103 return *this;
3104}
3105
3106
3107
3108template <int dim, int spacedim>
3109inline std::pair<unsigned int, unsigned int>
3111 const unsigned int index) const
3112{
3114 Assert(is_primitive(index),
3116 index)));
3118}
3119
3120
3121
3122template <int dim, int spacedim>
3123inline unsigned int
3125{
3126 return base_to_block_indices.size();
3127}
3128
3129
3130
3131template <int dim, int spacedim>
3132inline unsigned int
3134 const unsigned int index) const
3135{
3136 return static_cast<unsigned int>(base_to_block_indices.block_size(index));
3137}
3138
3139
3140
3141template <int dim, int spacedim>
3142inline unsigned int
3144 const unsigned int component,
3145 const unsigned int index) const
3146{
3147 AssertIndexRange(component, this->n_components());
3148 const std::vector<std::pair<unsigned int, unsigned int>>::const_iterator it =
3149 std::find(system_to_component_table.begin(),
3151 std::pair<unsigned int, unsigned int>(component, index));
3152
3154 ExcMessage("You are asking for the number of the shape function "
3155 "within a system element that corresponds to vector "
3156 "component " +
3157 Utilities::int_to_string(component) +
3158 " and within this to "
3159 "index " +
3161 ". But no such "
3162 "shape function exists."));
3163 return std::distance(system_to_component_table.begin(), it);
3164}
3165
3166
3167
3168template <int dim, int spacedim>
3169inline std::pair<unsigned int, unsigned int>
3171 const unsigned int index,
3172 const unsigned int face_no) const
3173{
3175 index,
3176 face_system_to_component_table[this->n_unique_faces() == 1 ? 0 : face_no]
3177 .size());
3178
3179 // in debug mode, check whether the
3180 // function is primitive, since
3181 // otherwise the result may have no
3182 // meaning
3183 //
3184 // since the primitivity tables are
3185 // all geared towards cell dof
3186 // indices, rather than face dof
3187 // indices, we have to work a
3188 // little bit...
3189 //
3190 // in 1d, the face index is equal
3191 // to the cell index
3192 Assert(is_primitive(this->face_to_cell_index(index, face_no)),
3194 index)));
3195
3197 0 :
3198 face_no][index];
3199}
3200
3201
3202
3203template <int dim, int spacedim>
3204inline std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
3206 const unsigned int index) const
3207{
3210}
3211
3212
3213
3214template <int dim, int spacedim>
3215inline std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
3217 const unsigned int index,
3218 const unsigned int face_no) const
3219{
3221 index,
3222 face_system_to_base_table[this->n_unique_faces() == 1 ? 0 : face_no]
3223 .size());
3224 return face_system_to_base_table[this->n_unique_faces() == 1 ? 0 : face_no]
3225 [index];
3226}
3227
3228
3229
3230template <int dim, int spacedim>
3233 const unsigned int index) const
3234{
3235 return base_to_block_indices.block_start(index);
3236}
3237
3238
3239
3240template <int dim, int spacedim>
3241inline std::pair<unsigned int, unsigned int>
3243 const unsigned int index) const
3244{
3246
3247 return component_to_base_table[index].first;
3248}
3249
3250
3251
3252template <int dim, int spacedim>
3253inline std::pair<unsigned int, unsigned int>
3255 const unsigned int index) const
3256{
3258}
3259
3260
3261
3262template <int dim, int spacedim>
3263inline std::pair<unsigned int, types::global_dof_index>
3265 const unsigned int index) const
3266{
3267 AssertIndexRange(index, this->n_dofs_per_cell());
3268 // The block is computed simply as
3269 // first block of this base plus
3270 // the index within the base blocks
3271 return std::pair<unsigned int, types::global_dof_index>(
3273 system_to_base_table[index].first.second,
3274 system_to_base_table[index].second);
3275}
3276
3277
3278
3279template <int dim, int spacedim>
3280inline bool
3282 const unsigned int index) const
3283{
3284 AssertIndexRange(index, this->n_dofs_per_cell());
3285 return restriction_is_additive_flags[index];
3286}
3287
3288
3289
3290template <int dim, int spacedim>
3291inline const ComponentMask &
3293{
3295 return nonzero_components[i];
3296}
3297
3298
3299
3300template <int dim, int spacedim>
3301inline unsigned int
3303{
3306}
3307
3308
3309
3310template <int dim, int spacedim>
3311inline bool
3313{
3314 return cached_primitivity;
3315}
3316
3317
3318
3319template <int dim, int spacedim>
3320inline bool
3321FiniteElement<dim, spacedim>::is_primitive(const unsigned int i) const
3322{
3324
3325 // return primitivity of a shape
3326 // function by checking whether it
3327 // has more than one non-zero
3328 // component or not. we could cache
3329 // this value in an array of bools,
3330 // but accessing a bit-vector (as
3331 // std::vector<bool> is) is
3332 // probably more expensive than
3333 // just comparing against 1
3334 //
3335 // for good measure, short circuit the test
3336 // if the entire FE is primitive
3337 return (is_primitive() || (n_nonzero_components_table[i] == 1));
3338}
3339
3340
3341
3342template <int dim, int spacedim>
3343inline GeometryPrimitive
3345 const unsigned int cell_dof_index) const
3346{
3347 AssertIndexRange(cell_dof_index, this->n_dofs_per_cell());
3348
3349 // just go through the usual cases, taking into account how DoFs
3350 // are enumerated on the reference cell
3351 if (cell_dof_index < this->get_first_line_index())
3353 else if (cell_dof_index < this->get_first_quad_index(0))
3355 else if (cell_dof_index < this->get_first_hex_index())
3357 else
3359}
3360
3361
3362
3364
3365#endif
size_type block_size(const unsigned int i) const
unsigned int size() const
size_type block_start(const unsigned int i) const
std::pair< unsigned int, size_type > global_to_local(const size_type i) const
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 override
unsigned int get_first_line_index() const
unsigned int n_dofs_per_cell() const
unsigned int get_first_quad_index(const unsigned int quad_no=0) const
unsigned int n_components() const
unsigned int n_unique_faces() const
unsigned int get_first_hex_index() const
virtual ~InternalDataBase()=default
InternalDataBase(const InternalDataBase &)=delete
UpdateFlags update_each
Definition: fe.h:719
virtual std::size_t memory_consumption() const
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
bool isotropic_prolongation_is_implemented() const
virtual std::string get_name() const =0
virtual Point< dim > unit_support_point(const unsigned int index) const
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
const std::vector< unsigned int > n_nonzero_components_table
Definition: fe.h:2599
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
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
const ComponentMask & get_nonzero_components(const unsigned int i) const
std::vector< std::vector< Point< dim - 1 > > > unit_face_support_points
Definition: fe.h:2458
bool prolongation_is_implemented() const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
ComponentMask component_mask(const BlockMask &block_mask) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
ComponentMask component_mask(const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
virtual std::unique_ptr< InternalDataBase > get_face_data(const UpdateFlags update_flags, const Mapping< dim, spacedim > &mapping, const hp::QCollection< dim - 1 > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index, const unsigned int face_no=0) const
bool is_primitive(const unsigned int i) const
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
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
bool has_face_support_points(const unsigned int face_no=0) const
const bool cached_primitivity
Definition: fe.h:2606
virtual const FiniteElement< dim, spacedim > & get_sub_fe(const unsigned int first_component, const unsigned int n_selected_components) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other, const unsigned int face_no=0) const
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
bool operator!=(const FiniteElement< dim, spacedim > &) const
ComponentMask component_mask(const FEValuesExtractors::Vector &vector) const
bool has_support_points() const
std::vector< std::vector< FullMatrix< double > > > restriction
Definition: fe.h:2413
bool has_generalized_support_points() const
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const =0
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
std::vector< Table< 2, int > > adjust_quad_dof_index_for_face_orientation_table
Definition: fe.h:2487
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const bool line_orientation) const
const std::vector< Point< dim > > & get_unit_support_points() const
static constexpr unsigned int space_dimension
Definition: fe.h:660
BlockIndices base_to_block_indices
Definition: fe.h:2551
const FiniteElement< dim, spacedim > & operator[](const unsigned int fe_index) const
bool is_primitive() const
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
const std::vector< Point< dim - 1 > > & get_unit_face_support_points(const unsigned int face_no=0) const
bool isotropic_restriction_is_implemented() const
std::vector< std::vector< Point< dim - 1 > > > generalized_face_support_points
Definition: fe.h:2470
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
Definition: fe.h:2502
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
bool restriction_is_additive(const unsigned int index) const
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
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
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
FiniteElement(const FiniteElement< dim, spacedim > &)=default
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix, const unsigned int face_no=0) const
unsigned int component_to_block_index(const unsigned int component) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const hp::QCollection< 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
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
Definition: fe.h:2538
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double > > &support_point_values, std::vector< double > &nodal_values) const
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
BlockMask block_mask(const FEValuesExtractors::Vector &vector) const
BlockMask block_mask(const ComponentMask &component_mask) const
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index(const unsigned int index) const
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
Definition: fe.h:2507
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index, const unsigned int face_no=0) const
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
unsigned int element_multiplicity(const unsigned int index) const
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
const std::vector< bool > restriction_is_additive_flags
Definition: fe.h:2581
virtual ~FiniteElement() override=default
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition: fe.h:2574
TableIndices< 2 > interface_constraints_size() const
std::vector< std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > > face_system_to_base_table
Definition: fe.h:2545
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2451
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
bool restriction_is_implemented() const
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
std::vector< std::vector< std::pair< unsigned int, unsigned int > > > face_system_to_component_table
Definition: fe.h:2519
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
FiniteElement(FiniteElement< dim, spacedim > &&)=default
unsigned int n_nonzero_components(const unsigned int i) const
FullMatrix< double > interface_constraints
Definition: fe.h:2439
const std::vector< Point< dim > > & get_generalized_support_points() const
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
BlockMask block_mask(const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const
unsigned int n_base_elements() const
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix, const unsigned int face_no=0) const
virtual std::size_t memory_consumption() const
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
GeometryPrimitive get_associated_geometry_primitive(const unsigned int cell_dof_index) const
std::vector< Point< dim > > generalized_support_points
Definition: fe.h:2464
unsigned int component_to_system_index(const unsigned int component, const unsigned int index) const
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const unsigned int face_no, const bool face_orientation, const bool face_flip, const bool face_rotation) const
const std::vector< ComponentMask > nonzero_components
Definition: fe.h:2590
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index, const unsigned int face_no=0) const
virtual bool hp_constraints_are_implemented() const
std::vector< std::vector< FullMatrix< double > > > prolongation
Definition: fe.h:2427
types::global_dof_index first_block_of_base(const unsigned int b) const
Abstract base class for mapping classes.
Definition: mapping.h:311
Definition: point.h:111
Definition: tensor.h:503
Definition: vector.h:109
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
UpdateFlags
Point< 2 > first
Definition: grid_out.cc:4603
#define DeclException0(Exception0)
Definition: exceptions.h:464
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1473
static ::ExceptionBase & ExcWrongInterfaceMatrixSize(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:532
static ::ExceptionBase & ExcUnitShapeValuesDoNotExist()
static ::ExceptionBase & ExcProjectionVoid()
#define AssertIndexRange(index, range)
Definition: exceptions.h:1732
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:487
static ::ExceptionBase & ExcFENotPrimitive()
static ::ExceptionBase & ExcEmbeddingVoid()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:509
static ::ExceptionBase & ExcInterpolationNotImplemented()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcFEHasNoSupportPoints()
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:473