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>
133 for (
unsigned int sid=0; sid<
side_fe_data.size(); sid++)
134 for (
unsigned int pid=0; pid<
side_fe_data[sid].size(); pid++)
140 template<
unsigned int spacedim>
141 template<
unsigned int DIM>
152 ASSERT_DBG(_fe.
n_components() == spacedim*spacedim).error(
"FETensor must have spacedim*spacedim components.");
161 if (fe_sys !=
nullptr)
163 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
185 template<
unsigned int spacedim>
186 template<
unsigned int DIM>
194 for (
unsigned int i=0; i<q.
size(); i++)
196 for (
unsigned int j=0; j<
n_dofs_; j++)
201 data->ref_shape_values[i][j] = trans(
shape_values.row(j));
206 for (
unsigned int i=0; i<q.
size(); i++)
208 for (
unsigned int j=0; j<
n_dofs_; j++)
214 data->ref_shape_grads[i][j] = grad;
222 template<
unsigned int spacedim>
231 template<
unsigned int spacedim>
240 template<
unsigned int spacedim>
242 const unsigned int point_no,
243 const unsigned int comp)
const 252 template<
unsigned int spacedim>
254 const unsigned int point_no,
255 const unsigned int comp)
const 264 template<
unsigned int spacedim>
271 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
272 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
277 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
278 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
283 template<
unsigned int spacedim>
292 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
293 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
296 for (
unsigned int c=0; c<spacedim; c++)
304 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
305 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
308 for (
unsigned int c=0; c<spacedim; c++)
315 template<
unsigned int spacedim>
324 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
325 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
328 for (
unsigned int c=0; c<spacedim; c++)
336 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
337 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
340 for (
unsigned int c=0; c<spacedim; c++)
347 template<
unsigned int spacedim>
356 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
357 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
360 for (
unsigned int c=0; c<spacedim; c++)
368 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
369 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
373 for (
unsigned int c=0; c<spacedim; c++)
380 template<
unsigned int spacedim>
389 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
390 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
393 for (
unsigned int c=0; c<spacedim*spacedim; c++)
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*spacedim; c++)
412 template<
unsigned int spacedim>
418 unsigned int comp_offset = 0;
432 unsigned int comp_offset = 0;
433 unsigned int shape_offset = 0;
437 for (
unsigned int i=0; i<fe_data.
n_points; i++)
442 comp_offset += fe_sys_n_space_components_[f];
451 unsigned int comp_offset = 0;
452 unsigned int shape_offset = 0;
456 for (
unsigned int i=0; i<fe_data.
n_points; i++)
461 comp_offset += fe_sys_n_space_components_[f];
469 template<
unsigned int spacedim>
492 ASSERT(
false).error(
"Not implemented.");
503 template<
unsigned int spacedim>
504 template<
unsigned int DIM>
510 if (DIM == 0)
return;
519 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
522 for (
unsigned int f=0; f<fe->
fe().size(); f++)
527 if ( q.dim() == DIM )
531 else if ( q.dim() + 1 == DIM )
534 for (
unsigned int sid = 0; sid < RefElement<DIM>::n_sides; sid++)
540 for (
unsigned int pid = 0; pid < RefElement<DIM>::n_side_permutations; pid++)
545 ASSERT_DBG(
false)(q.dim())(DIM).error(
"Dimension mismatch in FEValues::initialize().");
552 template<
unsigned int spacedim>
567 template<
unsigned int spacedim>
579 const unsigned int pid =
elm_values->side().element()->permutation_idx(sid);
593 fv[0].initialize(quadrature[0], *fe[0_d], flags);
594 fv[1].initialize(quadrature[1], *fe[1_d], flags);
595 fv[2].initialize(quadrature[2], *fe[2_d], flags);
596 fv[3].initialize(quadrature[3], *fe[3_d], flags);
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.
std::vector< std::vector< typename FEValues< spacedim >::FEInternalData * > > side_fe_data
Precomputed FE data (shape functions on reference element) for all sides and permuted quadrature poin...
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::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
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.
FEValues< spacedim >::FEInternalData * fe_data
Precomputed finite element data.
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.
FEInternalData * init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
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.
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.
ElementValues< spacedim > * elm_values
Auxiliary object for calculation of element-dependent data.
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.