49 template<
unsigned int dim,
unsigned int spacedim>
56 jacobians.resize(size);
59 determinants.resize(size);
63 JxW_values.resize(size);
66 inverse_jacobians.resize(size);
75 shape_gradients.resize(size,
vector<arma::vec::fixed<spacedim> >(n_comp));
82 normal_vectors.resize(size);
86 template<
unsigned int dim,
unsigned int spacedim>
92 switch (fv.
get_fe()->type_) {
106 ASSERT_DBG(fe_sys !=
nullptr).error(
"Mixed system must be represented by FESystem.");
111 unsigned int comp_offset = 0;
112 for (
auto fe : fe_sys->
fe())
132 comp_offset += fe->n_space_components(spacedim);
140 template<
unsigned int dim,
unsigned int spacedim>
142 : mapping(NULL), n_points_(0), fe(NULL), mapping_data(NULL), fe_data(NULL)
148 template<
unsigned int dim,
unsigned int spacedim>
156 template<
unsigned int dim,
unsigned int spacedim>
167 ASSERT_DBG(_fe.
n_components() == spacedim*spacedim).error(
"FETensor must have spacedim*spacedim components.");
183 template<
unsigned int dim,
unsigned int spacedim>
190 for (
unsigned int i=0; i<q->
size(); i++)
192 for (
unsigned int j=0; j<
fe->n_dofs(); j++)
194 for (
unsigned int c=0; c<
fe->n_components(); c++)
195 shape_values(j,c) =
fe->shape_value(j, q->
point<dim>(i).arma(), c);
197 data->ref_shape_values[i][j] = trans(shape_values.row(j));
202 for (
unsigned int i=0; i<q->
size(); i++)
204 for (
unsigned int j=0; j<
fe->n_dofs(); j++)
207 for (
unsigned int c=0; c<
fe->n_components(); c++)
208 grad.col(c) +=
fe->shape_grad(j, q->
point<dim>(i).arma(), c);
210 data->ref_shape_grads[i][j] = grad;
219 template<
unsigned int dim,
unsigned int spacedim>
228 template<
unsigned int dim,
unsigned int spacedim>
233 return data.shape_values[point_no][function_no];
237 template<
unsigned int dim,
unsigned int spacedim>
242 return data.shape_gradients[point_no][function_no];
246 template<
unsigned int dim,
unsigned int spacedim>
248 const unsigned int point_no,
249 const unsigned int comp)
const 258 template<
unsigned int dim,
unsigned int spacedim>
260 const unsigned int point_no,
261 const unsigned int comp)
const 270 template<
unsigned int dim,
unsigned int spacedim>
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 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
284 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
289 template<
unsigned int dim,
unsigned int spacedim>
297 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
298 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
301 for (
unsigned int c=0; c<spacedim; c++)
302 data.shape_values[i][j*spacedim+c] = fv_vec[c];
309 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
310 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
313 for (
unsigned int c=0; c<spacedim; c++)
314 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
320 template<
unsigned int dim,
unsigned int spacedim>
328 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
329 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
332 for (
unsigned int c=0; c<spacedim; c++)
333 data.shape_values[i][j*spacedim+c] = fv_vec[c];
340 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
341 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
344 for (
unsigned int c=0; c<spacedim; c++)
345 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
351 template<
unsigned int dim,
unsigned int spacedim>
359 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
360 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
363 for (
unsigned int c=0; c<spacedim; c++)
364 data.shape_values[i][j*spacedim+c] = fv_vec(c);
371 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
372 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
375 /
data.determinants[i];
376 for (
unsigned int c=0; c<spacedim; c++)
377 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
383 template<
unsigned int dim,
unsigned int spacedim>
391 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
392 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
395 for (
unsigned int c=0; c<spacedim*spacedim; c++)
396 data.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
403 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
404 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
407 for (
unsigned int c=0; c<spacedim*spacedim; c++)
408 data.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
414 template<
unsigned int dim,
unsigned int spacedim>
421 ASSERT_DBG(fe_sys !=
nullptr).error(
"Mixed system must be represented by FESystem.");
422 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
431 unsigned int n_space_components =
fe->n_space_components(spacedim);
437 unsigned int comp_offset = 0;
438 unsigned int shape_offset = 0;
439 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
441 unsigned int n_sub_space_components = fe_sys->
fe()[f]->n_space_components(spacedim);
444 for (
unsigned int i=0; i<fe_data.
n_points; i++)
445 for (
unsigned int n=0; n<fe_sys->
fe()[f]->n_dofs(); n++)
446 for (
unsigned int c=0; c<n_sub_space_components; c++)
447 data.shape_values[i][shape_offset+n_space_components*n+comp_offset+c] =
fe_values_vec[f]->
data.shape_values[i][n*n_sub_space_components+c];
449 comp_offset += n_sub_space_components;
450 shape_offset += fe_sys->
fe()[f]->n_dofs()*n_space_components;
458 unsigned int comp_offset = 0;
459 unsigned int shape_offset = 0;
460 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
462 unsigned int n_sub_space_components = fe_sys->
fe()[f]->n_space_components(spacedim);
465 for (
unsigned int i=0; i<fe_data.
n_points; i++)
466 for (
unsigned int n=0; n<fe_sys->
fe()[f]->n_dofs(); n++)
467 for (
unsigned int c=0; c<n_sub_space_components; c++)
468 data.shape_gradients[i][shape_offset+n_space_components*n+comp_offset+c] =
fe_values_vec[f]->
data.shape_gradients[i][n*n_sub_space_components+c];
470 comp_offset += n_sub_space_components;
471 shape_offset += fe_sys->
fe()[f]->n_dofs()*n_space_components;
478 template<
unsigned int dim,
unsigned int spacedim>
501 ASSERT(
false).error(
"Not implemented.");
512 template<
unsigned int dim,
unsigned int spacedim>
518 quadrature(&_quadrature)
521 this->
allocate(_mapping, _quadrature.
size(), _fe, _flags);
531 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
533 for (
auto fe_sub : fe->
fe())
540 template<
unsigned int dim,
unsigned int spacedim>
544 this->
data.present_cell = &cell;
547 this->
mapping->fill_fe_values(cell,
563 template<
unsigned int dim,
unsigned int spacedim>
574 this->
allocate(_mapping, _sub_quadrature.
size(), _fe, _flags);
576 for (
unsigned int sid = 0; sid < RefElement<dim>::n_sides; sid++)
578 for (
unsigned int pid = 0; pid < RefElement<dim>::n_side_permutations; pid++)
591 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
593 for (
auto fe_sub : fe->
fe())
600 template<
unsigned int dim,
unsigned int spacedim>
603 for (
unsigned int sid=0; sid<RefElement<dim>::n_sides; sid++)
605 for (
unsigned int pid=0; pid<RefElement<dim>::n_side_permutations; pid++)
614 template<
unsigned int dim,
unsigned int spacedim>
620 this->
data.present_cell = &cell;
626 this->
mapping->fill_fe_side_values(cell,
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell...
void fill_vec_contravariant_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for vectorial FE.
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe()
Transformed quadrature weight for cell sides.
#define ASSERT_EQ_DBG(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
void reinit(ElementAccessor< 3 > &cell, unsigned int sid)
Update cell-dependent data (gradients, Jacobians etc.)
void fill_scalar_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for scalar FE.
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...
Determinant of the Jacobian.
void fill_data(const FEInternalData &fe_data)
Computes the shape function values and gradients on the actual cell and fills the FEValues structure...
std::vector< std::vector< Quadrature > > side_quadrature
MappingInternalData * mapping_data
Precomputed mapping data.
Quadrature make_from_side(unsigned int sid, unsigned int pid)
void initialize(FEValuesBase &fv)
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
unsigned int side_perm_
Permutation index of current side.
void allocate(Mapping< dim, spacedim > &_mapping, unsigned int n_points, FiniteElement< dim > &_fe, UpdateFlags flags)
Allocates space for computed data.
FEValuesBase()
Default constructor.
Class FEValues calculates finite element data on the actual cells such as shape function values...
FEValuesData< dim, spacedim > data
Data computed by the mapping and finite element.
Class FESystem for compound finite elements.
#define ASSERT(expr)
Allow use shorter versions of macro names if these names is not used with external library...
FEInternalData * fe_data
Precomputed finite element data.
FiniteElement< dim > * get_fe() const
Returns the finite element in use.
unsigned int n_points()
Returns the number of quadrature points.
void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp)
Resize the data arrays.
unsigned int side_idx_
Current side on which FESideValues was recomputed.
Base class for quadrature rules on simplices in arbitrary dimensions.
unsigned int n_components_
Number of components of the FE.
virtual ~FESideValues()
Destructor.
Quadrature * quadrature
The quadrature rule used to calculate integrals.
Armor::vec< point_dim > point(const unsigned int i) const
Returns the ith quadrature point.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) override
Return the gradient of the function_no-th shape function at the point_no-th quadrature point...
Transformed quadrature points.
Abstract class for the mapping between reference and actual cell.
Compound finite element on dim dimensional simplex.
const unsigned int size() const
Returns number of quadrature points.
void fill_tensor_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for tensorial FE.
unsigned int n_sides() const
void fill_system_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for mixed system of FE.
void fill_vec_piola_data(const FEInternalData &fe_data)
Compute shape functions and gradients on the actual cell for Raviart-Thomas FE.
unsigned int n_components() const
Returns numer of components of the basis function.
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature &_quadrature, FiniteElement< dim > &_fe, UpdateFlags _flags)
Constructor.
Basic definitions of numerical quadrature rules.
Shape function gradients.
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
FEInternalData(unsigned int np, unsigned int nd)
Class Mapping calculates data related to the mapping of the reference cell to the actual cell...
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature &_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.
std::vector< std::vector< arma::mat > > ref_shape_grads
Precomputed gradients of basis functions at the quadrature points.
unsigned int n_points_
Number of integration points.
FEInternalData * side_fe_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
const Quadrature * sub_quadrature
Quadrature for the integration on the element sides.
double shape_value(const unsigned int function_no, const unsigned int point_no) override
Return the value of the function_no-th shape function at the point_no-th quadrature point...
Mapping< dim, spacedim > * mapping
The mapping from the reference cell to the actual cell.
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...
FiniteElement< dim > * fe
The used finite element.
Calculates finite element data on the actual cell.
unsigned int n_points
Number of quadrature points.
FEInternalData * init_fe_data(const Quadrature *q)
Precompute finite element data on reference element.
#define OLD_ASSERT_EQUAL(a, b)
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 description of finite elements.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Base class for FEValues and FESideValues.
unsigned int n_dofs
Number of dofs (shape functions).
unsigned int permutation_idx(unsigned int prm_idx) const
Return permutation_idx of given index.
void reinit(ElementAccessor< 3 > &cell)
Update cell-dependent data (gradients, Jacobians etc.)
Transformed quadrature weights.
FEType type_
Type of FiniteElement.
Calculates finite element data on a side.
#define ASSERT_LT_DBG(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.