37 template<
unsigned int dim,
unsigned int spacedim>
class FEValuesBase;
38 template<
unsigned int dim,
unsigned int spacedim>
class Mapping;
85 template<
unsigned int dim,
unsigned int spacedim>
96 void allocate(
unsigned int size,
UpdateFlags flags,
unsigned n_comp);
174 template<
unsigned int spacedim>
186 virtual double shape_value(
const unsigned int function_no,
const unsigned int point_no) = 0;
195 virtual arma::vec::fixed<spacedim> shape_grad(
const unsigned int function_no,
const unsigned int point_no) = 0;
203 virtual double JxW(
const unsigned int point_no) = 0;
210 virtual arma::vec::fixed<spacedim> normal_vector(
unsigned int point_no) = 0;
215 virtual unsigned int n_dofs() = 0;
222 template<
unsigned int dim,
unsigned int spacedim>
278 double shape_value(
const unsigned int function_no,
const unsigned int point_no);
288 arma::vec::fixed<spacedim> shape_grad(
const unsigned int function_no,
const unsigned int point_no);
299 double shape_value_component(
const unsigned int function_no,
300 const unsigned int point_no,
301 const unsigned int comp)
const;
312 arma::vec::fixed<spacedim> shape_grad_component(
const unsigned int function_no,
313 const unsigned int point_no,
314 const unsigned int comp)
const;
327 return data.determinants[point_no];
336 inline double JxW(
const unsigned int point_no)
339 return data.JxW_values[point_no];
347 inline arma::vec::fixed<spacedim>
point(
const unsigned int point_no)
350 return data.points[point_no];
371 return data.normal_vectors[point_no];
381 return views_cache_.scalars[i];
391 return views_cache_.vectors[i];
401 return views_cache_.tensors[i];
409 return quadrature->size();
531 template<
unsigned int dim,
unsigned int spacedim>
576 template<
unsigned int dim,
unsigned int spacedim>
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.
const FEValuesViews::Tensor< dim, spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
std::vector< double > JxW_values
Transformed quadrature weights.
FEInternalData * fe_data
Precomputed finite element data.
FiniteElement< dim > * get_fe() const
Returns the finite element in use.
vector< FEValuesViews::Tensor< dim, spacedim > > tensors
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.
vector< FEValuesViews::Vector< dim, spacedim > > vectors
unsigned int n_components_
Number of components of the FE.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Shape functions (for vectorial finite elements) evaluated at quadrature points.
vector< FEValuesViews::Scalar< dim, spacedim > > scalars
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.
virtual ~FEValuesSpaceBase()
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
FEInternalData(unsigned int np, unsigned int nd)
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
Mapping< dim, spacedim > * get_mapping() const
Returns the mapping in use.
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
Quadrature< dim > * get_quadrature() const
Returns the quadrature in use.
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
FiniteElement< dim > * fe
The used finite element.
ElementAccessor< 3 > * present_cell
Iterator to the last reinit-ed cell.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
Calculates finite element data on the actual cell.
unsigned int n_points
Number of quadrature points.
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...
std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
Structure for storing the precomputed finite element data.
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
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.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
unsigned int n_dofs
Number of dofs (shape functions).
Mapping data that can be precomputed on the actual cell.
const FEValuesViews::Scalar< dim, spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
Calculates finite element data on a side.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.