Flow123d  master-f44eb46
Public Member Functions | Private Member Functions | Private Attributes | List of all members
FESystem< dim > Class Template Reference

Compound finite element on dim dimensional simplex. More...

#include <fe_system.hh>

Inheritance diagram for FESystem< dim >:
Inheritance graph
[legend]
Collaboration diagram for FESystem< dim >:
Collaboration graph
[legend]

Public Member Functions

 FESystem (std::shared_ptr< FiniteElement< dim > > fe, FEType t)
 Constructor for FEVectorContravariant and FEVectorPiola. More...
 
 FESystem (const std::shared_ptr< FiniteElement< dim > > &fe, FEType t, unsigned int n)
 Constructor for FEVector, FETensor and FEMixedSystem. If t == FEVector then n must be the space dimension into which the reference cell will be mapped and that fe is scalar. If t == FETensor then n must be square of the space dimension into which the reference cell will be mapped and that fe is scalar. If t == FEMixedSystem, then n is the number of components. More...
 
 FESystem (std::vector< std::shared_ptr< FiniteElement< dim > > > fe)
 Constructor. FESystem for mixed elements. More...
 
std::vector< unsigned int > get_scalar_components () const
 
std::vector< unsigned int > get_vector_components () const
 
std::vector< unsigned int > get_tensor_components () const
 
UpdateFlags update_each (UpdateFlags flags) override
 Decides which additional quantities have to be computed for each cell. More...
 
virtual std::vector< arma::vec::fixed< dim+1 > > dof_points () const
 
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe () const
 
std::vector< unsigned int > fe_dofs (unsigned int fe_index)
 Return dof indices belonging to given sub-FE. More...
 
std::shared_ptr< FESystemFunctionSpacefunction_space () const
 Return function space casting to FESystemFunctionSpace type. More...
 
- Public Member Functions inherited from FiniteElement< dim >
 FiniteElement ()
 Constructor. More...
 
unsigned int n_dofs () const
 Returns the number of degrees of freedom needed by the finite element. More...
 
double shape_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 shape function at the point p on the reference element. More...
 
arma::vec::fixed< dim > shape_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 shape function at the point p on the reference element. More...
 
unsigned int n_components () const
 Returns numer of components of the basis function.
More...
 
const Dofdof (unsigned int i) const
 Returns i -th degree of freedom. More...
 
unsigned int n_space_components (unsigned int spacedim)
 Number of components of FE in a mapped space with dimension spacedim. More...
 
virtual ~FiniteElement ()
 Destructor. More...
 

Private Member Functions

void initialize ()
 Initialization of the internal structures from the vector of base FE. More...
 
void compute_node_matrix () override
 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...
 

Private Attributes

std::vector< std::shared_ptr< FiniteElement< dim > > > fe_
 Pointers to base FE objects. More...
 
std::vector< unsigned int > scalar_components_
 
std::vector< unsigned int > vector_components_
 
std::vector< unsigned int > tensor_components_
 

Additional Inherited Members

- Protected Member Functions inherited from FiniteElement< dim >
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...
 
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 inherited from FiniteElement< dim >
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...
 
std::shared_ptr< FunctionSpacefunction_space_
 Function space defining the FE. More...
 
std::vector< Dofdofs_
 Set of degrees of freedom (functionals) defining the FE. More...
 

Detailed Description

template<unsigned int dim>
class FESystem< dim >

Compound finite element on dim dimensional simplex.

This type of FE is used for vector-valued functions and for systems of equations.

Definition at line 101 of file fe_system.hh.

Constructor & Destructor Documentation

◆ FESystem() [1/3]

template<unsigned int dim>
FESystem< dim >::FESystem ( std::shared_ptr< FiniteElement< dim > >  fe,
FEType  t 
)

Constructor for FEVectorContravariant and FEVectorPiola.

Parameters
feBase finite element class (must be scalar).
tType (must be FEVectorContravariant or FEVectorPiola).

Definition at line 102 of file fe_system.cc.

◆ FESystem() [2/3]

template<unsigned int dim>
FESystem< dim >::FESystem ( const std::shared_ptr< FiniteElement< dim > > &  fe,
FEType  t,
unsigned int  n 
)

