19 #ifndef FINITE_ELEMENT_HH_ 20 #define FINITE_ELEMENT_HH_ 26 #include <boost/assign/list_of.hpp> 27 #include <boost/exception/info.hpp> 35 template<
unsigned int dim>
class FESystem;
36 template<
unsigned int dim,
unsigned int spacedim>
class FESideValues;
37 template<
unsigned int dim,
unsigned int spacedim>
class FEValues;
38 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
39 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesData;
40 template<
unsigned int dim>
class FE_P_disc;
77 Dof(
unsigned int dim_,
78 unsigned int n_face_idx_,
92 double evaluate(
const FS &function_space,
93 unsigned int basis_idx)
const;
129 unsigned int n_components)
131 : space_dim_(space_dim),
132 n_components_(n_components)
142 virtual double basis_value(
unsigned int basis_index,
144 unsigned int comp_index = 0
153 virtual const arma::vec basis_grad(
unsigned int basis_index,
155 unsigned int comp_index = 0
159 virtual unsigned int dim()
const = 0;
253 template<
unsigned int dim>
267 {
return dofs_.size(); }
278 double shape_value(
const unsigned int i,
279 const arma::vec::fixed<dim> &p,
280 const unsigned int comp = 0)
const;
291 arma::vec::fixed<dim> shape_grad(
const unsigned int i,
292 const arma::vec::fixed<dim> &p,
293 const unsigned int comp = 0)
const;
297 {
return function_space_->n_components(); }
300 inline const Dof &
dof(
unsigned int i)
const 304 unsigned int n_space_components(
unsigned int spacedim);
316 void init(
bool primitive =
true,
322 void setup_components();
338 virtual void compute_node_matrix();
345 {
return is_primitive_; }
352 {
return component_indices_[sys_idx]; }
359 {
return nonzero_components_[sys_idx]; }
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
arma::vec coords
Barycentric coordinates.
FunctionSpace(unsigned int space_dim, unsigned int n_components)
arma::mat node_matrix
Matrix that determines the coefficients of the raw basis functions from the values at the support poi...
unsigned int n_components_
Number of components of function values.
bool is_primitive_
Primitive FE is using componentwise shape functions, i.e. only one component is nonzero for each shap...
unsigned int dim
Association to n-face of given dimension (point, line, triangle, tetrahedron.
std::vector< std::vector< bool > > nonzero_components_
Footprints of nonzero components of shape functions.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
virtual ~FiniteElement()
Destructor.
unsigned int n_face_idx
Index of n-face to which the dof is associated.
unsigned int space_dim() const
Getter for space dimension.
Discontinuous Lagrangean finite element on dim dimensional simplex.
Compound finite element on dim dimensional simplex.
std::vector< Dof > dofs_
Set of degrees of freedom (functionals) defining the FE.
unsigned int n_components() const
Getter for number of components.
unsigned int n_components() const
Returns numer of components of the basis function.
bool is_primitive() const
Indicates whether the basis functions have one or more nonzero components (scalar FE spaces are alway...
unsigned int system_to_component_index(unsigned sys_idx) const
Returns the component index for vector valued finite elements.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
const Dof & dof(unsigned int i) const
Returns i -th degree of freedom.
double evaluate(const FS &function_space, unsigned int basis_idx) const
Evaulate dof for basis function of given function space.
Calculates finite element data on the actual cell.
Dof(unsigned int dim_, unsigned int n_face_idx_, arma::vec coords_, arma::vec coefs_, DofType type_)
unsigned int space_dim_
Space dimension of function arguments (i.e. 1, 2 or 3).
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
arma::vec coefs
Coefficients of linear combination of function value components.
std::vector< unsigned int > component_indices_
Indices of nonzero components of shape functions (for primitive FE).
const std::vector< bool > & get_nonzero_components(unsigned int sys_idx) const
Returns the mask of nonzero components for given basis function.
FEType type_
Type of FiniteElement.
Calculates finite element data on a side.