Flow123d
release_2.2.0-914-gf1a3a4f
|
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... | |
const unsigned int | n_dofs () const |
Returns the number of degrees of freedom needed by the finite element. More... | |
virtual double | basis_value (const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const |
Calculates the value of the comp-th component of the i-th raw basis function at the point p on the reference element. More... | |
virtual arma::vec::fixed< dim > | basis_grad (const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const |
Calculates the comp-th component of the gradient of the i-th raw basis function at the point p on the reference element. More... | |
virtual unsigned int | n_components () const |
Returns numer of components of the basis function. More... | |
Dof | dof (unsigned int i) const |
Returns i -th degree of freedom. More... | |
virtual | ~FiniteElement () |
Destructor. More... | |
Protected Member Functions | |
void | init (bool primitive=true, FEType type=FEScalar) |
Clears all internal structures. More... | |
void | setup_components () |
Initialize vectors with information about components of basis functions. More... | |
virtual FEInternalData * | initialize (const Quadrature< dim > &q) |
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 void | compute_node_matrix () |
Initializes the node_matrix for computing the coefficients of the shape functions in the raw basis of functions_space_ . This is done by evaluating the dofs_ for the basis function and by inverting the resulting matrix. More... | |
const bool | is_primitive () const |
Indicates whether the basis functions have one or more nonzero components (scalar FE spaces are always primitive). More... | |
unsigned int | system_to_component_index (unsigned sys_idx) const |
Returns the component index for vector valued finite elements. More... | |
const std::vector< bool > & | get_nonzero_components (unsigned int sys_idx) const |
Returns the mask of nonzero components for given basis function. More... | |
Protected Attributes | |
FEType | type_ |
Type of FiniteElement. More... | |
bool | is_primitive_ |
Primitive FE is using componentwise shape functions, i.e. only one component is nonzero for each shape function. More... | |
std::vector< unsigned int > | component_indices_ |
Indices of nonzero components of shape functions (for primitive FE). More... | |
std::vector< std::vector< bool > > | nonzero_components_ |
Footprints of nonzero components of shape functions. More... | |
arma::mat | node_matrix |
Matrix that determines the coefficients of the raw basis functions from the values at the support points. More... | |
FunctionSpace * | function_space_ |
Function space defining the FE. More... | |
std::vector< Dof > | dofs_ |
Set of degrees of freedom (functionals) defining the FE. More... | |
Friends | |
class | FESystem< dim, spacedim > |
class | FEValuesBase< dim, spacedim > |
class | FEValues< dim, spacedim > |
class | FESideValues< dim, spacedim > |
class | FE_P_disc< dim, spacedim > |
Abstract class for the description of a general finite element on a reference simplex in dim
dimensions.
The finite element is determined by a function_space_
and a set of dofs_
. Further it must be set whether the finite element is_primitive_
, which means that even if the functions in function_space_
are vector-valued, the basis functions have only one nonzero component (this is typical for tensor-product FE, e.g. vector-valued polynomial FE, but not for Raviart-Thomas FE).
Description of dofs:
The reference cell consists of lower dimensional entities (points, lines, triangles). Each dof is associated to one of these entities. 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.
Shape functions:
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 30 of file dofhandler.hh.
FiniteElement< dim, spacedim >::FiniteElement | ( | ) |
Constructor.
Definition at line 63 of file finite_element.cc.
|
virtual |
Destructor.
Definition at line 268 of file finite_element.cc.
|
virtual |
Calculates the comp-th
component of the gradient 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. |
comp | Number of vector component. |
Reimplemented in FESystem< dim, spacedim >.
Definition at line 153 of file finite_element.cc.
|
virtual |
Calculates the value of the comp-th
component 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. |
comp | Number of vector component. |
Reimplemented in FESystem< dim, spacedim >.
Definition at line 143 of file finite_element.cc.
|
inlineprotectedvirtual |
Initializes the node_matrix
for computing the coefficients of the shape functions in the raw basis of functions_space_
. This is done by evaluating the dofs_
for the basis function and by inverting the resulting matrix.
Reimplemented in FESystem< dim, spacedim >.
Definition at line 87 of file finite_element.cc.
|
inline |
Returns i
-th degree of freedom.
Definition at line 316 of file finite_element.hh.
|
inlineprotectedvirtual |
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 FESystem< dim, spacedim >.
Definition at line 193 of file finite_element.cc.
|
inlineprotected |
Returns the mask of nonzero components for given basis function.
sys_idx | Index of basis function. |
Definition at line 388 of file finite_element.hh.
|
protected |
Clears all internal structures.
Definition at line 70 of file finite_element.cc.
|
protectedvirtual |
Calculates the data on the reference cell.
q | Quadrature rule. |
flags | Update flags. |
Reimplemented in FESystem< dim, spacedim >.
Definition at line 100 of file finite_element.cc.
|
inlineprotected |
Indicates whether the basis functions have one or more nonzero components (scalar FE spaces are always primitive).
Definition at line 375 of file finite_element.hh.
|
inlinevirtual |
Returns numer of components of the basis function.
Reimplemented in FESystem< dim, spacedim >.
Definition at line 313 of file finite_element.hh.
|
inline |
Returns the number of degrees of freedom needed by the finite element.
Definition at line 286 of file finite_element.hh.
|
protected |
Initialize vectors with information about components of basis functions.
Definition at line 79 of file finite_element.cc.
|
inlineprotected |
Returns the component index for vector valued finite elements.
sys_idx | Index of shape function. |
Definition at line 381 of file finite_element.hh.
|
inlineprotectedvirtual |
Decides which additional quantities have to be computed for each cell.
flags | Computed update flags. |
Reimplemented in FESystem< dim, spacedim >.
Definition at line 164 of file finite_element.cc.
|
friend |
Definition at line 425 of file finite_element.hh.
|
friend |
Definition at line 424 of file finite_element.hh.
|
friend |
Definition at line 421 of file finite_element.hh.
|
friend |
Definition at line 423 of file finite_element.hh.
|
friend |
Definition at line 422 of file finite_element.hh.
|
protected |
Indices of nonzero components of shape functions (for primitive FE).
Definition at line 403 of file finite_element.hh.
|
protected |
Set of degrees of freedom (functionals) defining the FE.
Definition at line 418 of file finite_element.hh.
|
protected |
Function space defining the FE.
Definition at line 415 of file finite_element.hh.
|
protected |
Primitive FE is using componentwise shape functions, i.e. only one component is nonzero for each shape function.
Definition at line 400 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 412 of file finite_element.hh.
|
protected |
Footprints of nonzero components of shape functions.
Definition at line 406 of file finite_element.hh.
|
protected |
Type of FiniteElement.
Definition at line 394 of file finite_element.hh.