Flow123d
release_2.2.0-36-g163dc99
|
Abstract class for the description of a general finite element on a reference simplex in dim
dimensions.
More...
#include <dofhandler.hh>
Public Member Functions | |
FiniteElement () | |
Constructor. More... | |
void | init () |
Clears all internal structures. More... | |
const unsigned int | n_dofs () const |
Returns the number of degrees of freedom needed by the finite element. More... | |
const unsigned int | n_object_dofs (unsigned int object_dim, DofMultiplicity multiplicity) |
Returns the number of single dofs/dof pairs/triples/sextuples that lie on a single geometric entity of the dimension object_dim . More... | |
virtual double | basis_value (const unsigned int i, const arma::vec::fixed< dim > &p) const =0 |
Calculates the value of the i-th raw basis function at the point p on the reference element. More... | |
virtual arma::vec::fixed< dim > | basis_vector (const unsigned int i, const arma::vec::fixed< dim > &p) const =0 |
Variant of basis_value() for vectorial finite elements. More... | |
virtual arma::vec::fixed< dim > | basis_grad (const unsigned int i, const arma::vec::fixed< dim > &p) const =0 |
Calculates the gradient of the i-th raw basis function at the point p on the reference element. More... | |
virtual arma::mat::fixed< dim, dim > | basis_grad_vector (const unsigned int i, const arma::vec::fixed< dim > &p) const =0 |
Variant of basis_grad() for vectorial finite elements. More... | |
virtual void | compute_node_matrix () |
Initializes the node_matrix for computing the coefficients of the raw basis functions from values at support points. More... | |
virtual FEInternalData * | initialize (const Quadrature< dim > &q, UpdateFlags flags) |
Calculates the data on the reference cell. More... | |
virtual UpdateFlags | update_each (UpdateFlags flags) |
Decides which additional quantities have to be computed for each cell. More... | |
virtual void | fill_fe_values (const Quadrature< dim > &q, FEInternalData &data, FEValuesData< dim, spacedim > &fv_data) |
Computes the shape function values and gradients on the actual cell and fills the FEValues structure. More... | |
virtual const unsigned int | polynomial_order () const |
Returns the maximum degree of space of polynomials contained in the finite element space. More... | |
const bool | is_scalar () const |
Indicates whether the finite element function space is scalar or vectorial. More... | |
const std::vector< arma::vec::fixed< dim > > & | get_generalized_support_points () |
Returns either the generalized support points (if they are defined) or the unit support points. More... | |
virtual | ~FiniteElement () |
Destructor. More... | |
Protected Attributes | |
unsigned int | number_of_dofs |
Total number of degrees of freedom at one finite element. More... | |
unsigned int | number_of_single_dofs [dim+1] |
Number of single dofs at one geometrical entity of the given dimension (point, line, triangle, tetrahedron). More... | |
unsigned int | number_of_pairs [dim+1] |
Number of pairs of dofs at one geometrical entity of the given dimension (applicable to lines and triangles). More... | |
unsigned int | number_of_triples [dim+1] |
Number of triples of dofs associated to one triangle. More... | |
unsigned int | number_of_sextuples [dim+1] |
Number of sextuples of dofs associated to one triangle. More... | |
unsigned int | order |
Polynomial order - to be possibly used in hp methods. More... | |
bool | is_scalar_fe |
Indicator of scalar versus vectorial finite element. More... | |
arma::mat | node_matrix |
Matrix that determines the coefficients of the raw basis functions from the values at the support points. More... | |
std::vector< arma::vec::fixed< dim > > | unit_support_points |
Support points for Lagrangean finite elements. More... | |
std::vector< arma::vec::fixed< dim > > | generalized_support_points |
Support points for non-Lagrangean finite elements. More... | |
Abstract class for the description of a general finite element on a reference simplex in dim
dimensions.
Description of dofs:
The reference cell consists of lower dimensional entities (points, lines, triangles). Each dof is associated to one of these entities. This means that if the entity is shared by 2 or more neighbouring cells in the mesh then this dof is shared by the finite elements on all of these cells. If a dof is associated to the cell itself then it is not shared with neighbouring cells. The ordering of nodes in the entity may not be appropriate for the finite elements on the neighbouring cells, hence we need to describe how the order of dofs changes with the relative configuration of the entity with respect to the actual cell. For this reason we define the dof multiplicity which allows to group the dofs as described in DofMultiplicity.
Support points:
Sometimes it is convenient to describe the function space using a basis (called the raw basis) that is different from the set of shape functions for the finite element (the actual basis). For this reason we define the support points which play the role of nodal functionals associated to the particular dofs. To convert between the two bases one can use the node_matrix
, which is constructed by the method compute_node_matrix(). In the case of non-Lagrangean finite elements the dofs are not associated to nodal functionals but e.g. to derivatives or averages. For that reason we distinguish between the unit support points which are uniquely associated to the dofs and the generalized support points that are auxiliary for the calculation of the dof functionals.
Definition at line 29 of file dofhandler.hh.
FiniteElement< dim, spacedim >::FiniteElement | ( | ) |
Constructor.
Definition at line 33 of file finite_element.cc.
|
virtual |
Destructor.
Definition at line 175 of file finite_element.cc.
|
pure virtual |
Calculates the gradient of the i-th
raw basis function at the point p
on the reference element.
The gradient components are relative to the reference cell coordinate system.
i | Number of the basis function. |
p | Point of evaluation. |
Implemented in FE_P_disc< degree, dim, spacedim >, FE_P_disc< 1, 2, 3 >, FE_P_disc< 1, 3, 3 >, FE_P_disc< 1, 1, 3 >, FE_P_disc< 0, dim, 3 >, FE_P< degree, dim, spacedim >, FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
|
pure virtual |
Variant of basis_grad() for vectorial finite elements.
Calculates the gradient of the i-th
vector-valued raw basis function at the point p
on the reference element. The gradient components are relative to the reference cell coordinate system.
i | Number of the basis function. |
p | Point of evaluation. |
Implemented in FE_P_disc< degree, dim, spacedim >, FE_P_disc< 1, 2, 3 >, FE_P_disc< 1, 3, 3 >, FE_P_disc< 1, 1, 3 >, FE_P_disc< 0, dim, 3 >, FE_P< degree, dim, spacedim >, FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
|
pure virtual |
Calculates the value of the i-th
raw basis function at the point p
on the reference element.
i | Number of the basis function. |
p | Point of evaluation. |
Implemented in FE_P_disc< degree, dim, spacedim >, FE_P_disc< 1, 2, 3 >, FE_P_disc< 1, 3, 3 >, FE_P_disc< 1, 1, 3 >, FE_P_disc< 0, dim, 3 >, FE_P< degree, dim, spacedim >, FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
|
pure virtual |
Variant of basis_value() for vectorial finite elements.
Calculates the value i-th
vector-valued raw basis function at the point p
on the reference element.
i | Number of the basis function. |
p | Point of evaluation. |
Implemented in FE_P_disc< degree, dim, spacedim >, FE_P_disc< 1, 2, 3 >, FE_P_disc< 1, 3, 3 >, FE_P_disc< 1, 1, 3 >, FE_P_disc< 0, dim, 3 >, FE_P< degree, dim, spacedim >, FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
|
inlinevirtual |
Initializes the node_matrix
for computing the coefficients of the raw basis functions from values at support points.
The method is implemented for the case of Langrangean finite element. In other cases it may be reimplemented.
Reimplemented in FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
Definition at line 80 of file finite_element.cc.
|
inlinevirtual |
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
q | Quadrature rule. |
data | Precomputed finite element data. |
fv_data | Data to be computed. |
Reimplemented in FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
Definition at line 138 of file finite_element.cc.
const vector< arma::vec::fixed< dim > > & FiniteElement< dim, spacedim >::get_generalized_support_points | ( | ) |
Returns either the generalized support points (if they are defined) or the unit support points.
Definition at line 161 of file finite_element.cc.
void FiniteElement< dim, spacedim >::init | ( | ) |
Clears all internal structures.
Definition at line 39 of file finite_element.cc.
|
virtual |
Calculates the data on the reference cell.
q | Quadrature rule. |
flags | Update flags. |
Reimplemented in FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
Definition at line 95 of file finite_element.cc.
|
inline |
Indicates whether the finite element function space is scalar or vectorial.
Definition at line 269 of file finite_element.hh.
|
inline |
Returns the number of degrees of freedom needed by the finite element.
Definition at line 53 of file finite_element.cc.
|
inline |
Returns the number of single dofs/dof pairs/triples/sextuples that lie on a single geometric entity of the dimension object_dim
.
object_dim | Dimension of the geometric entity. |
multiplicity | Multiplicity of dofs. |
Definition at line 59 of file finite_element.cc.
|
inlinevirtual |
Returns the maximum degree of space of polynomials contained in the finite element space.
For possible use in hp methods.
Definition at line 261 of file finite_element.hh.
|
inlinevirtual |
Decides which additional quantities have to be computed for each cell.
flags | Computed update flags. |
Reimplemented in FE_RT0< dim, spacedim >, and FE_RT0< dim, 3 >.
Definition at line 127 of file finite_element.cc.
|
protected |
Support points for non-Lagrangean finite elements.
In case of non-Lagrangean finite elements the meaning of the support points is different, hence we denote the structure as generalized_support_points
.
Definition at line 347 of file finite_element.hh.
|
protected |
Indicator of scalar versus vectorial finite element.
Definition at line 321 of file finite_element.hh.
|
protected |
Matrix that determines the coefficients of the raw basis functions from the values at the support points.
Definition at line 327 of file finite_element.hh.
|
protected |
Total number of degrees of freedom at one finite element.
Definition at line 289 of file finite_element.hh.
|
protected |
Number of pairs of dofs at one geometrical entity of the given dimension (applicable to lines and triangles).
Definition at line 301 of file finite_element.hh.
|
protected |
Number of sextuples of dofs associated to one triangle.
Definition at line 311 of file finite_element.hh.
|
protected |
Number of single dofs at one geometrical entity of the given dimension (point, line, triangle, tetrahedron).
Definition at line 295 of file finite_element.hh.
|
protected |
Number of triples of dofs associated to one triangle.
Definition at line 306 of file finite_element.hh.
|
protected |
Polynomial order - to be possibly used in hp methods.
Definition at line 316 of file finite_element.hh.
|
protected |
Support points for Lagrangean finite elements.
Support points are points in the reference element where function values determine the dofs. In case of Lagrangean finite elements the dof values are precisely the function values at unit_support_points
.
Definition at line 338 of file finite_element.hh.