32 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
33 template<
unsigned int dim,
unsigned int spacedim>
class Mapping;
45 template<
unsigned int dim,
unsigned int spacedim>
134 template<
unsigned int spacedim>
146 virtual double shape_value(
const unsigned int function_no,
const unsigned int point_no) = 0;
155 virtual arma::vec::fixed<spacedim> shape_grad(
const unsigned int function_no,
const unsigned int point_no) = 0;
163 virtual double JxW(
const unsigned int point_no) = 0;
170 virtual arma::vec::fixed<spacedim> normal_vector(
unsigned int point_no) = 0;
175 virtual unsigned int n_dofs() = 0;
182 template<
unsigned int dim,
unsigned int spacedim>
226 inline double shape_value(
const unsigned int function_no,
const unsigned int point_no)
228 return data.shape_values[point_no][function_no];
239 inline arma::vec::fixed<spacedim>
shape_grad(
const unsigned int function_no,
const unsigned int point_no)
241 return trans(data.shape_gradients[point_no].row(function_no));
254 inline arma::vec::fixed<spacedim>
shape_vector(
const unsigned int function_no,
const unsigned int point_no)
256 return data.shape_vectors[point_no][function_no];
268 inline arma::mat::fixed<spacedim,spacedim>
shape_grad_vector(
const unsigned int function_no,
const unsigned int point_no)
270 return data.shape_grad_vectors[point_no][function_no];
283 return data.determinants[point_no];
292 inline double JxW(
const unsigned int point_no)
294 return data.JxW_values[point_no];
302 inline arma::vec::fixed<spacedim>
point(
const unsigned int point_no)
304 return data.points[point_no];
324 return data.normal_vectors[point_no];
332 return quadrature->size();
415 template<
unsigned int dim,
unsigned int spacedim>
460 template<
unsigned int dim,
unsigned int spacedim>
std::vector< arma::vec > shape_values
Shape functions evaluated at the quadrature points.
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point...
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
unsigned int n_dofs()
Returns the number of shape functions.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
MappingInternalData * mapping_data
Precomputed mapping data.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
FiniteElement< dim, spacedim > * get_fe() const
Returns the finite element in use.
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
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.
FEInternalData * fe_data
Precomputed finite element data.
unsigned int n_points()
Returns the number of quadrature points.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
Base class for quadrature rules on simplices in arbitrary dimensions.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Normal vectors to the element at the quadrature points lying on a side.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
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.
Abstract class for the mapping between reference and actual cell.
FiniteElement< dim, spacedim > * fe
The used finite element.
std::vector< arma::mat > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
virtual ~FEValuesSpaceBase()
ElementFullIter * present_cell
Iterator to the last reinit-ed cell.
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
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...
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...
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
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...
Mapping< dim, spacedim > * get_mapping() const
Returns the mapping in use.
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.
Quadrature< dim > * get_quadrature() const
Returns the quadrature in use.
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Calculates finite element data on the actual cell.
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc...
Structure for storing the precomputed finite element data.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
Abstract base class with certain methods independent of the template parameter dim.
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...
Mapping data that can be precomputed on the actual cell.
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Calculates finite element data on a side.
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.