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.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 1998 - 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_h
17#define dealii_fe_h
18
19#include <deal.II/base/config.h>
20
23#include <deal.II/fe/fe_base.h>
26#include <deal.II/fe/mapping.h>
27
29
30#include <memory>
31
32
34
35template <int dim, int spacedim>
36class FEValuesBase;
37template <int dim, int spacedim>
38class FEValues;
39template <int dim, int spacedim>
40class FEFaceValues;
41template <int dim, int spacedim>
42class FESubfaceValues;
43template <int dim, int spacedim>
44class FESystem;
45
644template <int dim, int spacedim = dim>
645class FiniteElement : public Subscriptor, public FiniteElementData<dim>
646{
647public:
651 static const unsigned int space_dimension = spacedim;
652
678 {
679 public:
685
689 virtual ~InternalDataBase() = default;
690
695
711
715 virtual std::size_t
717 };
718
719public:
763 const std::vector<bool> & restriction_is_additive_flags,
764 const std::vector<ComponentMask> &nonzero_components);
765
770
775
780 virtual ~FiniteElement() override = default;
781
790 std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>, unsigned int>
791 operator^(const unsigned int multiplicity) const;
792
804 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
805 clone() const = 0;
806
817 virtual std::string
818 get_name() const = 0;
819
845 operator[](const unsigned int fe_index) const;
846
872 virtual double
873 shape_value(const unsigned int i, const Point<dim> &p) const;
874
881 virtual double
882 shape_value_component(const unsigned int i,
883 const Point<dim> & p,
884 const unsigned int component) const;
885
907 virtual Tensor<1, dim>
908 shape_grad(const unsigned int i, const Point<dim> &p) const;
909
916 virtual Tensor<1, dim>
917 shape_grad_component(const unsigned int i,
918 const Point<dim> & p,
919 const unsigned int component) const;
920
942 virtual Tensor<2, dim>
943 shape_grad_grad(const unsigned int i, const Point<dim> &p) const;
944
951 virtual Tensor<2, dim>
952 shape_grad_grad_component(const unsigned int i,
953 const Point<dim> & p,
954 const unsigned int component) const;
955
977 virtual Tensor<3, dim>
978 shape_3rd_derivative(const unsigned int i, const Point<dim> &p) const;
979
986 virtual Tensor<3, dim>
987 shape_3rd_derivative_component(const unsigned int i,
988 const Point<dim> & p,
989 const unsigned int component) const;
990
1012 virtual Tensor<4, dim>
1013 shape_4th_derivative(const unsigned int i, const Point<dim> &p) const;
1014
1021 virtual Tensor<4, dim>
1023 const Point<dim> & p,
1024 const unsigned int component) const;
1035 virtual bool
1036 has_support_on_face(const unsigned int shape_index,
1037 const unsigned int face_index) const;
1038
1040
1060 virtual const FullMatrix<double> &
1061 get_restriction_matrix(const unsigned int child,
1062 const RefinementCase<dim> &refinement_case =
1064
1094 virtual const FullMatrix<double> &
1095 get_prolongation_matrix(const unsigned int child,
1096 const RefinementCase<dim> &refinement_case =
1098
1120 bool
1122
1138 bool
1140
1162 bool
1164
1180 bool
1182
1183
1192 bool
1193 restriction_is_additive(const unsigned int index) const;
1194
1206 const FullMatrix<double> &
1207 constraints(const ::internal::SubfaceCase<dim> &subface_case =
1209
1225 bool
1227 const ::internal::SubfaceCase<dim> &subface_case =
1229
1230
1252 virtual bool
1254
1255
1267 virtual void
1269 FullMatrix<double> & matrix) const;
1271
1289 virtual void
1292 const unsigned int face_no = 0) const;
1293
1294
1306 virtual void
1308 const unsigned int subface,
1310 const unsigned int face_no = 0) const;
1312
1313
1334 virtual std::vector<std::pair<unsigned int, unsigned int>>
1336
1341 virtual std::vector<std::pair<unsigned int, unsigned int>>
1343
1348 virtual std::vector<std::pair<unsigned int, unsigned int>>
1350 const unsigned int face_no = 0) const;
1351
1369 const unsigned int codim = 0) const;
1370
1372
1403 virtual bool
1405
1410 bool
1412
1448 std::pair<unsigned int, unsigned int>
1449 system_to_component_index(const unsigned int index) const;
1450
1460 unsigned int
1461 component_to_system_index(const unsigned int component,
1462 const unsigned int index) const;
1463
1473 std::pair<unsigned int, unsigned int>
1474 face_system_to_component_index(const unsigned int index,
1475 const unsigned int face_no = 0) const;
1476
1485 unsigned int
1487 const unsigned int face_no,
1488 const bool face_orientation,
1489 const bool face_flip,
1490 const bool face_rotation) const;
1491
1546 virtual unsigned int
1547 face_to_cell_index(const unsigned int face_dof_index,
1548 const unsigned int face,
1549 const bool face_orientation = true,
1550 const bool face_flip = false,
1551 const bool face_rotation = false) const;
1552
1560 unsigned int
1562 const bool line_orientation) const;
1563
1580 const ComponentMask &
1581 get_nonzero_components(const unsigned int i) const;
1582
1593 unsigned int
1594 n_nonzero_components(const unsigned int i) const;
1595
1604 bool
1606
1617 bool
1618 is_primitive(const unsigned int i) const;
1619
1631 unsigned int
1633
1638 virtual const FiniteElement<dim, spacedim> &
1639 base_element(const unsigned int index) const;
1640
1647 unsigned int
1648 element_multiplicity(const unsigned int index) const;
1649
1723 get_sub_fe(const ComponentMask &mask) const;
1724
1733 virtual const FiniteElement<dim, spacedim> &
1734 get_sub_fe(const unsigned int first_component,
1735 const unsigned int n_selected_components) const;
1736
1759 std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1760 system_to_base_index(const unsigned int index) const;
1761
1770 std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1771 face_system_to_base_index(const unsigned int index,
1772 const unsigned int face_no = 0) const;
1773
1779 first_block_of_base(const unsigned int b) const;
1780
1793 std::pair<unsigned int, unsigned int>
1794 component_to_base_index(const unsigned int component) const;
1795
1796
1801 std::pair<unsigned int, unsigned int>
1802 block_to_base_index(const unsigned int block) const;
1803
1807 std::pair<unsigned int, types::global_dof_index>
1808 system_to_block_index(const unsigned int component) const;
1809
1813 unsigned int
1814 component_to_block_index(const unsigned int component) const;
1815
1817
1837
1852
1868 const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
1869
1886
1907 BlockMask
1909
1926 BlockMask
1928
1946 BlockMask
1948
1971 BlockMask
1973
1989 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
1991
1993
2029 const std::vector<Point<dim>> &
2031
2049 bool
2051
2064 virtual Point<dim>
2065 unit_support_point(const unsigned int index) const;
2066
2093 const std::vector<Point<dim - 1>> &
2094 get_unit_face_support_points(const unsigned int face_no = 0) const;
2095
2104 bool
2105 has_face_support_points(const unsigned int face_no = 0) const;
2106
2111 virtual Point<dim - 1>
2112 unit_face_support_point(const unsigned int index,
2113 const unsigned int face_no = 0) const;
2114
2127 const std::vector<Point<dim>> &
2129
2139 bool
2141
2186 get_associated_geometry_primitive(const unsigned int cell_dof_index) const;
2187
2188
2266 virtual void
2268 const std::vector<Vector<double>> &support_point_values,
2269 std::vector<double> & nodal_values) const;
2270
2272
2281 virtual std::size_t
2283
2290 int,
2291 << "The shape function with index " << arg1
2292 << " is not primitive, i.e. it is vector-valued and "
2293 << "has more than one non-zero vector component. This "
2294 << "function cannot be called for these shape functions. "
2295 << "Maybe you want to use the same function with the "
2296 << "_component suffix?");
2310 "You are trying to access the values or derivatives of shape functions "
2311 "on the reference cell of an element that does not define its shape "
2312 "functions through mapping from the reference cell. Consequently, "
2313 "you cannot ask for shape function values or derivatives on the "
2314 "reference cell.");
2315
2323 "You are trying to access the support points of a finite "
2324 "element that either has no support points at all, or for "
2325 "which the corresponding tables have not been implemented.");
2326
2334 "You are trying to access the matrices that describe how "
2335 "to embed a finite element function on one cell into the "
2336 "finite element space on one of its children (i.e., the "
2337 "'embedding' or 'prolongation' matrices). However, the "
2338 "current finite element can either not define this sort of "
2339 "operation, or it has not yet been implemented.");
2340
2349 "You are trying to access the matrices that describe how "
2350 "to restrict a finite element function from the children "
2351 "of one cell to the finite element space defined on their "
2352 "parent (i.e., the 'restriction' or 'projection' matrices). "
2353 "However, the current finite element can either not define "
2354 "this sort of operation, or it has not yet been "
2355 "implemented.");
2356
2362 int,
2363 int,
2364 << "The interface matrix has a size of " << arg1 << "x" << arg2
2365 << ", which is not reasonable for the current element "
2366 "in the present dimension.");
2372
2373protected:
2387 void
2389 const bool isotropic_restriction_only = false,
2390 const bool isotropic_prolongation_only = false);
2391
2404 std::vector<std::vector<FullMatrix<double>>> restriction;
2405
2418 std::vector<std::vector<FullMatrix<double>>> prolongation;
2419
2431
2442 std::vector<Point<dim>> unit_support_points;
2443
2449 std::vector<std::vector<Point<dim - 1>>> unit_face_support_points;
2450
2455 std::vector<Point<dim>> generalized_support_points;
2456
2461 std::vector<std::vector<Point<dim - 1>>> generalized_face_support_points;
2462
2479
2494
2498 std::vector<std::pair<unsigned int, unsigned int>> system_to_component_table;
2499
2509 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
2511
2528 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>
2530
2534 std::vector<
2535 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>>
2537
2543
2564 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>
2566
2572 const std::vector<bool> restriction_is_additive_flags;
2573
2581 const std::vector<ComponentMask> nonzero_components;
2582
2590 const std::vector<unsigned int> n_nonzero_components_table;
2591
2598
2612
2618 static std::vector<unsigned int>
2620 const std::vector<ComponentMask> &nonzero_components);
2621
2642 virtual UpdateFlags
2643 requires_update_flags(const UpdateFlags update_flags) const = 0;
2644
2721 virtual std::unique_ptr<InternalDataBase>
2722 get_data(const UpdateFlags update_flags,
2723 const Mapping<dim, spacedim> &mapping,
2724 const Quadrature<dim> & quadrature,
2725 ::internal::FEValuesImplementation::
2726 FiniteElementRelatedData<dim, spacedim> &output_data) const = 0;
2727
2769 virtual std::unique_ptr<InternalDataBase>
2770 get_face_data(const UpdateFlags update_flags,
2771 const Mapping<dim, spacedim> & mapping,
2772 const hp::QCollection<dim - 1> &quadrature,
2773 ::internal::FEValuesImplementation::
2774 FiniteElementRelatedData<dim, spacedim> &output_data) const;
2775
2779 virtual std::unique_ptr<InternalDataBase>
2781 const UpdateFlags update_flags,
2782 const Mapping<dim, spacedim> &mapping,
2783 const Quadrature<dim - 1> & quadrature,
2785 &output_data) const;
2786
2828 virtual std::unique_ptr<InternalDataBase>
2830 const UpdateFlags update_flags,
2831 const Mapping<dim, spacedim> &mapping,
2832 const Quadrature<dim - 1> & quadrature,
2834 spacedim>
2835 &output_data) const;
2836
2917 virtual void
2920 const CellSimilarity::Similarity cell_similarity,
2921 const Quadrature<dim> & quadrature,
2922 const Mapping<dim, spacedim> & mapping,
2923 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
2924 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
2925 spacedim>
2926 & mapping_data,
2927 const InternalDataBase &fe_internal,
2929 spacedim>
2930 &output_data) const = 0;
2931
2974 virtual void
2977 const unsigned int face_no,
2978 const hp::QCollection<dim - 1> & quadrature,
2979 const Mapping<dim, spacedim> & mapping,
2980 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
2981 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
2982 spacedim>
2983 & mapping_data,
2984 const InternalDataBase &fe_internal,
2986 spacedim>
2987 &output_data) const;
2988
2992 virtual void
2995 const unsigned int face_no,
2996 const Quadrature<dim - 1> & quadrature,
2997 const Mapping<dim, spacedim> & mapping,
2998 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3000 & mapping_data,
3001 const InternalDataBase &fe_internal,
3003 &output_data) const;
3004
3050 virtual void
3053 const unsigned int face_no,
3054 const unsigned int sub_no,
3055 const Quadrature<dim - 1> & quadrature,
3056 const Mapping<dim, spacedim> & mapping,
3057 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3058 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
3059 spacedim>
3060 & mapping_data,
3061 const InternalDataBase &fe_internal,
3063 spacedim>
3064 &output_data) const = 0;
3065
3066 friend class InternalDataBase;
3067 friend class FEValuesBase<dim, spacedim>;
3068 friend class FEValues<dim, spacedim>;
3069 friend class FEFaceValues<dim, spacedim>;
3070 friend class FESubfaceValues<dim, spacedim>;
3071 friend class FESystem<dim, spacedim>;
3072
3073 // explicitly check for sensible template arguments, but not on windows
3074 // because MSVC creates bogus warnings during normal compilation
3075#ifndef DEAL_II_MSVC
3076 static_assert(dim <= spacedim,
3077 "The dimension <dim> of a FiniteElement must be less than or "
3078 "equal to the space dimension <spacedim> in which it lives.");
3079#endif
3080};
3081
3082
3083//----------------------------------------------------------------------//
3084
3085
3086template <int dim, int spacedim>
3088 operator[](const unsigned int fe_index) const
3089{
3090 (void)fe_index;
3091 Assert(fe_index == 0,
3092 ExcMessage("A fe_index of zero is the only index allowed here"));
3093 return *this;
3094}
3095
3096
3097
3098template <int dim, int spacedim>
3099inline std::pair<unsigned int, unsigned int>
3101 const unsigned int index) const
3102{
3104 Assert(is_primitive(index),
3106 index)));
3107 return system_to_component_table[index];
3108}
3109
3110
3111
3112template <int dim, int spacedim>
3113inline unsigned int
3115{
3116 return base_to_block_indices.size();
3117}
3118
3119
3120
3121template <int dim, int spacedim>
3122inline unsigned int
3124 const unsigned int index) const
3125{
3126 return static_cast<unsigned int>(base_to_block_indices.block_size(index));
3127}
3128
3129
3130
3131template <int dim, int spacedim>
3132inline unsigned int
3134 const unsigned int component,
3135 const unsigned int index) const
3136{
3137 AssertIndexRange(component, this->n_components());
3138 const std::vector<std::pair<unsigned int, unsigned int>>::const_iterator it =
3139 std::find(system_to_component_table.begin(),
3141 std::pair<unsigned int, unsigned int>(component, index));
3142
3144 ExcMessage("You are asking for the number of the shape function "
3145 "within a system element that corresponds to vector "
3146 "component " +
3147 Utilities::int_to_string(component) +
3148 " and within this to "
3149 "index " +
3151 ". But no such "
3152 "shape function exists."));
3153 return std::distance(system_to_component_table.begin(), it);
3154}
3155
3156
3157
3158template <int dim, int spacedim>
3159inline std::pair<unsigned int, unsigned int>
3161 const unsigned int index,
3162 const unsigned int face_no) const
3163{
3165 index,
3166 face_system_to_component_table[this->n_unique_faces() == 1 ? 0 : face_no]
3167 .size());
3168
3169 // in debug mode, check whether the
3170 // function is primitive, since
3171 // otherwise the result may have no
3172 // meaning
3173 //
3174 // since the primitivity tables are
3175 // all geared towards cell dof
3176 // indices, rather than face dof
3177 // indices, we have to work a
3178 // little bit...
3179 //
3180 // in 1d, the face index is equal
3181 // to the cell index
3182 Assert(is_primitive(this->face_to_cell_index(index, face_no)),
3184 index)));
3185
3187 0 :
3188 face_no][index];
3189}
3190
3191
3192
3193template <int dim, int spacedim>
3194inline std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
3196 const unsigned int index) const
3197{
3199 return system_to_base_table[index];
3200}
3201
3202
3203
3204template <int dim, int spacedim>
3205inline std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
3207 const unsigned int index,
3208 const unsigned int face_no) const
3209{
3211 index,
3212 face_system_to_base_table[this->n_unique_faces() == 1 ? 0 : face_no]
3213 .size());
3214 return face_system_to_base_table[this->n_unique_faces() == 1 ? 0 : face_no]
3215 [index];
3216}
3217
3218
3219
3220template <int dim, int spacedim>
3223 const unsigned int index) const
3224{
3225 return base_to_block_indices.block_start(index);
3226}
3227
3228
3229
3230template <int dim, int spacedim>
3231inline std::pair<unsigned int, unsigned int>
3233 const unsigned int index) const
3234{
3236
3237 return component_to_base_table[index].first;
3238}
3239
3240
3241
3242template <int dim, int spacedim>
3243inline std::pair<unsigned int, unsigned int>
3245 const unsigned int index) const
3246{
3248}
3249
3250
3251
3252template <int dim, int spacedim>
3253inline std::pair<unsigned int, types::global_dof_index>
3255 const unsigned int index) const
3256{
3257 AssertIndexRange(index, this->n_dofs_per_cell());
3258 // The block is computed simply as
3259 // first block of this base plus
3260 // the index within the base blocks
3261 return std::pair<unsigned int, types::global_dof_index>(
3263 system_to_base_table[index].first.second,
3264 system_to_base_table[index].second);
3265}
3266
3267
3268
3269template <int dim, int spacedim>
3270inline bool
3272 const unsigned int index) const
3273{
3274 AssertIndexRange(index, this->n_dofs_per_cell());
3275 return restriction_is_additive_flags[index];
3276}
3277
3278
3279
3280template <int dim, int spacedim>
3281inline const ComponentMask &
3283{
3285 return nonzero_components[i];
3286}
3287
3288
3289
3290template <int dim, int spacedim>
3291inline unsigned int
3293{
3296}
3297
3298
3299
3300template <int dim, int spacedim>
3301inline bool
3303{
3304 return cached_primitivity;
3305}
3306
3307
3308
3309template <int dim, int spacedim>
3310inline bool
3311FiniteElement<dim, spacedim>::is_primitive(const unsigned int i) const
3312{
3314
3315 // return primitivity of a shape
3316 // function by checking whether it
3317 // has more than one non-zero
3318 // component or not. we could cache
3319 // this value in an array of bools,
3320 // but accessing a bit-vector (as
3321 // std::vector<bool> is) is
3322 // probably more expensive than
3323 // just comparing against 1
3324 //
3325 // for good measure, short circuit the test
3326 // if the entire FE is primitive
3327 return (is_primitive() || (n_nonzero_components_table[i] == 1));
3328}
3329
3330
3331
3332template <int dim, int spacedim>
3333inline GeometryPrimitive
3335 const unsigned int cell_dof_index) const
3336{
3337 AssertIndexRange(cell_dof_index, this->n_dofs_per_cell());
3338
3339 // just go through the usual cases, taking into account how DoFs
3340 // are enumerated on the reference cell
3341 if (cell_dof_index < this->get_first_line_index())
3343 else if (cell_dof_index < this->get_first_quad_index(0))
3345 else if (cell_dof_index < this->get_first_hex_index())
3347 else
3349}
3350
3351
3352
3354
3355#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
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:710
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:2590
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:2449
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:2597
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:2404
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:2478
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
static const unsigned int space_dimension
Definition: fe.h:651
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
BlockIndices base_to_block_indices
Definition: fe.h:2542
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:2461
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:2493
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:2529
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:2498
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:2572
virtual ~FiniteElement() override=default
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition: fe.h:2565
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:2536
std::vector< Point< dim > > unit_support_points
Definition: fe.h:2442
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:2510
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:2430
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:2455
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:2581
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:2418
types::global_dof_index first_block_of_base(const unsigned int b) const
Definition: point.h:111
Definition: vector.h:110
#define DEAL_II_DEPRECATED
Definition: config.h:162
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:402
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:403
UpdateFlags
Point< 2 > first
Definition: grid_out.cc:4587
#define DeclException0(Exception0)
Definition: exceptions.h:470
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
#define Assert(cond, exc)
Definition: exceptions.h:1465
static ::ExceptionBase & ExcWrongInterfaceMatrixSize(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition: exceptions.h:538
static ::ExceptionBase & ExcUnitShapeValuesDoNotExist()
static ::ExceptionBase & ExcProjectionVoid()
#define AssertIndexRange(index, range)
Definition: exceptions.h:1690
#define DeclExceptionMsg(Exception, defaulttext)
Definition: exceptions.h:493
static ::ExceptionBase & ExcFENotPrimitive()
static ::ExceptionBase & ExcEmbeddingVoid()
#define DeclException1(Exception1, type1, outsequence)
Definition: exceptions.h:515
static ::ExceptionBase & ExcInterpolationNotImplemented()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcFEHasNoSupportPoints()
@ matrix
Contents is actually a matrix.
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition: utilities.cc:473