deal.II version GIT relicensing-1927-g3de9220933 2024-10-03 08:40:00+00:00
\(\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\}}\)
Loading...
Searching...
No Matches
fe.h
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 1998 - 2024 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#ifndef dealii_fe_h
16#define dealii_fe_h
17
18#include <deal.II/base/config.h>
19
22#include <deal.II/fe/fe_data.h>
25#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 virtual double
855 shape_value(const unsigned int i, const Point<dim> &p) const;
856
863 virtual double
864 shape_value_component(const unsigned int i,
865 const Point<dim> &p,
866 const unsigned int component) const;
867
889 virtual Tensor<1, dim>
890 shape_grad(const unsigned int i, const Point<dim> &p) const;
891
898 virtual Tensor<1, dim>
899 shape_grad_component(const unsigned int i,
900 const Point<dim> &p,
901 const unsigned int component) const;
902
924 virtual Tensor<2, dim>
925 shape_grad_grad(const unsigned int i, const Point<dim> &p) const;
926
933 virtual Tensor<2, dim>
934 shape_grad_grad_component(const unsigned int i,
935 const Point<dim> &p,
936 const unsigned int component) const;
937
959 virtual Tensor<3, dim>
960 shape_3rd_derivative(const unsigned int i, const Point<dim> &p) const;
961
969 virtual Tensor<3, dim>
970 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>
996 shape_4th_derivative(const unsigned int i, const Point<dim> &p) const;
997
1005 virtual Tensor<4, dim>
1007 const Point<dim> &p,
1008 const unsigned int component) const;
1019 virtual bool
1020 has_support_on_face(const unsigned int shape_index,
1021 const unsigned int face_index) const;
1022
1059 virtual const FullMatrix<double> &
1060 get_restriction_matrix(const unsigned int child,
1061 const RefinementCase<dim> &refinement_case =
1063
1093 virtual const FullMatrix<double> &
1094 get_prolongation_matrix(const unsigned int child,
1095 const RefinementCase<dim> &refinement_case =
1097
1119 bool
1121
1137 bool
1139
1161 bool
1163
1179 bool
1181
1182
1191 bool
1192 restriction_is_additive(const unsigned int index) const;
1193
1205 const FullMatrix<double> &
1206 constraints(const ::internal::SubfaceCase<dim> &subface_case =
1208
1224 bool
1226 const ::internal::SubfaceCase<dim> &subface_case =
1228
1229
1251 virtual bool
1253
1254
1266 virtual void
1268 FullMatrix<double> &matrix) const;
1288 virtual void
1290 FullMatrix<double> &matrix,
1291 const unsigned int face_no = 0) const;
1292
1293
1305 virtual void
1307 const unsigned int subface,
1308 FullMatrix<double> &matrix,
1309 const unsigned int face_no = 0) const;
1333 virtual std::vector<std::pair<unsigned int, unsigned int>>
1335
1340 virtual std::vector<std::pair<unsigned int, unsigned int>>
1342
1347 virtual std::vector<std::pair<unsigned int, unsigned int>>
1349 const unsigned int face_no = 0) const;
1350
1368 const unsigned int codim = 0) const;
1369
1402 virtual bool
1404
1409 bool
1411
1447 std::pair<unsigned int, unsigned int>
1448 system_to_component_index(const unsigned int index) const;
1449
1459 unsigned int
1460 component_to_system_index(const unsigned int component,
1461 const unsigned int index) const;
1462
1472 std::pair<unsigned int, unsigned int>
1473 face_system_to_component_index(const unsigned int index,
1474 const unsigned int face_no = 0) const;
1475
1482 unsigned int
1484 const unsigned int index,
1485 const unsigned int face_no,
1486 const unsigned char combined_orientation) const;
1487
1535 virtual unsigned int
1537 const unsigned int face_dof_index,
1538 const unsigned int face,
1539 const unsigned char combined_orientation =
1541
1551 unsigned int
1553 const unsigned int index,
1554 const unsigned char combined_orientation) const;
1555
1572 const ComponentMask &
1573 get_nonzero_components(const unsigned int i) const;
1574
1585 unsigned int
1586 n_nonzero_components(const unsigned int i) const;
1587
1596 bool
1598
1609 bool
1610 is_primitive(const unsigned int i) const;
1611
1623 unsigned int
1625
1630 virtual const FiniteElement<dim, spacedim> &
1631 base_element(const unsigned int index) const;
1632
1639 unsigned int
1640 element_multiplicity(const unsigned int index) const;
1641
1715 get_sub_fe(const ComponentMask &mask) const;
1716
1725 virtual const FiniteElement<dim, spacedim> &
1726 get_sub_fe(const unsigned int first_component,
1727 const unsigned int n_selected_components) const;
1728
1751 std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1752 system_to_base_index(const unsigned int index) const;
1753
1762 std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
1763 face_system_to_base_index(const unsigned int index,
1764 const unsigned int face_no = 0) const;
1765
1771 first_block_of_base(const unsigned int b) const;
1772
1785 std::pair<unsigned int, unsigned int>
1786 component_to_base_index(const unsigned int component) const;
1787
1788
1793 std::pair<unsigned int, unsigned int>
1794 block_to_base_index(const unsigned int block) const;
1795
1799 std::pair<unsigned int, types::global_dof_index>
1800 system_to_block_index(const unsigned int component) const;
1801
1805 unsigned int
1806 component_to_block_index(const unsigned int component) const;
1807
1829
1844
1860 const FEValuesExtractors::SymmetricTensor<2> &sym_tensor) const;
1861
1878
1899 BlockMask
1901
1918 BlockMask
1920
1938 BlockMask
1940
1963 BlockMask
1965
1981 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
1983
2021 const std::vector<Point<dim>> &
2023
2050 bool
2052
2065 virtual Point<dim>
2066 unit_support_point(const unsigned int index) const;
2067
2094 const std::vector<Point<dim - 1>> &
2095 get_unit_face_support_points(const unsigned int face_no = 0) const;
2096
2108 bool
2109 has_face_support_points(const unsigned int face_no = 0) const;
2110
2115 virtual Point<dim - 1>
2116 unit_face_support_point(const unsigned int index,
2117 const unsigned int face_no = 0) const;
2118
2146 const std::vector<Point<dim>> &
2148
2163 bool
2165
2212 get_associated_geometry_primitive(const unsigned int cell_dof_index) const;
2213
2214
2292 virtual void
2294 const std::vector<Vector<double>> &support_point_values,
2295 std::vector<double> &nodal_values) const;
2296
2307 virtual std::size_t
2309
2316 int,
2317 << "The shape function with index " << arg1
2318 << " is not primitive, i.e. it is vector-valued and "
2319 << "has more than one non-zero vector component. This "
2320 << "function cannot be called for these shape functions. "
2321 << "Maybe you want to use the same function with the "
2322 << "_component suffix?");
2336 "You are trying to access the values or derivatives of shape functions "
2337 "on the reference cell of an element that does not define its shape "
2338 "functions through mapping from the reference cell. Consequently, "
2339 "you cannot ask for shape function values or derivatives on the "
2340 "reference cell.");
2341
2349 "You are trying to access the support points of a finite "
2350 "element that either has no support points at all, or for "
2351 "which the corresponding tables have not been implemented.");
2352
2360 "You are trying to access the matrices that describe how "
2361 "to embed a finite element function on one cell into the "
2362 "finite element space on one of its children (i.e., the "
2363 "'embedding' or 'prolongation' matrices). However, the "
2364 "current finite element can either not define this sort of "
2365 "operation, or it has not yet been implemented.");
2366
2375 "You are trying to access the matrices that describe how "
2376 "to restrict a finite element function from the children "
2377 "of one cell to the finite element space defined on their "
2378 "parent (i.e., the 'restriction' or 'projection' matrices). "
2379 "However, the current finite element can either not define "
2380 "this sort of operation, or it has not yet been "
2381 "implemented.");
2382
2388 int,
2389 int,
2390 << "The interface matrix has a size of " << arg1 << 'x' << arg2
2391 << ", which is not reasonable for the current element "
2392 "in the present dimension.");
2398
2399protected:
2413 void
2415 const bool isotropic_restriction_only = false,
2416 const bool isotropic_prolongation_only = false);
2417
2430 std::vector<std::vector<FullMatrix<double>>> restriction;
2431
2444 std::vector<std::vector<FullMatrix<double>>> prolongation;
2445
2457
2468 std::vector<Point<dim>> unit_support_points;
2469
2475 std::vector<std::vector<Point<dim - 1>>> unit_face_support_points;
2476
2481 std::vector<Point<dim>> generalized_support_points;
2482
2487 std::vector<std::vector<Point<dim - 1>>> generalized_face_support_points;
2488
2505
2518
2522 std::vector<std::pair<unsigned int, unsigned int>> system_to_component_table;
2523
2533 std::vector<std::vector<std::pair<unsigned int, unsigned int>>>
2535
2552 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>
2554
2558 std::vector<
2559 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>>
2561
2567
2588 std::vector<std::pair<std::pair<unsigned int, unsigned int>, unsigned int>>
2590
2596 const std::vector<bool> restriction_is_additive_flags;
2597
2605 const std::vector<ComponentMask> nonzero_components;
2606
2614 const std::vector<unsigned int> n_nonzero_components_table;
2615
2622
2637
2643 static std::vector<unsigned int>
2645 const std::vector<ComponentMask> &nonzero_components);
2646
2667 virtual UpdateFlags
2668 requires_update_flags(const UpdateFlags update_flags) const = 0;
2669
2747 virtual std::unique_ptr<InternalDataBase>
2748 get_data(const UpdateFlags update_flags,
2749 const Mapping<dim, spacedim> &mapping,
2750 const Quadrature<dim> &quadrature,
2751 ::internal::FEValuesImplementation::
2752 FiniteElementRelatedData<dim, spacedim> &output_data) const = 0;
2753
2796 virtual std::unique_ptr<InternalDataBase>
2797 get_face_data(const UpdateFlags update_flags,
2798 const Mapping<dim, spacedim> &mapping,
2799 const hp::QCollection<dim - 1> &quadrature,
2800 ::internal::FEValuesImplementation::
2801 FiniteElementRelatedData<dim, spacedim> &output_data) const;
2802
2806 virtual std::unique_ptr<InternalDataBase>
2808 const UpdateFlags update_flags,
2809 const Mapping<dim, spacedim> &mapping,
2810 const Quadrature<dim - 1> &quadrature,
2812 &output_data) const;
2813
2857 virtual std::unique_ptr<InternalDataBase>
2859 const UpdateFlags update_flags,
2860 const Mapping<dim, spacedim> &mapping,
2861 const Quadrature<dim - 1> &quadrature,
2863 spacedim>
2864 &output_data) const;
2865
2946 virtual void
2949 const CellSimilarity::Similarity cell_similarity,
2950 const Quadrature<dim> &quadrature,
2951 const Mapping<dim, spacedim> &mapping,
2952 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
2954 &mapping_data,
2955 const InternalDataBase &fe_internal,
2957 spacedim>
2958 &output_data) const = 0;
2959
3002 virtual void
3005 const unsigned int face_no,
3006 const hp::QCollection<dim - 1> &quadrature,
3007 const Mapping<dim, spacedim> &mapping,
3008 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3010 &mapping_data,
3011 const InternalDataBase &fe_internal,
3013 spacedim>
3014 &output_data) const;
3015
3019 virtual void
3022 const unsigned int face_no,
3023 const Quadrature<dim - 1> &quadrature,
3024 const Mapping<dim, spacedim> &mapping,
3025 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3027 &mapping_data,
3028 const InternalDataBase &fe_internal,
3030 &output_data) const;
3031
3077 virtual void
3080 const unsigned int face_no,
3081 const unsigned int sub_no,
3082 const Quadrature<dim - 1> &quadrature,
3083 const Mapping<dim, spacedim> &mapping,
3084 const typename Mapping<dim, spacedim>::InternalDataBase &mapping_internal,
3086 &mapping_data,
3087 const InternalDataBase &fe_internal,
3089 spacedim>
3090 &output_data) const = 0;
3091
3092 friend class InternalDataBase;
3093 friend class FEValuesBase<dim, spacedim>;
3094 friend class FEValues<dim, spacedim>;
3095 friend class FEFaceValues<dim, spacedim>;
3096 friend class FESubfaceValues<dim, spacedim>;
3097 friend class NonMatching::FEImmersedSurfaceValues<dim>;
3098 friend class FESystem<dim, spacedim>;
3099
3100 // explicitly check for sensible template arguments, but not on windows
3101 // because MSVC creates bogus warnings during normal compilation
3102#ifndef DEAL_II_MSVC
3103 static_assert(dim <= spacedim,
3104 "The dimension <dim> of a FiniteElement must be less than or "
3105 "equal to the space dimension <spacedim> in which it lives.");
3106#endif
3107};
3108
3109
3110//----------------------------------------------------------------------//
3111#ifndef DOXYGEN
3112
3113template <int dim, int spacedim>
3114inline std::pair<unsigned int, unsigned int>
3116 const unsigned int index) const
3117{
3119 Assert(is_primitive(index),
3121 index)));
3122 return system_to_component_table[index];
3123}
3124
3125
3126
3127template <int dim, int spacedim>
3128inline unsigned int
3130{
3131 return base_to_block_indices.size();
3132}
3133
3134
3135
3136template <int dim, int spacedim>
3137inline unsigned int
3139 const unsigned int index) const
3140{
3141 return static_cast<unsigned int>(base_to_block_indices.block_size(index));
3142}
3143
3144
3145
3146template <int dim, int spacedim>
3147inline unsigned int
3149 const unsigned int component,
3150 const unsigned int index) const
3151{
3152 AssertIndexRange(component, this->n_components());
3153 const std::vector<std::pair<unsigned int, unsigned int>>::const_iterator it =
3154 std::find(system_to_component_table.begin(),
3156 std::pair<unsigned int, unsigned int>(component, index));
3157
3159 ExcMessage("You are asking for the number of the shape function "
3160 "within a system element that corresponds to vector "
3161 "component " +
3162 Utilities::int_to_string(component) +
3163 " and within this to "
3164 "index " +
3166 ". But no such "
3167 "shape function exists."));
3168 return std::distance(system_to_component_table.begin(), it);
3169}
3170
3171
3172
3173template <int dim, int spacedim>
3174inline std::pair<unsigned int, unsigned int>
3176 const unsigned int index,
3177 const unsigned int face_no) const
3178{
3180 index,
3181 face_system_to_component_table[this->n_unique_faces() == 1 ? 0 : face_no]
3182 .size());
3183
3184 // in debug mode, check whether the
3185 // function is primitive, since
3186 // otherwise the result may have no
3187 // meaning
3188 //
3189 // since the primitivity tables are
3190 // all geared towards cell dof
3191 // indices, rather than face dof
3192 // indices, we have to work a
3193 // little bit...
3194 //
3195 // in 1d, the face index is equal
3196 // to the cell index
3197 Assert(is_primitive(this->face_to_cell_index(index, face_no)),
3199 index)));
3200
3202 0 :
3203 face_no][index];
3204}
3205
3206
3207
3208template <int dim, int spacedim>
3209inline std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
3211 const unsigned int index) const
3212{
3215}
3216
3217
3218
3219template <int dim, int spacedim>
3220inline std::pair<std::pair<unsigned int, unsigned int>, unsigned int>
3222 const unsigned int index,
3223 const unsigned int face_no) const
3224{
3226 index,
3227 face_system_to_base_table[this->n_unique_faces() == 1 ? 0 : face_no]
3228 .size());
3229 return face_system_to_base_table[this->n_unique_faces() == 1 ? 0 : face_no]
3230 [index];
3231}
3232
3233
3234
3235template <int dim, int spacedim>
3238 const unsigned int index) const
3239{
3240 return base_to_block_indices.block_start(index);
3241}
3242
3243
3244
3245template <int dim, int spacedim>
3246inline std::pair<unsigned int, unsigned int>
3248 const unsigned int index) const
3249{
3251
3252 return component_to_base_table[index].first;
3253}
3254
3255
3256
3257template <int dim, int spacedim>
3258inline std::pair<unsigned int, unsigned int>
3260 const unsigned int index) const
3261{
3263}
3264
3265
3266
3267template <int dim, int spacedim>
3268inline std::pair<unsigned int, types::global_dof_index>
3270 const unsigned int index) const
3271{
3272 AssertIndexRange(index, this->n_dofs_per_cell());
3273 // The block is computed simply as
3274 // first block of this base plus
3275 // the index within the base blocks
3276 return std::pair<unsigned int, types::global_dof_index>(
3278 system_to_base_table[index].first.second,
3279 system_to_base_table[index].second);
3280}
3281
3282
3283
3284template <int dim, int spacedim>
3285inline bool
3287 const unsigned int index) const
3288{
3289 AssertIndexRange(index, this->n_dofs_per_cell());
3290 return restriction_is_additive_flags[index];
3291}
3292
3293
3294
3295template <int dim, int spacedim>
3296inline const ComponentMask &
3298{
3300 return nonzero_components[i];
3301}
3302
3303
3304
3305template <int dim, int spacedim>
3306inline unsigned int
3308{
3310 return n_nonzero_components_table[i];
3311}
3312
3313
3314
3315template <int dim, int spacedim>
3316inline bool
3318{
3319 return cached_primitivity;
3320}
3321
3322
3323
3324template <int dim, int spacedim>
3325inline bool
3326FiniteElement<dim, spacedim>::is_primitive(const unsigned int i) const
3327{
3329
3330 // return primitivity of a shape
3331 // function by checking whether it
3332 // has more than one non-zero
3333 // component or not. we could cache
3334 // this value in an array of bools,
3335 // but accessing a bit-vector (as
3336 // std::vector<bool> is) is
3337 // probably more expensive than
3338 // just comparing against 1
3339 //
3340 // for good measure, short circuit the test
3341 // if the entire FE is primitive
3342 return (is_primitive() || (n_nonzero_components_table[i] == 1));
3343}
3344
3345
3346
3347template <int dim, int spacedim>
3348inline GeometryPrimitive
3350 const unsigned int cell_dof_index) const
3351{
3352 AssertIndexRange(cell_dof_index, this->n_dofs_per_cell());
3353
3354 // just go through the usual cases, taking into account how DoFs
3355 // are enumerated on the reference cell
3356 if (cell_dof_index < this->get_first_line_index())
3358 else if (cell_dof_index < this->get_first_quad_index(0))
3360 else if (cell_dof_index < this->get_first_hex_index())
3362 else
3364}
3365
3366#endif
3367
3369
3370#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 unsigned char combined_orientation=ReferenceCell::default_combined_face_orientation()) 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
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:2614
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
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< std::vector< Point< dim - 1 > > > unit_face_support_points
Definition fe.h:2475
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:2621
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
virtual unsigned int face_to_cell_index(const unsigned int face_dof_index, const unsigned int face, const unsigned char combined_orientation=ReferenceCell::default_combined_face_orientation()) const
std::vector< std::vector< FullMatrix< double > > > restriction
Definition fe.h:2430
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:2504
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
const std::vector< Point< dim > > & get_unit_support_points() 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
static constexpr unsigned int space_dimension
Definition fe.h:660
BlockIndices base_to_block_indices
Definition fe.h:2566
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:2487
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:2517
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 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
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
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
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:2553
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
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const unsigned int face_no, const unsigned char combined_orientation) const
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
Definition fe.h:2522
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:2596
virtual ~FiniteElement() override=default
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const unsigned char combined_orientation) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
Definition fe.h:2589
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:2560
std::vector< Point< dim > > unit_support_points
Definition fe.h:2468
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:2534
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:2456
const std::vector< Point< dim > > & get_generalized_support_points() const
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:2481
unsigned int component_to_system_index(const unsigned int component, const unsigned int index) const
const std::vector< ComponentMask > nonzero_components
Definition fe.h:2605
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:2444
types::global_dof_index first_block_of_base(const unsigned int b) const
Abstract base class for mapping classes.
Definition mapping.h:318
Definition point.h:111
static constexpr unsigned char default_combined_face_orientation()
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:498
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:499
Point< 2 > first
Definition grid_out.cc:4623
#define DeclException0(Exception0)
Definition exceptions.h:466
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
#define Assert(cond, exc)
static ::ExceptionBase & ExcWrongInterfaceMatrixSize(int arg1, int arg2)
#define DeclException2(Exception2, type1, type2, outsequence)
Definition exceptions.h:534
static ::ExceptionBase & ExcUnitShapeValuesDoNotExist()
static ::ExceptionBase & ExcProjectionVoid()
#define AssertIndexRange(index, range)
#define DeclExceptionMsg(Exception, defaulttext)
Definition exceptions.h:489
static ::ExceptionBase & ExcFENotPrimitive()
static ::ExceptionBase & ExcEmbeddingVoid()
#define DeclException1(Exception1, type1, outsequence)
Definition exceptions.h:511
static ::ExceptionBase & ExcInterpolationNotImplemented()
static ::ExceptionBase & ExcMessage(std::string arg1)
static ::ExceptionBase & ExcFEHasNoSupportPoints()
UpdateFlags
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
Definition utilities.cc:470