41 template<
class FV,
unsigned int dim>
class MapScalar;
42 template<
class FV,
unsigned int dim>
class MapPiola;
44 template<
class FV,
unsigned int dim>
class MapVector;
45 template<
class FV,
unsigned int dim>
class MapTensor;
46 template<
class FV,
unsigned int dim>
class MapSystem;
48 template<
unsigned int spcedim>
class FEValues;
70 unsigned int first_component_idx,
71 unsigned int ncomps = 1);
99 template<
class FV,
unsigned int spacedim = 3>
110 template<
unsigned int DIM>
127 template<
unsigned int DIM>
139 template<
unsigned int DIM>
159 inline unsigned int dim()
const
206 template<
class MapType>
220 template<
unsigned int DIM>
266 friend class MapPiola<FV, spacedim>;
287 template<
unsigned int spacedim = 3>
297 template<
unsigned int DIM>
332 inline double shape_value(
const unsigned int function_no,
const unsigned int point_no)
const
347 inline arma::vec::fixed<spacedim>
shape_grad(
const unsigned int function_no,
const unsigned int point_no)
const
364 const unsigned int point_no,
365 const unsigned int comp)
const
383 const unsigned int point_no,
384 const unsigned int comp)
const;
389 inline void set_shape_value(
unsigned int i_point,
unsigned int i_func_comp,
double val)
397 inline void set_shape_gradient(
unsigned int i_point,
unsigned int i_func_comp, arma::vec::fixed<spacedim> val)
422 inline double JxW(
const unsigned int point_no)
435 inline arma::vec::fixed<spacedim>
point(
const unsigned int point_no)
498 template<
unsigned int spacedim = 3>
516 inline void get_cell(
const unsigned int patch_cell_idx) {
521 inline void get_side(
unsigned int patch_cell_idx,
unsigned int side_idx) {
532 inline double shape_value(
const unsigned int function_no,
const unsigned int point_no)
const
578 inline arma::vec::fixed<spacedim>
shape_grad(
const unsigned int function_no,
const unsigned int point_no)
const
624 const unsigned int point_no,
625 const unsigned int comp)
const
643 const unsigned int point_no,
644 const unsigned int comp)
const;
650 inline void set_shape_value(
unsigned int i_point,
unsigned int i_func_comp,
double val)
658 inline void set_shape_gradient(
unsigned int i_point,
unsigned int i_func_comp, arma::vec::fixed<spacedim> val)
683 inline double JxW(
const unsigned int point_no)
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
Base point accessor class.
unsigned int elem_patch_idx() const
unsigned int eval_point_idx() const
Return index in EvalPoints object.
Class for computation of data on cell and side.
Structure for storing the precomputed finite element data.
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
unsigned int n_points
Number of quadrature points.
FEInternalData(unsigned int np, unsigned int nd)
unsigned int n_dofs
Number of dofs (shape functions).
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
void allocate(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
unsigned int n_points_
Number of integration points.
FV * fv_
Helper object, we need its for ViewsCache initialization.
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
std::vector< shared_ptr< FEInternalData > > side_fe_data_
Precomputed FE data (shape functions on reference element) for all side quadrature points.
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
FEType fe_type_
Type of finite element (scalar, vector, tensor).
unsigned int n_dofs_
Number of finite element dofs.
const FEValuesViews::Scalar< FV, spacedim > & scalar_view(unsigned int i) const
Accessor to scalar values of multicomponent FE.
unsigned int n_points() const
Returns the number of quadrature points.
virtual void init_fe_val_vec()=0
Initialize fe_values_vec only in PatchFEValues_TEMP.
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
virtual void allocate_in(unsigned int)=0
Initialize vectors declared separately in descendants.
void fill_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure.
std::vector< FV > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
const FEValuesViews::Tensor< FV, spacedim > & tensor_view(unsigned int i) const
Accessor to tensor values of multicomponent FE.
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
unsigned int n_components_
Number of components of the FE.
unsigned int dim() const
Return dimension of reference space.
FEValuesBase()
Default constructor with postponed initialization.
void fill_data_specialized(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure....
unsigned int n_dofs() const
Returns the number of shape functions.
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
virtual void initialize_in(Quadrature &, unsigned int)=0
Initialize ElementValues separately in descendants.
std::shared_ptr< FEInternalData > fe_data_
Precomputed finite element data.
const FEValuesViews::Vector< FV, spacedim > & vector_view(unsigned int i) const
Accessor to vector values of multicomponent FE.
unsigned int dim_
Dimension of reference space.
Calculates finite element data on the actual cell.
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
~FEValues()
Correct deallocation of objects created by 'initialize' methods.
void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
FEValues()
Default constructor with postponed initialization.
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.
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
const Armor::array & point_list() const
Return coordinates of all quadrature points in the actual cell system.
FEValues(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
std::vector< std::vector< double > > shape_values_
Shape functions evaluated at the quadrature points.
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
void initialize_in(Quadrature &q, unsigned int dim) override
Implement FEValuesBase::initialize_in.
void init_fe_val_vec() override
Implement FEValuesBase::initialize_in.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed< spacedim > val)
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.
std::shared_ptr< ElementValues< spacedim > > elm_values_
Auxiliary object for calculation of element-dependent data.
arma::vec::fixed< spacedim > point(const unsigned int point_no)
Return coordinates of the quadrature point in the actual cell system.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients_
void allocate_in(unsigned int q_dim) override
Implement FEValuesBase::allocate_in.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Helper class allows update values and gradients of FEValues of FEVectorContravariant type.
Helper class allows update values and gradients of FEValues of FEVectorPiola type.
Helper class allows update values and gradients of FEValues of FEScalar type.
Helper class allows update values and gradients of FEValues of FEMixedSystem type.
Helper class allows update values and gradients of FEValues of FETensor type.
Helper class allows update values and gradients of FEValues of FEVector type.
std::shared_ptr< ElementValues< spacedim > > elm_values_
Auxiliary object for calculation of element-dependent data.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients_
std::vector< std::vector< double > > shape_values_
Shape functions evaluated at the quadrature points.
void get_cell(const unsigned int patch_cell_idx)
Set element that is selected for processing. Element is given by index on patch.
PatchFEValues_TEMP(unsigned int max_size=0)
Constructor, set maximal number of elements on patch.
double shape_value(const unsigned int function_no, const BulkPoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
arma::vec::fixed< spacedim > normal_vector(unsigned int point_no)
Returns the normal vector to a side at given quadrature point.
unsigned int used_size() const
arma::vec::fixed< spacedim > normal_vector(const SidePoint &p)
Returns the normal vector to a side at given quadrature point.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const SidePoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
void allocate_in(unsigned int q_dim) override
Implement FEValuesBase::allocate_in.
double JxW(const unsigned int point_no)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
unsigned int max_n_elem_
Maximal number of elements on patch.
void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
MeshObjectType object_type_
Distinguishes using of PatchFEValues_TEMP for storing data of elements or sides.
std::map< unsigned int, unsigned int > element_patch_map_
Map of element patch indexes to element_data_.
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.
void get_side(unsigned int patch_cell_idx, unsigned int side_idx)
Set element and side that are selected for processing. Element is given by index on patch.
double shape_value(const unsigned int function_no, const SidePoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
double determinant(const unsigned int point_no)
Return the relative volume change of the cell (Jacobian determinant).
unsigned int max_size() const
std::vector< ElementFEData > element_data_
Data of elements / sides on patch.
unsigned int patch_data_idx_
Patch index of processed element / side.
void initialize_in(Quadrature &q, unsigned int dim) override
Implement FEValuesBase::initialize_in.
double JxW(const SidePoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
void init_fe_val_vec() override
Implement FEValuesBase::initialize_in.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const BulkPoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
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.
double JxW(const BulkPoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
unsigned int used_size_
Number of elements / sides on patch. Must be less or equal to size of element_data vector.
void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed< spacedim > val)
void resize(unsigned int max_size)
Set size of ElementFEData. Important: Use only during the initialization of FESystem !
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
void reinit(PatchElementsList patch_elements)
Reinit data.
std::array< QGauss, 4 > array
Base class for quadrature rules on simplices in arbitrary dimensions.
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
virtual unsigned int side_idx() const =0
Intermediate step in implementation of PatcFEValues.
unsigned int local_point_idx() const
Return local index in quadrature. Temporary method - intermediate step in implementation of PatcFEVal...
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
std::list< std::pair< ElementAccessor< 3 >, unsigned int > > PatchElementsList
Definitions of particular quadrature rules on simplices.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
vector< FEValuesViews::Scalar< FV, spacedim > > scalars
vector< FEValuesViews::Vector< FV, spacedim > > vectors
void initialize(const FV &fv, const FiniteElement< DIM > &fe)
vector< FEValuesViews::Tensor< FV, spacedim > > tensors
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.