Reference documentation for deal.II version 9.0.0
|
#include <deal.II/fe/fe_poly_tensor.h>
Classes | |
class | InternalData |
Public Member Functions | |
FE_PolyTensor (const unsigned int degree, const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components) | |
virtual UpdateFlags | requires_update_flags (const UpdateFlags update_flags) const |
virtual double | shape_value (const unsigned int i, const Point< dim > &p) const |
virtual double | shape_value_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 1, dim > | shape_grad (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 1, dim > | shape_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 2, dim > | shape_grad_grad (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 2, dim > | shape_grad_grad_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
Public Member Functions inherited from FiniteElement< dim, spacedim > | |
FiniteElement (const FiniteElementData< dim > &fe_data, const std::vector< bool > &restriction_is_additive_flags, const std::vector< ComponentMask > &nonzero_components) | |
FiniteElement (FiniteElement< dim, spacedim > &&)=default | |
FiniteElement (const FiniteElement< dim, spacedim > &)=default | |
virtual | ~FiniteElement ()=default |
std::pair< std::unique_ptr< FiniteElement< dim, spacedim > >, unsigned int > | operator^ (const unsigned int multiplicity) const |
virtual std::unique_ptr< FiniteElement< dim, spacedim > > | clone () const =0 |
virtual std::string | get_name () const =0 |
const FiniteElement< dim, spacedim > & | operator[] (const unsigned int fe_index) const |
virtual bool | operator== (const FiniteElement< dim, spacedim > &fe) const |
bool | operator!= (const FiniteElement< dim, spacedim > &) const |
virtual std::size_t | memory_consumption () const |
virtual Tensor< 3, dim > | shape_3rd_derivative (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 3, dim > | shape_3rd_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual Tensor< 4, dim > | shape_4th_derivative (const unsigned int i, const Point< dim > &p) const |
virtual Tensor< 4, dim > | shape_4th_derivative_component (const unsigned int i, const Point< dim > &p, const unsigned int component) const |
virtual bool | has_support_on_face (const unsigned int shape_index, const unsigned int face_index) const |
virtual const FullMatrix< double > & | get_restriction_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const |
virtual const FullMatrix< double > & | get_prolongation_matrix (const unsigned int child, const RefinementCase< dim > &refinement_case=RefinementCase< dim >::isotropic_refinement) const |
bool | prolongation_is_implemented () const |
bool | isotropic_prolongation_is_implemented () const |
bool | restriction_is_implemented () const |
bool | isotropic_restriction_is_implemented () const |
bool | restriction_is_additive (const unsigned int index) const |
const FullMatrix< double > & | constraints (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const |
bool | constraints_are_implemented (const ::internal::SubfaceCase< dim > &subface_case=::internal::SubfaceCase< dim >::case_isotropic) const |
virtual bool | hp_constraints_are_implemented () const |
virtual void | get_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const |
virtual void | get_face_interpolation_matrix (const FiniteElement< dim, spacedim > &source, FullMatrix< double > &matrix) const |
virtual void | get_subface_interpolation_matrix (const FiniteElement< dim, spacedim > &source, const unsigned int subface, FullMatrix< double > &matrix) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_vertex_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_line_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const |
virtual std::vector< std::pair< unsigned int, unsigned int > > | hp_quad_dof_identities (const FiniteElement< dim, spacedim > &fe_other) const |
virtual FiniteElementDomination::Domination | compare_for_face_domination (const FiniteElement< dim, spacedim > &fe_other) const |
std::pair< unsigned int, unsigned int > | system_to_component_index (const unsigned int index) const |
unsigned int | component_to_system_index (const unsigned int component, const unsigned int index) const |
std::pair< unsigned int, unsigned int > | face_system_to_component_index (const unsigned int index) const |
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 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 |
unsigned int | adjust_line_dof_index_for_line_orientation (const unsigned int index, const bool line_orientation) const |
const ComponentMask & | get_nonzero_components (const unsigned int i) const |
unsigned int | n_nonzero_components (const unsigned int i) const |
bool | is_primitive () const |
bool | is_primitive (const unsigned int i) const |
unsigned int | n_base_elements () const |
virtual const FiniteElement< dim, spacedim > & | base_element (const unsigned int index) const |
unsigned int | element_multiplicity (const unsigned int index) const |
const FiniteElement< dim, spacedim > & | get_sub_fe (const ComponentMask &mask) const |
virtual const FiniteElement< dim, spacedim > & | get_sub_fe (const unsigned int first_component, const unsigned int n_selected_components) const |
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > | system_to_base_index (const unsigned int index) const |
std::pair< std::pair< unsigned int, unsigned int >, unsigned int > | face_system_to_base_index (const unsigned int index) const |
types::global_dof_index | first_block_of_base (const unsigned int b) const |
std::pair< unsigned int, unsigned int > | component_to_base_index (const unsigned int component) const |
std::pair< unsigned int, unsigned int > | block_to_base_index (const unsigned int block) const |
std::pair< unsigned int, types::global_dof_index > | system_to_block_index (const unsigned int component) const |
unsigned int | component_to_block_index (const unsigned int component) const |
ComponentMask | component_mask (const FEValuesExtractors::Scalar &scalar) const |
ComponentMask | component_mask (const FEValuesExtractors::Vector &vector) const |
ComponentMask | component_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const |
ComponentMask | component_mask (const BlockMask &block_mask) const |
BlockMask | block_mask (const FEValuesExtractors::Scalar &scalar) const |
BlockMask | block_mask (const FEValuesExtractors::Vector &vector) const |
BlockMask | block_mask (const FEValuesExtractors::SymmetricTensor< 2 > &sym_tensor) const |
BlockMask | block_mask (const ComponentMask &component_mask) const |
virtual std::pair< Table< 2, bool >, std::vector< unsigned int > > | get_constant_modes () const |
const std::vector< Point< dim > > & | get_unit_support_points () const |
bool | has_support_points () const |
virtual Point< dim > | unit_support_point (const unsigned int index) const |
const std::vector< Point< dim-1 > > & | get_unit_face_support_points () const |
bool | has_face_support_points () const |
virtual Point< dim-1 > | unit_face_support_point (const unsigned int index) const |
const std::vector< Point< dim > > & | get_generalized_support_points () const |
bool | has_generalized_support_points () const |
const std::vector< Point< dim-1 > > & | get_generalized_face_support_points () const |
bool | has_generalized_face_support_points () const |
GeometryPrimitive | get_associated_geometry_primitive (const unsigned int cell_dof_index) 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 |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Public Member Functions inherited from FiniteElementData< dim > | |
FiniteElementData (const std::vector< unsigned int > &dofs_per_object, const unsigned int n_components, const unsigned int degree, const Conformity conformity=unknown, const BlockIndices &block_indices=BlockIndices()) | |
unsigned int | n_dofs_per_vertex () const |
unsigned int | n_dofs_per_line () const |
unsigned int | n_dofs_per_quad () const |
unsigned int | n_dofs_per_hex () const |
unsigned int | n_dofs_per_face () const |
unsigned int | n_dofs_per_cell () const |
template<int structdim> | |
unsigned int | n_dofs_per_object () const |
unsigned int | n_components () const |
unsigned int | n_blocks () const |
const BlockIndices & | block_indices () const |
unsigned int | tensor_degree () const |
bool | conforms (const Conformity) const |
bool | operator== (const FiniteElementData &) const |
Protected Member Functions | |
virtual std::unique_ptr< typename FiniteElement< dim, spacedim >::InternalDataBase > | get_data (const UpdateFlags update_flags, const Mapping< dim, spacedim > &, const Quadrature< dim > &quadrature, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &) const |
Protected Member Functions inherited from FiniteElement< dim, spacedim > | |
void | reinit_restriction_and_prolongation_matrices (const bool isotropic_restriction_only=false, const bool isotropic_prolongation_only=false) |
TableIndices< 2 > | interface_constraints_size () 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 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 |
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 |
virtual void | fill_fe_face_values (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const unsigned int face_no, const Quadrature< dim-1 > &quadrature, const Mapping< dim, spacedim > &mapping, const typename Mapping< dim, spacedim >::InternalDataBase &mapping_internal, const ::internal::FEValuesImplementation::MappingRelatedData< dim, spacedim > &mapping_data, const InternalDataBase &fe_internal, ::internal::FEValuesImplementation::FiniteElementRelatedData< dim, spacedim > &output_data) const =0 |
virtual 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 |
Protected Attributes | |
MappingType | mapping_type |
PolynomialType | poly_space |
FullMatrix< double > | inverse_node_matrix |
Threads::Mutex | cache_mutex |
Point< dim > | cached_point |
std::vector< Tensor< 1, dim > > | cached_values |
std::vector< Tensor< 2, dim > > | cached_grads |
std::vector< Tensor< 3, dim > > | cached_grad_grads |
Protected Attributes inherited from FiniteElement< dim, spacedim > | |
std::vector< std::vector< FullMatrix< double > > > | restriction |
std::vector< std::vector< FullMatrix< double > > > | prolongation |
FullMatrix< double > | interface_constraints |
std::vector< Point< dim > > | unit_support_points |
std::vector< Point< dim-1 > > | unit_face_support_points |
std::vector< Point< dim > > | generalized_support_points |
std::vector< Point< dim-1 > > | generalized_face_support_points |
Table< 2, int > | adjust_quad_dof_index_for_face_orientation_table |
std::vector< int > | adjust_line_dof_index_for_line_orientation_table |
std::vector< std::pair< unsigned int, unsigned int > > | system_to_component_table |
std::vector< std::pair< unsigned int, unsigned int > > | face_system_to_component_table |
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | system_to_base_table |
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | face_system_to_base_table |
BlockIndices | base_to_block_indices |
std::vector< std::pair< std::pair< unsigned int, unsigned int >, unsigned int > > | component_to_base_table |
const std::vector< bool > | restriction_is_additive_flags |
const std::vector< ComponentMask > | nonzero_components |
const std::vector< unsigned int > | n_nonzero_components_table |
const bool | cached_primitivity |
Additional Inherited Members | |
Public Types inherited from FiniteElementData< dim > | |
enum | Conformity { unknown = 0x00, L2 = 0x01, Hcurl = 0x02, Hdiv = 0x04, H1 = Hcurl | Hdiv, H2 = 0x0e } |
Static Public Member Functions inherited from FiniteElement< dim, spacedim > | |
static ::ExceptionBase & | ExcShapeFunctionNotPrimitive (int arg1) |
static ::ExceptionBase & | ExcFENotPrimitive () |
static ::ExceptionBase & | ExcUnitShapeValuesDoNotExist () |
static ::ExceptionBase & | ExcFEHasNoSupportPoints () |
static ::ExceptionBase & | ExcEmbeddingVoid () |
static ::ExceptionBase & | ExcProjectionVoid () |
static ::ExceptionBase & | ExcWrongInterfaceMatrixSize (int arg1, int arg2) |
static ::ExceptionBase & | ExcInterpolationNotImplemented () |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
Public Attributes inherited from FiniteElementData< dim > | |
const unsigned int | dofs_per_vertex |
const unsigned int | dofs_per_line |
const unsigned int | dofs_per_quad |
const unsigned int | dofs_per_hex |
const unsigned int | first_line_index |
const unsigned int | first_quad_index |
const unsigned int | first_hex_index |
const unsigned int | first_face_line_index |
const unsigned int | first_face_quad_index |
const unsigned int | dofs_per_face |
const unsigned int | dofs_per_cell |
const unsigned int | components |
const unsigned int | degree |
const Conformity | conforming_space |
const BlockIndices | block_indices_data |
Static Public Attributes inherited from FiniteElement< dim, spacedim > | |
static const unsigned int | space_dimension = spacedim |
Static Public Attributes inherited from FiniteElementData< dim > | |
static const unsigned int | dimension = dim |
Static Protected Member Functions inherited from FiniteElement< dim, spacedim > | |
static std::vector< unsigned int > | compute_n_nonzero_components (const std::vector< ComponentMask > &nonzero_components) |
This class gives a unified framework for the implementation of FiniteElement classes based on Tensor valued polynomial spaces like PolynomialsBDM and PolynomialsRaviartThomas.
In essence, what this class requires is that a derived class describes to it a (vector-valued) polynomial space in which every polynomial has exactly dim
vector components. The polynomial space is described through the PolynomialType
template argument, which needs to provide a function of the following signature:
For more information on the template parameter spacedim
, see the documentation for the class Triangulation.
This class is not a fully implemented FiniteElement class, but implements some common features of vector valued elements based on vector valued polynomial classes. What's missing here in particular is information on the topological location of the node values, and derived classes need to provide this information.
Similarly, in many cases, node functionals depend on the shape of the mesh cell, since they evaluate normal or tangential components on the faces. In order to allow for a set of transformations, the variable mapping_type has been introduced. It should needs be set in the constructor of a derived class.
Any derived class must decide on the polynomial space to use. This polynomial space should be implemented simply as a set of vector valued polynomials like PolynomialsBDM and PolynomialsRaviartThomas. In order to facilitate this implementation, which basis the polynomial space chooses is not of importance to the current class – as described next, this class handles the transformation from the basis chosen by the polynomial space template argument to the basis we want to use for finite element computations internally.
In most cases, the basis used by the class that describes the polynomial space, \(\{\tilde\varphi_j(\hat{\mathbf x})\}\), does not match the one we want to use for the finite element description, \(\{\varphi_j(\hat{\mathbf x})\}\). Rather, we need to express the finite element shape functions as a linear combination of the basis provided by the polynomial space:
\begin{align*} \varphi_j = \sum_k c_{jk} \tilde\varphi_j. \end{align*}
These expansion coefficients \(c_{jk}\) are typically computed in the constructors of derived classes. To facilitate this, this class at first (unless told otherwise, see below), assumes that the shape functions should be exactly the ones provided by the polynomial space. In the constructor of the derived class, one then typically has code of the form
The FETools::compute_node_matrix() function explains in more detail what exactly it computes, and how; in any case, the result is that inverse_node_matrix
now contains the expansion coefficients \(c_{jk}\), and the fact that this block of code now sets the matrix to a non-zero size indicates to the functions of the current class that it should from then on use the expanded basis, \(\{\varphi_j(\hat{\mathbf x})\}\), and no longer the original, "raw" basis \(\{\tilde\varphi_j(\hat{\mathbf x})\}\) when asked for values or derivatives of shape functions.
In order for this scheme to work, it is important to ensure that the size of the inverse_node_matrix
be zero at the time when FETools::compute_node_matrix() is called; thus, the call to this function cannot be inlined into the last line – the result of the call really does need to be stored in the temporary object M
.
In most cases, vector valued basis functions must be transformed when mapped from the reference cell to the actual grid cell. These transformations can be selected from the set MappingType and stored in mapping_type. Therefore, each constructor should contain a line like:
(in case no mapping is required) or using whatever value among the ones defined in MappingType is appropriate for the element you are implementing.
Definition at line 140 of file fe_poly_tensor.h.
FE_PolyTensor< PolynomialType, dim, spacedim >::FE_PolyTensor | ( | const unsigned int | degree, |
const FiniteElementData< dim > & | fe_data, | ||
const std::vector< bool > & | restriction_is_additive_flags, | ||
const std::vector< ComponentMask > & | nonzero_components | ||
) |
Constructor.
degree:
constructor argument for poly. May be different from fe_data.degree
. Definition at line 145 of file fe_poly_tensor.cc.
|
virtual |
Given a set of update flags, compute which other quantities also need to be computed in order to satisfy the request by the given flags. Then return the combination of the original set of flags and those just computed.
As an example, if update_flags
contains update_gradients a finite element class will typically require the computation of the inverse of the Jacobian matrix in order to rotate the gradient of shape functions on the reference cell to the real cell. It would then return not just update_gradients, but also update_covariant_transformation, the flag that makes the mapping class produce the inverse of the Jacobian matrix.
An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module.
Implements FiniteElement< dim, spacedim >.
Definition at line 1735 of file fe_poly_tensor.cc.
|
virtual |
Compute the (scalar) value of shape function i
at the given quadrature point p
. Since the elements represented by this class are vector valued, there is no such scalar value and the function therefore throws an exception.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 173 of file fe_poly_tensor.cc.
|
virtual |
Just like for shape_value(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the value of the component-th
vector component of the ith
shape function at point p
.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 185 of file fe_poly_tensor.cc.
|
virtual |
Compute the gradient of (scalar) shape function i
at the given quadrature point p
. Since the elements represented by this class are vector valued, there is no such scalar value and the function therefore throws an exception.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 217 of file fe_poly_tensor.cc.
|
virtual |
Just like for shape_grad(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the gradient of the component-th
vector component of the ith
shape function at point p
.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 229 of file fe_poly_tensor.cc.
|
virtual |
Compute the Hessian of (scalar) shape function i
at the given quadrature point p
. Since the elements represented by this class are vector valued, there is no such scalar value and the function therefore throws an exception.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 263 of file fe_poly_tensor.cc.
|
virtual |
Just like for shape_grad_grad(), but this function will be called when the shape function has more than one non-zero vector component. In that case, this function should return the gradient of the component-th
vector component of the ith
shape function at point p
.
Reimplemented from FiniteElement< dim, spacedim >.
Definition at line 275 of file fe_poly_tensor.cc.
|
inlineprotectedvirtual |
Create an internal data object and return a pointer to it of which the caller of this function then assumes ownership. This object will then be passed to the FiniteElement::fill_fe_values() every time the finite element shape functions and their derivatives are evaluated on a concrete cell. The object created here is therefore used by derived classes as a place for scratch objects that are used in evaluating shape functions, as well as to store information that can be pre-computed once and re-used on every cell (e.g., for evaluating the values and gradients of shape functions on the reference cell, for later re-use when transforming these values to a concrete cell).
This function is the first one called in the process of initializing a FEValues object for a given mapping and finite element object. The returned object will later be passed to FiniteElement::fill_fe_values() for a concrete cell, which will itself place its output into an object of type internal::FEValuesImplementation::FiniteElementRelatedData. Since there may be data that can already be computed in its final form on the reference cell, this function also receives a reference to the internal::FEValuesImplementation::FiniteElementRelatedData object as its last argument. This output argument is guaranteed to always be the same one when used with the InternalDataBase object returned by this function. In other words, the subdivision of scratch data and final data in the returned object and the output_data
object is as follows: If data can be pre- computed on the reference cell in the exact form in which it will later be needed on a concrete cell, then this function should already emplace it in the output_data
object. An example are the values of shape functions at quadrature points for the usual Lagrange elements which on a concrete cell are identical to the ones on the reference cell. On the other hand, if some data can be pre-computed to make computations on a concrete cell cheaper, then it should be put into the returned object for later re-use in a derive class's implementation of FiniteElement::fill_fe_values(). An example are the gradients of shape functions on the reference cell for Lagrange elements: to compute the gradients of the shape functions on a concrete cell, one has to multiply the gradients on the reference cell by the inverse of the Jacobian of the mapping; consequently, we cannot already compute the gradients on a concrete cell at the time the current function is called, but we can at least pre-compute the gradients on the reference cell, and store it in the object returned.
An extensive discussion of the interaction between this function and FEValues can be found in the How Mapping, FiniteElement, and FEValues work together documentation module. See also the documentation of the InternalDataBase class.
[in] | update_flags | A set of UpdateFlags values that describe what kind of information the FEValues object requests the finite element to compute. This set of flags may also include information that the finite element can not compute, e.g., flags that pertain to data produced by the mapping. An implementation of this function needs to set up all data fields in the returned object that are necessary to produce the finite- element related data specified by these flags, and may already pre- compute part of this information as discussed above. Elements may want to store these update flags (or a subset of these flags) in InternalDataBase::update_each so they know at the time when FiniteElement::fill_fe_values() is called what they are supposed to compute |
[in] | mapping | A reference to the mapping used for computing values and derivatives of shape functions. |
[in] | quadrature | A reference to the object that describes where the shape functions should be evaluated. |
[out] | output_data | A reference to the object that FEValues will use in conjunction with the object returned here and where an implementation of FiniteElement::fill_fe_values() will place the requested information. This allows the current function to already pre-compute pieces of information that can be computed on the reference cell, as discussed above. FEValues guarantees that this output object and the object returned by the current function will always be used together. |
Implements FiniteElement< dim, spacedim >.
Definition at line 213 of file fe_poly_tensor.h.
|
protected |
The mapping type to be used to map shape functions from the reference cell to the mesh cell.
Definition at line 206 of file fe_poly_tensor.h.
|
protected |
The polynomial space. Its type is given by the template parameter PolynomialType.
Definition at line 420 of file fe_poly_tensor.h.
|
protected |
The inverse of the matrix aij of node values Ni applied to polynomial pj. This matrix is used to convert polynomials in the "raw" basis provided in poly_space to the basis dual to the node functionals on the reference cell.
This object is not filled by FE_PolyTensor, but is a chance for a derived class to allow for reorganization of the basis functions. If it is left empty, the basis in poly_space is used.
Definition at line 433 of file fe_poly_tensor.h.
|
mutableprotected |
A mutex to be used to guard access to the variables below.
Definition at line 438 of file fe_poly_tensor.h.
|
mutableprotected |
If a shape function is computed at a single point, we must compute all of them to apply inverse_node_matrix. In order to avoid too much overhead, we cache the point and the function values for the next evaluation.
Definition at line 445 of file fe_poly_tensor.h.
|
mutableprotected |
Cached shape function values after call to shape_value_component().
Definition at line 450 of file fe_poly_tensor.h.
|
mutableprotected |
Cached shape function gradients after call to shape_grad_component().
Definition at line 455 of file fe_poly_tensor.h.
|
mutableprotected |
Cached second derivatives of shape functions after call to shape_grad_grad_component().
Definition at line 461 of file fe_poly_tensor.h.