|
Reference documentation for deal.II version 9.2.0
|
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\)
\(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\)
\(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\)
\(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Go to the documentation of this file.
35 template <
int dim,
int spacedim>
37 template <
int dim,
int spacedim>
39 template <
int dim,
int spacedim>
41 template <
int dim,
int spacedim>
43 template <
int dim,
int spacedim>
647 template <
int dim,
int spacedim = dim>
795 std::pair<std::unique_ptr<FiniteElement<dim, spacedim>>,
unsigned int>
796 operator^(
const unsigned int multiplicity)
const;
809 virtual std::unique_ptr<FiniteElement<dim, spacedim>>
847 operator[](
const unsigned int fe_index)
const;
886 const unsigned int component)
const;
921 const unsigned int component)
const;
956 const unsigned int component)
const;
991 const unsigned int component)
const;
1026 const unsigned int component)
const;
1039 const unsigned int face_index)
const;
1209 constraints(const ::internal::SubfaceCase<dim> &subface_case =
1229 const ::internal::SubfaceCase<dim> &subface_case =
1309 const unsigned int subface,
1334 virtual std::vector<std::pair<unsigned int, unsigned int>>
1341 virtual std::vector<std::pair<unsigned int, unsigned int>>
1348 virtual std::vector<std::pair<unsigned int, unsigned int>>
1383 const unsigned int codim = 0)
const;
1462 std::pair<unsigned int, unsigned int>
1476 const unsigned int index)
const;
1487 std::pair<unsigned int, unsigned int>
1500 const bool face_orientation,
1501 const bool face_flip,
1502 const bool face_rotation)
const;
1558 virtual unsigned int
1560 const unsigned int face,
1561 const bool face_orientation =
true,
1562 const bool face_flip =
false,
1563 const bool face_rotation =
false)
const;
1574 const bool line_orientation)
const;
1746 get_sub_fe(
const unsigned int first_component,
1747 const unsigned int n_selected_components)
const;
1771 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1782 std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
1804 std::pair<unsigned int, unsigned int>
1812 std::pair<unsigned int, unsigned int>
1818 std::pair<unsigned int, types::global_dof_index>
2000 virtual std::pair<Table<2, bool>, std::vector<unsigned int>>
2040 const std::vector<Point<dim>> &
2104 const std::vector<
Point<dim - 1>> &
2122 virtual Point<dim - 1>
2137 const std::vector<Point<dim>> &
2279 std::vector<double> & nodal_values)
const;
2301 <<
"The shape function with index " << arg1
2302 <<
" is not primitive, i.e. it is vector-valued and "
2303 <<
"has more than one non-zero vector component. This "
2304 <<
"function cannot be called for these shape functions. "
2305 <<
"Maybe you want to use the same function with the "
2306 <<
"_component suffix?");
2320 "You are trying to access the values or derivatives of shape functions "
2321 "on the reference cell of an element that does not define its shape "
2322 "functions through mapping from the reference cell. Consequently, "
2323 "you cannot ask for shape function values or derivatives on the "
2333 "You are trying to access the support points of a finite "
2334 "element that either has no support points at all, or for "
2335 "which the corresponding tables have not been implemented.");
2344 "You are trying to access the matrices that describe how "
2345 "to embed a finite element function on one cell into the "
2346 "finite element space on one of its children (i.e., the "
2347 "'embedding' or 'prolongation' matrices). However, the "
2348 "current finite element can either not define this sort of "
2349 "operation, or it has not yet been implemented.");
2359 "You are trying to access the matrices that describe how "
2360 "to restrict a finite element function from the children "
2361 "of one cell to the finite element space defined on their "
2362 "parent (i.e., the 'restriction' or 'projection' matrices). "
2363 "However, the current finite element can either not define "
2364 "this sort of operation, or it has not yet been "
2374 <<
"The interface matrix has a size of " << arg1 <<
"x" << arg2
2375 <<
", which is not reasonable for the current element "
2376 "in the present dimension.");
2399 const bool isotropic_restriction_only =
false,
2400 const bool isotropic_prolongation_only =
false);
2519 std::vector<std::pair<unsigned int, unsigned int>>
2538 std::vector<std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>>
2544 std::vector<std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>>
2573 std::vector<std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>>
2627 static std::vector<unsigned int>
2730 virtual std::unique_ptr<InternalDataBase>
2735 FiniteElementRelatedData<dim, spacedim> &output_data)
const = 0;
2778 virtual std::unique_ptr<InternalDataBase>
2783 FiniteElementRelatedData<dim, spacedim> &output_data)
const;
2826 virtual std::unique_ptr<InternalDataBase>
2833 &output_data)
const;
2922 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
2928 &output_data)
const = 0;
2975 const unsigned int face_no,
2979 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
2985 &output_data)
const = 0;
3035 const unsigned int face_no,
3036 const unsigned int sub_no,
3040 const ::internal::FEValuesImplementation::MappingRelatedData<dim,
3046 &output_data)
const = 0;
3057 #ifndef DEAL_II_MSVC
3058 static_assert(dim <= spacedim,
3059 "The dimension <dim> of a FiniteElement must be less than or "
3060 "equal to the space dimension <spacedim> in which it lives.");
3068 template <
int dim,
int spacedim>
3074 ExcMessage(
"A fe_index of zero is the only index allowed here"));
3080 template <
int dim,
int spacedim>
3081 inline std::pair<unsigned int, unsigned int>
3083 const unsigned int index)
const
3086 Assert(is_primitive(index),
3089 return system_to_component_table[index];
3094 template <
int dim,
int spacedim>
3098 return base_to_block_indices.size();
3103 template <
int dim,
int spacedim>
3106 const unsigned int index)
const
3108 return static_cast<unsigned int>(base_to_block_indices.block_size(index));
3113 template <
int dim,
int spacedim>
3116 const unsigned int component,
3117 const unsigned int index)
const
3120 const std::vector<std::pair<unsigned int, unsigned int>>::const_iterator it =
3121 std::find(system_to_component_table.begin(),
3122 system_to_component_table.end(),
3123 std::pair<unsigned int, unsigned int>(component, index));
3125 Assert(it != system_to_component_table.end(),
3126 ExcMessage(
"You are asking for the number of the shape function "
3127 "within a system element that corresponds to vector "
3130 " and within this to "
3134 "shape function exists."));
3135 return std::distance(system_to_component_table.begin(), it);
3140 template <
int dim,
int spacedim>
3141 inline std::pair<unsigned int, unsigned int>
3143 const unsigned int index)
const
3164 return face_system_to_component_table[index];
3169 template <
int dim,
int spacedim>
3170 inline std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
3172 const unsigned int index)
const
3175 return system_to_base_table[index];
3180 template <
int dim,
int spacedim>
3181 inline std::pair<std::pair<unsigned int, unsigned int>,
unsigned int>
3183 const unsigned int index)
const
3186 return face_system_to_base_table[index];
3191 template <
int dim,
int spacedim>
3194 const unsigned int index)
const
3196 return base_to_block_indices.block_start(index);
3201 template <
int dim,
int spacedim>
3202 inline std::pair<unsigned int, unsigned int>
3204 const unsigned int index)
const
3208 return component_to_base_table[index].first;
3213 template <
int dim,
int spacedim>
3214 inline std::pair<unsigned int, unsigned int>
3216 const unsigned int index)
const
3218 return base_to_block_indices.global_to_local(index);
3223 template <
int dim,
int spacedim>
3224 inline std::pair<unsigned int, types::global_dof_index>
3226 const unsigned int index)
const
3232 return std::pair<unsigned int, types::global_dof_index>(
3233 first_block_of_base(system_to_base_table[index].
first.first) +
3234 system_to_base_table[index].first.second,
3235 system_to_base_table[index].second);
3240 template <
int dim,
int spacedim>
3243 const unsigned int index)
const
3246 return restriction_is_additive_flags[index];
3251 template <
int dim,
int spacedim>
3256 return nonzero_components[i];
3261 template <
int dim,
int spacedim>
3266 return n_nonzero_components_table[i];
3271 template <
int dim,
int spacedim>
3275 return cached_primitivity;
3280 template <
int dim,
int spacedim>
3298 return (is_primitive() || (n_nonzero_components_table[i] == 1));
3303 template <
int dim,
int spacedim>
3306 const unsigned int cell_dof_index)
const
3312 if (cell_dof_index < this->first_line_index)
3314 else if (cell_dof_index < this->first_quad_index)
3316 else if (cell_dof_index < this->first_hex_index)
virtual Tensor< 2, dim > shape_grad_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
std::vector< int > adjust_line_dof_index_for_line_orientation_table
bool restriction_is_implemented() const
#define DeclExceptionMsg(Exception, defaulttext)
unsigned int adjust_line_dof_index_for_line_orientation(const unsigned int index, const bool line_orientation) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_line_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
unsigned int component_to_block_index(const unsigned int component) const
static ::ExceptionBase & ExcProjectionVoid()
const FiniteElement< dim, spacedim > & operator[](const unsigned int fe_index) const
bool isotropic_prolongation_is_implemented() const
TableIndices< 2 > interface_constraints_size() const
std::vector< std::pair< unsigned int, unsigned int > > system_to_component_table
virtual std::size_t memory_consumption() const
std::string int_to_string(const unsigned int value, const unsigned int digits=numbers::invalid_unsigned_int)
virtual ~InternalDataBase()=default
ComponentMask component_mask(const FEValuesExtractors::Scalar &scalar) const
static std::vector< unsigned int > compute_n_nonzero_components(const std::vector< ComponentMask > &nonzero_components)
unsigned int element_multiplicity(const unsigned int index) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_quad_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > component_to_base_table
const std::vector< unsigned int > n_nonzero_components_table
FullMatrix< double > interface_constraints
bool is_primitive() const
BlockIndices base_to_block_indices
friend class InternalDataBase
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
#define AssertIndexRange(index, range)
static ::ExceptionBase & ExcInterpolationNotImplemented()
unsigned int n_base_elements() const
unsigned int n_components(const DoFHandler< dim, spacedim > &dh)
virtual void fill_fe_subface_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const unsigned int sub_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
static ::ExceptionBase & ExcFENotPrimitive()
FiniteElement(const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components)
bool has_face_support_points() const
bool has_support_points() const
static ::ExceptionBase & ExcFEHasNoSupportPoints()
virtual Tensor< 4, dim > shape_4th_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const ComponentMask & get_nonzero_components(const unsigned int i) const
virtual Tensor< 1, dim > shape_grad_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
bool has_generalized_support_points() const
std::pair< unsigned int, unsigned int > block_to_base_index(const unsigned int block) const
Table< 2, int > adjust_quad_dof_index_for_face_orientation_table
static ::ExceptionBase & ExcWrongInterfaceMatrixSize(int arg1, int arg2)
Abstract base class for mapping classes.
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > system_to_base_index(const unsigned int index) const
virtual std::vector< std::pair< unsigned int, unsigned int > > hp_vertex_dof_identities(const FiniteElement< dim, spacedim > &fe_other) const
virtual bool hp_constraints_are_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
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
static ::ExceptionBase & ExcMessage(std::string arg1)
bool constraints_are_implemented(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
std::pair< unsigned int, unsigned int > system_to_component_index(const unsigned int index) const
virtual Tensor< 1, dim > shape_grad(const unsigned int i, const Point< dim > &p) const
virtual std::size_t memory_consumption() const
virtual Tensor< 3, dim > shape_3rd_derivative_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
virtual const FiniteElement< dim, spacedim > & base_element(const unsigned int index) const
const std::vector< ComponentMask > nonzero_components
#define DEAL_II_NAMESPACE_OPEN
std::vector< std::vector< FullMatrix< double > > > prolongation
SymmetricTensor< 2, dim, Number > b(const Tensor< 2, dim, Number > &F)
@ matrix
Contents is actually a matrix.
const FullMatrix< double > & constraints(const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const
virtual UpdateFlags requires_update_flags(const UpdateFlags update_flags) const =0
#define DEAL_II_DEPRECATED
#define DeclException1(Exception1, type1, outsequence)
std::pair< unsigned int, unsigned int > component_to_base_index(const unsigned int component) const
virtual FiniteElementDomination::Domination compare_for_domination(const FiniteElement< dim, spacedim > &fe_other, const unsigned int codim=0) const
static ::ExceptionBase & ExcShapeFunctionNotPrimitive(int arg1)
virtual Point< dim - 1 > unit_face_support_point(const unsigned int index) const
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > get_constant_modes() const
virtual Tensor< 3, dim > shape_3rd_derivative(const unsigned int i, const Point< dim > &p) const
GeometryPrimitive get_associated_geometry_primitive(const unsigned int cell_dof_index) const
const std::vector< Point< dim > > & get_unit_support_points() const
unsigned int component_to_system_index(const unsigned int component, const unsigned int index) const
bool operator!=(const FiniteElement< dim, spacedim > &) const
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
static ::ExceptionBase & ExcUnitShapeValuesDoNotExist()
virtual bool has_support_on_face(const unsigned int shape_index, const unsigned int face_index) const
virtual double shape_value(const unsigned int i, const Point< dim > &p) const
virtual void fill_fe_face_values(const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim - 1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0
static ::ExceptionBase & ExcEmbeddingVoid()
bool prolongation_is_implemented() const
#define Assert(cond, exc)
BlockMask block_mask(const FEValuesExtractors::Scalar &scalar) const
bool restriction_is_additive(const unsigned int index) const
types::global_dof_index first_block_of_base(const unsigned int b) const
std::vector< Point< dim - 1 > > unit_face_support_points
std::vector< std::vector< FullMatrix< double > > > restriction
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
virtual std::unique_ptr< FiniteElement< dim, spacedim > > clone() const =0
#define DeclException0(Exception0)
virtual const FullMatrix< double > & get_prolongation_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
std::vector< Point< dim > > generalized_support_points
const bool cached_primitivity
virtual Tensor< 2, dim > shape_grad_grad(const unsigned int i, const Point< dim > &p) const
static const unsigned int space_dimension
virtual void get_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual std::string get_name() const =0
const std::vector< Point< dim > > & get_generalized_support_points() const
virtual Point< dim > unit_support_point(const unsigned int index) const
bool isotropic_restriction_is_implemented() const
std::pair< unsigned int, types::global_dof_index > system_to_block_index(const unsigned int component) const
const std::vector< bool > restriction_is_additive_flags
virtual void get_subface_interpolation_matrix(const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const
unsigned int n_nonzero_components(const unsigned int i) const
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > face_system_to_base_index(const unsigned int index) 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
void reinit_restriction_and_prolongation_matrices(const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false)
std::pair< unsigned int, unsigned int > face_system_to_component_index(const unsigned int index) const
#define DEAL_II_NAMESPACE_CLOSE
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > system_to_base_table
const std::vector< Point< dim - 1 > > & get_unit_face_support_points() const
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > operator^(const unsigned int multiplicity) const
virtual ~FiniteElement() override=default
#define DeclException2(Exception2, type1, type2, outsequence)
virtual void get_face_interpolation_matrix(const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const
virtual bool operator==(const FiniteElement< dim, spacedim > &fe) const
virtual const FullMatrix< double > & get_restriction_matrix(const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const
virtual Tensor< 4, dim > shape_4th_derivative(const unsigned int i, const Point< dim > &p) const
std::vector< Point< dim - 1 > > generalized_face_support_points
virtual FiniteElementDomination::Domination compare_for_face_domination(const FiniteElement< dim, spacedim > &fe_other) const final
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > face_system_to_base_table
unsigned int adjust_quad_dof_index_for_face_orientation(const unsigned int index, const bool face_orientation, const bool face_flip, const bool face_rotation) const
virtual void convert_generalized_support_point_values_to_dof_values(const std::vector< Vector< double >> &support_point_values, std::vector< double > &nodal_values) const
std::vector< std::pair< unsigned int, unsigned int > > face_system_to_component_table
virtual double shape_value_component(const unsigned int i, const Point< dim > &p, const unsigned int component) const
const FiniteElement< dim, spacedim > & get_sub_fe(const ComponentMask &mask) const
std::vector< Point< dim > > unit_support_points