Constructor for FEVector, FETensor and FEMixedSystem. If t == FEVector then n must be the space dimension into which the reference cell will be mapped and that fe is scalar. If t == FETensor then n must be square of the space dimension into which the reference cell will be mapped and that fe is scalar. If t == FEMixedSystem, then n is the number of components.

Parameters
feBase finite element class (must be scalar if t is FEVector or FETensor).
tType of FESystem (must be either FEVector, FETensor or FEMixedSystem).
nMultiplicity (number of components).

Definition at line 115 of file fe_system.cc.

◆ FESystem() [3/3]

template<unsigned int dim>
FESystem< dim >::FESystem ( std::vector< std::shared_ptr< FiniteElement< dim > > >  fe)

Constructor. FESystem for mixed elements.

Parameters
feBase finite element classes.

Definition at line 129 of file fe_system.cc.

Member Function Documentation

◆ compute_node_matrix()

template<unsigned int dim>
void FESystem< dim >::compute_node_matrix
overrideprivatevirtual

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 from FiniteElement< dim >.

Definition at line 246 of file fe_system.cc.

◆ dof_points()

template<unsigned int dim>
vector< arma::vec::fixed< dim+1 > > FESystem< dim >::dof_points
virtual

Get barycentric coordinates of the points on the reference element associated with the dofs. Used in BDDC for unknown reason.

Reimplemented from FiniteElement< dim >.

Definition at line 233 of file fe_system.cc.

◆ fe()

template<unsigned int dim>
const std::vector<std::shared_ptr<FiniteElement<dim> > >& FESystem< dim >::fe ( ) const
inline

Definition at line 147 of file fe_system.hh.

Here is the caller graph for this function:

◆ fe_dofs()

template<unsigned int dim>
std::vector< unsigned int > FESystem< dim >::fe_dofs ( unsigned int  fe_index)

Return dof indices belonging to given sub-FE.

Definition at line 267 of file fe_system.cc.

Here is the caller graph for this function:

◆ function_space()

template<unsigned int dim>
std::shared_ptr<FESystemFunctionSpace> FESystem< dim >::function_space ( ) const
inline

Return function space casting to FESystemFunctionSpace type.

Definition at line 154 of file fe_system.hh.

◆ get_scalar_components()

template<unsigned int dim>
std::vector<unsigned int> FESystem< dim >::get_scalar_components ( ) const
inline

Definition at line 131 of file fe_system.hh.

Here is the caller graph for this function:

◆ get_tensor_components()

template<unsigned int dim>
std::vector<unsigned int> FESystem< dim >::get_tensor_components ( ) const
inline

Definition at line 137 of file fe_system.hh.

Here is the caller graph for this function:

◆ get_vector_components()

template<unsigned int dim>
std::vector<unsigned int> FESystem< dim >::get_vector_components ( ) const
inline

Definition at line 134 of file fe_system.hh.

Here is the caller graph for this function:

◆ initialize()

template<unsigned int dim>
void FESystem< dim >::initialize
private

Initialization of the internal structures from the vector of base FE.

Definition at line 140 of file fe_system.cc.

◆ update_each()

template<unsigned int dim>
UpdateFlags FESystem< dim >::update_each ( UpdateFlags  flags)
inlineoverridevirtual

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

Parameters
flagsComputed update flags.

Reimplemented from FiniteElement< dim >.

Definition at line 222 of file fe_system.cc.

Member Data Documentation

◆ fe_

template<unsigned int dim>
std::vector<std::shared_ptr<FiniteElement<dim> > > FESystem< dim >::fe_
private

Pointers to base FE objects.

Definition at line 167 of file fe_system.hh.

◆ scalar_components_

template<unsigned int dim>
std::vector<unsigned int> FESystem< dim >::scalar_components_
private

Definition at line 169 of file fe_system.hh.

◆ tensor_components_

template<unsigned int dim>
std::vector<unsigned int> FESystem< dim >::tensor_components_
private

Definition at line 171 of file fe_system.hh.

◆ vector_components_

template<unsigned int dim>
std::vector<unsigned int> FESystem< dim >::vector_components_
private

Definition at line 170 of file fe_system.hh.


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