39 template<
unsigned int spacedim>
49 template<
unsigned int spacedim>
52 unsigned int first_component_idx,
56 for (
unsigned int ip=0; ip<
n_points; ip++)
57 for (
unsigned int id=0;
id<dof_indices.size();
id++)
66 template<
unsigned int spacedim>
67 template<
unsigned int DIM>
87 ASSERT_DBG(fe_sys !=
nullptr).error(
"Mixed system must be represented by FESystem.");
92 unsigned int comp_offset = 0;
93 for (
auto fe : fe_sys->
fe())
121 template<
unsigned int spacedim>
129 template<
unsigned int spacedim>
134 template<
unsigned int spacedim>
135 template<
unsigned int DIM>
141 if (DIM == 0)
return;
150 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
153 for (
unsigned int f=0; f<fe->
fe().size(); f++)
158 if ( q.dim() == DIM )
162 else if ( q.dim() + 1 == DIM )
165 for (
unsigned int sid = 0; sid < RefElement<DIM>::n_sides; sid++)
171 for (
unsigned int pid = 0; pid < RefElement<DIM>::n_side_permutations; pid++)
176 ASSERT_DBG(
false)(q.dim())(DIM).error(
"Dimension mismatch in FEValues::initialize().");
181 template<
unsigned int spacedim>
182 template<
unsigned int DIM>
193 ASSERT_DBG(_fe.
n_components() == spacedim*spacedim).error(
"FETensor must have spacedim*spacedim components.");
206 if (fe_sys !=
nullptr)
208 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
230 template<
unsigned int spacedim>
231 template<
unsigned int DIM>
236 std::shared_ptr<FEInternalData> data = std::make_shared<FEInternalData>(q.
size(),
n_dofs_);
239 for (
unsigned int i=0; i<q.
size(); i++)
241 for (
unsigned int j=0; j<
n_dofs_; j++)
246 data->ref_shape_values[i][j] = trans(
shape_values.row(j));
251 for (
unsigned int i=0; i<q.
size(); i++)
253 for (
unsigned int j=0; j<
n_dofs_; j++)
259 data->ref_shape_grads[i][j] = grad;
267 template<
unsigned int spacedim>
276 template<
unsigned int spacedim>
285 template<
unsigned int spacedim>
287 const unsigned int point_no,
288 const unsigned int comp)
const 297 template<
unsigned int spacedim>
299 const unsigned int point_no,
300 const unsigned int comp)
const 309 template<
unsigned int spacedim>
316 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
317 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
322 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
323 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
328 template<
unsigned int spacedim>
337 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
338 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
341 for (
unsigned int c=0; c<spacedim; c++)
349 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
350 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
353 for (
unsigned int c=0; c<spacedim; c++)
360 template<
unsigned int spacedim>
369 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
370 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
373 for (
unsigned int c=0; c<spacedim; c++)
381 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
382 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
385 for (
unsigned int c=0; c<spacedim; c++)
392 template<
unsigned int spacedim>
401 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
402 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
405 for (
unsigned int c=0; c<spacedim; c++)
413 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
414 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
418 for (
unsigned int c=0; c<spacedim; c++)
425 template<
unsigned int spacedim>
434 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
435 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
438 for (
unsigned int c=0; c<spacedim*spacedim; c++)
446 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
447 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
450 for (
unsigned int c=0; c<spacedim*spacedim; c++)
457 template<
unsigned int spacedim>
463 unsigned int comp_offset = 0;
477 unsigned int comp_offset = 0;
478 unsigned int shape_offset = 0;
482 for (
unsigned int i=0; i<fe_data.
n_points; i++)
487 comp_offset += fe_sys_n_space_components_[f];
496 unsigned int comp_offset = 0;
497 unsigned int shape_offset = 0;
501 for (
unsigned int i=0; i<fe_data.
n_points; i++)
506 comp_offset += fe_sys_n_space_components_[f];
514 template<
unsigned int spacedim>
537 ASSERT(
false).error(
"Not implemented.");
547 template<
unsigned int spacedim>
562 template<
unsigned int spacedim>
574 const unsigned int pid =
elm_values->side().element()->permutation_idx(sid);
588 fv[0].initialize(quadrature[0], *fe[0_d], flags);
589 fv[1].initialize(quadrature[1], *fe[1_d], flags);
590 fv[2].initialize(quadrature[2], *fe[2_d], flags);
591 fv[3].initialize(quadrature[3], *fe[3_d], flags);
std::shared_ptr< ElementValues< spacedim > > elm_values
Auxiliary object for calculation of element-dependent data.
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell...
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
unsigned int n_dofs() const
Returns the number of degrees of freedom needed by the finite element.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
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...
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
unsigned int side_idx() const
Returns local index of the side on the element.
~FEValues()
Correct deallocation of objects created by 'initialize' methods.
FEType fe_type_
Type of finite element (scalar, vector, tensor).
std::shared_ptr< FEInternalData > fe_data
Precomputed finite element data.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
arma::vec::fixed< dim > shape_grad(const unsigned int i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the comp-th component of the gradient of the i-th shape function at the point p on the ref...
static UpdateFlags update_each(UpdateFlags flags)
Determines which additional quantities have to be computed.
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...
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...
void fill_system_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for mixed system of FE.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
double determinant(const unsigned int point_no) const
Return the relative volume change of the cell (Jacobian determinant).
Class FEValues calculates finite element data on the actual cells such as shape function values...
Class for computation of data on cell and side.
arma::mat inverse_jacobian(const unsigned int point_no) const
Return inverse Jacobian matrix at point point_no.
Class FESystem for compound finite elements.
void fill_vec_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
void initialize(const FEValues &fv, const FiniteElement< DIM > &fe)
void reinit(const ElementAccessor< spacedim > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
void allocate(unsigned int n_points, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
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...
Base class for quadrature rules on simplices in arbitrary dimensions.
arma::mat jacobian(const unsigned int point_no) const
Return Jacobian matrix at point point_no.
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
ArmaMat< double, N, M > mat
unsigned int n_points
Number of quadrature points.
unsigned int dim() const
Returns dimension of the side, that is dimension of the element minus one.
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 i, const arma::vec::fixed< dim > &p, const unsigned int comp=0) const
Calculates the value of the comp-th component of the i-th shape function at the point p on the refere...
Compound finite element on dim dimensional simplex.
virtual UpdateFlags update_each(UpdateFlags flags)
Decides which additional quantities have to be computed for each cell.
std::vector< unsigned int > fe_dofs(unsigned int fe_index)
Return dof indices belonging to given sub-FE.
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
unsigned int n_points_
Number of integration points.
unsigned int n_components() const
Returns numer of components of the basis function.
Basic definitions of numerical quadrature rules.
Shape function gradients.
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
unsigned int n_space_components(unsigned int spacedim)
Number of components of FE in a mapped space with dimension spacedim.
unsigned int dim_
Dimension of reference space.
unsigned int n_dofs
Number of dofs (shape functions).
Structure for storing the precomputed finite element data.
void fill_vec_piola_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
std::array< QGauss, 4 > array
int LongIdx
Define type that represents indices of large arrays (elements, nodes, dofs etc.)
FEInternalData(unsigned int np, unsigned int nd)
void fill_vec_contravariant_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
std::vector< std::vector< shared_ptr< FEInternalData > > > side_fe_data
Precomputed FE data (shape functions on reference element) for all sides and permuted quadrature poin...
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
Calculates finite element data on the actual cell.
unsigned int n_dofs_
Number of finite element dofs.
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
Abstract class for description of finite elements.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int n_components_
Number of components of the FE.
std::vector< FEValues< spacedim > > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
unsigned int size() const
Returns number of quadrature points.
void fill_tensor_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
FEValues()
Default constructor with postponed initialization.
std::vector< FEValues< 3 > > mixed_fe_values(QGauss::array &quadrature, MixedPtr< FiniteElement > fe, UpdateFlags flags)
void fill_scalar_data(const ElementValues< spacedim > &elm_values, const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
unsigned int n_points() const
Returns the number of quadrature points.
std::vector< std::vector< arma::vec > > ref_shape_values
Precomputed values of basis functions at the quadrature points.
std::vector< std::vector< double > > shape_values
Shape functions evaluated at the quadrature points.
FEType type_
Type of FiniteElement.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.