Flow123d  release_2.2.0-914-gf1a3a4f
Public Member Functions | Protected Member Functions | Protected Attributes | Friends | List of all members
FiniteElement< dim, spacedim > Class Template Reference

Abstract class for the description of a general finite element on a reference simplex in dim dimensions. More...

#include <dofhandler.hh>

Inheritance diagram for FiniteElement< dim, spacedim >:
Inheritance graph
[legend]
Collaboration diagram for FiniteElement< dim, spacedim >:
Collaboration graph
[legend]

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 FEInternalDatainitialize (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...
 
FunctionSpacefunction_space_
 Function space defining the FE. More...
 
std::vector< Dofdofs_
 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 >
 

Detailed Description

template<unsigned int dim, unsigned int spacedim>
class FiniteElement< 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.

Constructor & Destructor Documentation

template<unsigned int dim, unsigned int spacedim>
FiniteElement< dim, spacedim >::FiniteElement ( )

Constructor.

Definition at line 63 of file finite_element.cc.

template<unsigned int dim, unsigned int spacedim>
FiniteElement< dim, spacedim >::~FiniteElement ( )
virtual

Destructor.

Definition at line 268 of file finite_element.cc.

Member Function Documentation

template<unsigned int dim, unsigned int spacedim>
arma::vec::fixed< dim > FiniteElement< dim, spacedim >::basis_grad ( const unsigned int  i,
const arma::vec::fixed< dim > &  p,
const unsigned int  comp = 0 
) const
virtual

Calculates the comp-th component of the gradient of the i-th raw basis function at the point p on the reference element.

Parameters
iNumber of the basis function.
pPoint of evaluation.
compNumber of vector component.

Reimplemented in FESystem< dim, spacedim >.

Definition at line 153 of file finite_element.cc.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
double FiniteElement< dim, spacedim >::basis_value ( const unsigned int  i,
const arma::vec::fixed< dim > &  p,
const unsigned int  comp = 0 
) const
virtual

Calculates the value of the comp-th component of the i-th raw basis function at the point p on the reference element.

Parameters
iNumber of the basis function.
pPoint of evaluation.
compNumber of vector component.

Reimplemented in FESystem< dim, spacedim >.

Definition at line 143 of file finite_element.cc.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
void FiniteElement< dim, spacedim >::compute_node_matrix ( )
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.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
Dof FiniteElement< dim, spacedim >::dof ( unsigned int  i) const
inline

Returns i -th degree of freedom.

Definition at line 316 of file finite_element.hh.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
void FiniteElement< dim, spacedim >::fill_fe_values ( const Quadrature< dim > &  q,
FEInternalData data,
FEValuesData< dim, spacedim > &  fv_data 
)
inlineprotectedvirtual

Computes the shape function values and gradients on the actual cell and fills the FEValues structure.

Parameters
qQuadrature rule.
dataPrecomputed finite element data.
fv_dataData to be computed.

Reimplemented in FESystem< dim, spacedim >.

Definition at line 193 of file finite_element.cc.

template<unsigned int dim, unsigned int spacedim>
const std::vector<bool>& FiniteElement< dim, spacedim >::get_nonzero_components ( unsigned int  sys_idx) const
inlineprotected

Returns the mask of nonzero components for given basis function.

Parameters
sys_idxIndex of basis function.

Definition at line 388 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
void FiniteElement< dim, spacedim >::init ( bool  primitive = true,
FEType  type = FEScalar 
)
protected

Clears all internal structures.

Definition at line 70 of file finite_element.cc.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
FEInternalData * FiniteElement< dim, spacedim >::initialize ( const Quadrature< dim > &  q)
protectedvirtual

Calculates the data on the reference cell.

Parameters
qQuadrature rule.
flagsUpdate flags.

Reimplemented in FESystem< dim, spacedim >.

Definition at line 100 of file finite_element.cc.

template<unsigned int dim, unsigned int spacedim>
const bool FiniteElement< dim, spacedim >::is_primitive ( ) const
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.

template<unsigned int dim, unsigned int spacedim>
virtual unsigned int FiniteElement< dim, spacedim >::n_components ( ) const
inlinevirtual

Returns numer of components of the basis function.

Reimplemented in FESystem< dim, spacedim >.

Definition at line 313 of file finite_element.hh.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
const unsigned int FiniteElement< dim, spacedim >::n_dofs ( ) const
inline

Returns the number of degrees of freedom needed by the finite element.

Definition at line 286 of file finite_element.hh.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
void FiniteElement< dim, spacedim >::setup_components ( )
protected

Initialize vectors with information about components of basis functions.

Definition at line 79 of file finite_element.cc.

Here is the caller graph for this function:

template<unsigned int dim, unsigned int spacedim>
unsigned int FiniteElement< dim, spacedim >::system_to_component_index ( unsigned  sys_idx) const
inlineprotected

Returns the component index for vector valued finite elements.

Parameters
sys_idxIndex of shape function.

Definition at line 381 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
UpdateFlags FiniteElement< dim, spacedim >::update_each ( UpdateFlags  flags)
inlineprotectedvirtual

Decides which additional quantities have to be computed for each cell.

Parameters
flagsComputed update flags.

Reimplemented in FESystem< dim, spacedim >.

Definition at line 164 of file finite_element.cc.

Friends And Related Function Documentation

template<unsigned int dim, unsigned int spacedim>
friend class FE_P_disc< dim, spacedim >
friend

Definition at line 425 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
friend class FESideValues< dim, spacedim >
friend

Definition at line 424 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
friend class FESystem< dim, spacedim >
friend

Definition at line 421 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
friend class FEValues< dim, spacedim >
friend

Definition at line 423 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
friend class FEValuesBase< dim, spacedim >
friend

Definition at line 422 of file finite_element.hh.

Member Data Documentation

template<unsigned int dim, unsigned int spacedim>
std::vector<unsigned int> FiniteElement< dim, spacedim >::component_indices_
protected

Indices of nonzero components of shape functions (for primitive FE).

Definition at line 403 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
std::vector<Dof> FiniteElement< dim, spacedim >::dofs_
protected

Set of degrees of freedom (functionals) defining the FE.

Definition at line 418 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
FunctionSpace* FiniteElement< dim, spacedim >::function_space_
protected

Function space defining the FE.

Definition at line 415 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
bool FiniteElement< dim, spacedim >::is_primitive_
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.

template<unsigned int dim, unsigned int spacedim>
arma::mat FiniteElement< dim, spacedim >::node_matrix
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.

template<unsigned int dim, unsigned int spacedim>
std::vector<std::vector<bool> > FiniteElement< dim, spacedim >::nonzero_components_
protected

Footprints of nonzero components of shape functions.

Definition at line 406 of file finite_element.hh.

template<unsigned int dim, unsigned int spacedim>
FEType FiniteElement< dim, spacedim >::type_
protected

Type of FiniteElement.

Definition at line 394 of file finite_element.hh.


The documentation for this class was generated from the following files: