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>
189 for (
unsigned int i=0; i<q->
size(); i++)
191 for (
unsigned int j=0; j<
fe->n_dofs(); j++)
193 for (
unsigned int c=0; c<
fe->n_components(); c++)
194 shape_values(j,c) =
fe->shape_value(j, q->
point(i), c);
201 for (
unsigned int i=0; i<q->
size(); i++)
203 for (
unsigned int j=0; j<
fe->n_dofs(); j++)
206 for (
unsigned int c=0; c<
fe->n_components(); c++)
207 grad.col(c) +=
fe->shape_grad(j, q->
point(i), c);
218 template<
unsigned int dim,
unsigned int spacedim>
227 template<
unsigned int dim,
unsigned int spacedim>
232 return data.shape_values[point_no][function_no];
236 template<
unsigned int dim,
unsigned int spacedim>
241 return data.shape_gradients[point_no][function_no];
245 template<
unsigned int dim,
unsigned int spacedim>
247 const unsigned int point_no,
248 const unsigned int comp)
const 257 template<
unsigned int dim,
unsigned int spacedim>
259 const unsigned int point_no,
260 const unsigned int comp)
const 269 template<
unsigned int dim,
unsigned int spacedim>
276 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
277 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
282 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
283 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
288 template<
unsigned int dim,
unsigned int spacedim>
296 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
297 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
300 for (
unsigned int c=0; c<spacedim; c++)
301 data.shape_values[i][j*spacedim+c] = fv_vec[c];
308 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
309 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
312 for (
unsigned int c=0; c<spacedim; c++)
313 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
319 template<
unsigned int dim,
unsigned int spacedim>
327 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
328 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
331 for (
unsigned int c=0; c<spacedim; c++)
332 data.shape_values[i][j*spacedim+c] = fv_vec[c];
339 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
340 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
343 for (
unsigned int c=0; c<spacedim; c++)
344 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
350 template<
unsigned int dim,
unsigned int spacedim>
358 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
359 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
362 for (
unsigned int c=0; c<spacedim; c++)
363 data.shape_values[i][j*spacedim+c] = fv_vec(c);
370 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
371 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
374 /
data.determinants[i];
375 for (
unsigned int c=0; c<spacedim; c++)
376 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
382 template<
unsigned int dim,
unsigned int spacedim>
390 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
391 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
394 for (
unsigned int c=0; c<spacedim*spacedim; c++)
395 data.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
402 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
403 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
406 for (
unsigned int c=0; c<spacedim*spacedim; c++)
407 data.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
413 template<
unsigned int dim,
unsigned int spacedim>
420 ASSERT_DBG(fe_sys !=
nullptr).error(
"Mixed system must be represented by FESystem.");
421 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
430 unsigned int n_space_components =
fe->n_space_components(spacedim);
436 unsigned int comp_offset = 0;
437 unsigned int shape_offset = 0;
438 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
440 unsigned int n_sub_space_components = fe_sys->
fe()[f]->n_space_components(spacedim);
443 for (
unsigned int i=0; i<fe_data.
n_points; i++)
444 for (
unsigned int n=0; n<fe_sys->
fe()[f]->n_dofs(); n++)
445 for (
unsigned int c=0; c<n_sub_space_components; c++)
446 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];
448 comp_offset += n_sub_space_components;
449 shape_offset += fe_sys->
fe()[f]->n_dofs()*n_space_components;
457 unsigned int comp_offset = 0;
458 unsigned int shape_offset = 0;
459 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
461 unsigned int n_sub_space_components = fe_sys->
fe()[f]->n_space_components(spacedim);
464 for (
unsigned int i=0; i<fe_data.
n_points; i++)
465 for (
unsigned int n=0; n<fe_sys->
fe()[f]->n_dofs(); n++)
466 for (
unsigned int c=0; c<n_sub_space_components; c++)
467 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];
469 comp_offset += n_sub_space_components;
470 shape_offset += fe_sys->
fe()[f]->n_dofs()*n_space_components;
477 template<
unsigned int dim,
unsigned int spacedim>
500 ASSERT(
false).error(
"Not implemented.");
511 template<
unsigned int dim,
unsigned int spacedim>
517 quadrature(&_quadrature)
519 this->
allocate(_mapping, _quadrature.
size(), _fe, _flags);
529 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
531 for (
auto fe_sub : fe->
fe())
538 template<
unsigned int dim,
unsigned int spacedim>
542 this->
data.present_cell = &cell;
545 this->
mapping->fill_fe_values(cell,
561 template<
unsigned int dim,
unsigned int spacedim>
569 this->
allocate(_mapping, _sub_quadrature.
size(), _fe, _flags);
571 for (
unsigned int sid = 0; sid < RefElement<dim>::n_sides; sid++)
573 for (
unsigned int pid = 0; pid < RefElement<dim>::n_side_permutations; pid++)
586 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
588 for (
auto fe_sub : fe->
fe())
595 template<
unsigned int dim,
unsigned int spacedim>
598 for (
unsigned int sid=0; sid<RefElement<dim>::n_sides; sid++)
600 for (
unsigned int pid=0; pid<RefElement<dim>::n_side_permutations; pid++)
609 template<
unsigned int dim,
unsigned int spacedim>
615 this->
data.present_cell = &cell;
621 this->
mapping->fill_fe_side_values(cell,
arma::vec::fixed< dim > point(const unsigned int i) const
Returns the ith quadrature point.
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
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...
MappingInternalData * mapping_data
Precomputed mapping data.
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.
FEValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags _flags)
Constructor.
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.
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.
Basic definitions of numerical quadrature rules.
const Quadrature< dim-1 > * sub_quadrature
Quadrature for the integration on the element sides.
Shape function gradients.
MappingInternalData * side_mapping_data[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
Quadrature< dim > side_quadrature[RefElement< dim >::n_sides][RefElement< dim >::n_side_permutations]
FEInternalData * init_fe_data(const Quadrature< dim > *q)
Precompute finite element data on reference element.
FEInternalData(unsigned int np, unsigned int nd)
Class Mapping calculates data related to the mapping of the reference cell to the actual cell...
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]
FESideValues(Mapping< dim, spacedim > &_mapping, Quadrature< dim-1 > &_sub_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Constructor.
ViewsCache views_cache_
Auxiliary storage of FEValuesViews accessors.
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.
const unsigned int size() const
Returns number of quadrature points.
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.
#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.