43 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
44 template<
unsigned int dim,
unsigned int spacedim>
class Mapping;
56 template<
unsigned int dim,
unsigned int spacedim>
145 template<
unsigned int spacedim>
156 virtual const double shape_value(
const unsigned int function_no,
const unsigned int point_no) = 0;
165 virtual const arma::vec::fixed<spacedim>
shape_grad(
const unsigned int function_no,
const unsigned int point_no) = 0;
173 virtual const double JxW(
const unsigned int point_no) = 0;
180 virtual const arma::vec::fixed<spacedim>
normal_vector(
unsigned int point_no) = 0;
185 virtual const unsigned int n_dofs() = 0;
192 template<
unsigned int dim,
unsigned int spacedim>
236 const double shape_value(
const unsigned int function_no,
const unsigned int point_no);
245 const arma::vec::fixed<spacedim>
shape_grad(
const unsigned int function_no,
const unsigned int point_no);
256 const arma::vec::fixed<spacedim>
shape_vector(
const unsigned int function_no,
const unsigned int point_no);
267 const arma::mat::fixed<spacedim,spacedim>
shape_grad_vector(
const unsigned int function_no,
const unsigned int point_no);
277 const double determinant(
const unsigned int point_no);
285 const double JxW(
const unsigned int point_no);
292 const arma::vec::fixed<spacedim>
point(
const unsigned int point_no);
305 const arma::vec::fixed<spacedim>
normal_vector(
unsigned int point_no);
315 const unsigned int n_dofs();
369 template<
unsigned int dim,
unsigned int spacedim>
414 template<
unsigned int dim,
unsigned int spacedim>
virtual const double JxW(const unsigned int point_no)=0
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
std::vector< arma::vec > shape_values
Shape functions evaluated at the quadrature points.
void allocate(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim, spacedim > &_fe, UpdateFlags flags)
Allocates space for computed data.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
Abstract class for the mapping between reference and actual cell.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Base class for FEValues and FESideValues.
const arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
Base class for quadrature rules on simplices in arbitrary dimensions.
const double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
std::vector< std::vector< arma::mat::fixed< spacedim, spacedim > > > shape_grad_vectors
Gradients of shape functions (for vectorial finite elements).
std::vector< double > JxW_values
Transformed quadrature weights.
const double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
FEValuesBase()
Default constructor.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
virtual ~FESideValues()
Destructor.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Normal vectors to the element at the quadrature points lying on a side.
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
const unsigned int n_dofs()
Returns the number of shape functions.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.
const Quadrature< dim > * get_quadrature() const
Returns the quadrature in use.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
virtual const arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)=0
Returns the normal vector to a side at given quadrature point.
MappingInternalData * mapping_data
Precomputed mapping data.
FEInternalData * fe_data
Precomputed finite element data.
std::vector< arma::mat > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
ElementFullIter * present_cell
Iterator to the last reinit-ed cell.
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
const vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
const arma::mat::fixed< spacedim, spacedim > shape_grad_vector(const unsigned int function_no, const unsigned int point_no)
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
virtual const arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no)=0
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Quadrature< dim > side_quadrature[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
const arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no)
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
virtual const unsigned int n_dofs()=0
Returns the number of shape functions.
void allocate(unsigned int size, UpdateFlags flags, bool is_scalar=true)
Resize the data arrays.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_vectors
Shape functions (for vectorial finite elements) evaluated at quadrature points.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
FEInternalData * side_fe_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
virtual const double shape_value(const unsigned int function_no, const unsigned int point_no)=0
Return the value of the function_no-th shape function at the point_no-th quadrature point...
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim, spacedim > &_fe, UpdateFlags _flags)
Constructor.
void reinit(ElementFullIter &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
void reinit(ElementFullIter &cell)
Update cell-dependent data (gradients, Jacobians etc.)
const arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
FiniteElement< dim, spacedim > * fe
The used finite element.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Calculates finite element data on the actual cell.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
const double shape_value(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Structure for storing the precomputed finite element data.
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim-1 > &_sub_quadrature, FiniteElement< dim, spacedim > &_fe, UpdateFlags flags)
Constructor.
const unsigned int n_points()
Returns the number of quadrature points.
const arma::vec::fixed< spacedim > shape_vector(const unsigned int function_no, const unsigned int point_no)
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Abstract base class with certain methods independent of the template parameter dim.
Mapping data that can be precomputed on the actual cell.
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Calculates finite element data on a side.