Flow123d
release_3.0.0-968-gc87a28e79
|
Go to the documentation of this file.
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>
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);
300 const unsigned int point_no,
301 const unsigned int comp)
const;
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];
531 template<
unsigned int dim,
unsigned int spacedim>
576 template<
unsigned int dim,
unsigned int spacedim>
void fill_vec_piola_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
void fill_tensor_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
double shape_value_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
void fill_system_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for mixed system of FE.
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
virtual 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.
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags _flags)
Constructor.
unsigned int n_dofs
Number of dofs (shape functions).
virtual ~FESideValues()
Destructor.
vector< FEValuesViews::Vector< dim, spacedim > > vectors
Mapping data that can be precomputed on the actual cell.
Structure for storing the precomputed finite element data.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
unsigned int n_dofs()
Returns the number of shape functions.
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
FiniteElement< dim > * fe
The used finite element.
void fill_scalar_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
FiniteElement< dim > * get_fe() const
Returns the finite element in use.
Quadrature< dim > * get_quadrature() const
Returns the quadrature in use.
void fill_data(const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
Calculates finite element data on the actual cell.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
Gradients of shape functions evaluated at the quadrature points.
Quadrature< dim > side_quadrature[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Abstract class for the mapping between reference and actual cell.
std::vector< std::shared_ptr< FEValuesBase< dim, spacedim > > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
FEValuesBase()
Default constructor.
Mapping< dim, spacedim > * get_mapping() const
Returns the mapping in use.
virtual ~FEValuesSpaceBase()
FEInternalData * side_fe_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim-1 > &_sub_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Constructor.
void fill_vec_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
const FEValuesViews::Vector< dim, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
unsigned int n_points
Number of quadrature points.
std::vector< arma::mat::fixed< dim, spacedim > > inverse_jacobians
Inverse Jacobians at the quadrature points.
Abstract base class with certain methods independent of the template parameter dim.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
void allocate(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Allocates space for computed data.
Base class for FEValues and FESideValues.
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
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.
std::vector< double > determinants
Determinants of Jacobians at quadrature points.
arma::vec::fixed< spacedim > shape_grad_component(const unsigned int function_no, const unsigned int point_no, const unsigned int comp) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
std::vector< arma::mat::fixed< spacedim, dim > > jacobians
Jacobians of the mapping at the quadrature points.
void initialize(FEValuesBase &fv)
vector< arma::vec::fixed< spacedim > > & point_list()
Return coordinates of all quadrature points in the actual cell system.
virtual arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)=0
Returns the normal vector to a side at given quadrature point.
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
virtual 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.
MappingInternalData * mapping_data
Precomputed mapping data.
Class FEValuesData holds the arrays of data computed by Mapping and FiniteElement.
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
std::vector< double > JxW_values
Transformed quadrature weights.
FEInternalData * init_fe_data(const Quadrature< dim > *q)
Precompute finite element data on reference element.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp)
Resize the data arrays.
FEInternalData * fe_data
Precomputed finite element data.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
std::vector< arma::vec::fixed< spacedim > > normal_vectors
Normal vectors to the element at the quadrature points lying on a side.
void fill_vec_contravariant_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
vector< FEValuesViews::Tensor< dim, spacedim > > tensors
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.
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
Calculates finite element data on a side.
unsigned int n_components_
Number of components of the FE.
FEInternalData(unsigned int np, unsigned int nd)
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
Base class for quadrature rules on simplices in arbitrary dimensions.
vector< FEValuesViews::Scalar< dim, spacedim > > scalars
virtual double JxW(const unsigned int point_no)=0
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
unsigned int n_points()
Returns the number of quadrature points.
ElementAccessor< 3 > * present_cell
Iterator to the last reinit-ed cell.
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
virtual unsigned int n_dofs()=0
Returns the number of shape functions.
const FEValuesViews::Scalar< dim, spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
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.
std::vector< arma::vec::fixed< spacedim > > points
Coordinates of quadrature points in the actual cell coordinate system.