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 OLD_ASSERT(fe !=
nullptr,
"Mixed system must be represented by FESystem.");
122 template<
unsigned int dim,
unsigned int spacedim>
124 : mapping(NULL), quadrature(NULL), fe(NULL), mapping_data(NULL), fe_data(NULL)
130 template<
unsigned int dim,
unsigned int spacedim>
138 template<
unsigned int dim,
unsigned int spacedim>
149 ASSERT_DBG(_fe.
n_components() == spacedim*spacedim).error(
"FETensor must have spacedim*spacedim components.");
179 template<
unsigned int dim,
unsigned int spacedim>
184 arma::mat shape_values(
fe->n_dofs(),
fe->n_components());
185 for (
unsigned int i=0; i<q->
size(); i++)
187 for (
unsigned int j=0; j<
fe->n_dofs(); j++)
189 for (
unsigned int c=0; c<
fe->n_components(); c++)
190 shape_values(j,c) =
fe->shape_value(j, q->
point(i), c);
196 arma::mat grad(dim,
fe->n_components());
197 for (
unsigned int i=0; i<q->
size(); i++)
199 for (
unsigned int j=0; j<
fe->n_dofs(); j++)
202 for (
unsigned int c=0; c<
fe->n_components(); c++)
203 grad.col(c) +=
fe->shape_grad(j, q->
point(i), c);
214 template<
unsigned int dim,
unsigned int spacedim>
223 template<
unsigned int dim,
unsigned int spacedim>
228 return data.shape_values[point_no][function_no];
232 template<
unsigned int dim,
unsigned int spacedim>
237 return data.shape_gradients[point_no][function_no];
241 template<
unsigned int dim,
unsigned int spacedim>
243 const unsigned int point_no,
244 const unsigned int comp)
const 253 template<
unsigned int dim,
unsigned int spacedim>
255 const unsigned int point_no,
256 const unsigned int comp)
const 265 template<
unsigned int dim,
unsigned int spacedim>
272 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
273 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
278 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
279 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
284 template<
unsigned int dim,
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++)
297 data.shape_values[i][j*spacedim+c] = fv_vec[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++)
309 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
315 template<
unsigned int dim,
unsigned int spacedim>
323 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
324 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
327 for (
unsigned int c=0; c<spacedim; c++)
328 data.shape_values[i][j*spacedim+c] = fv_vec[c];
335 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
336 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
339 for (
unsigned int c=0; c<spacedim; c++)
340 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
346 template<
unsigned int dim,
unsigned int spacedim>
354 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
355 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
358 for (
unsigned int c=0; c<spacedim; c++)
359 data.shape_values[i][j*spacedim+c] = fv_vec(c);
366 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
367 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
370 /
data.determinants[i];
371 for (
unsigned int c=0; c<spacedim; c++)
372 data.shape_gradients[i][j*spacedim+c] = grads.col(c);
378 template<
unsigned int dim,
unsigned int spacedim>
386 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
387 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
390 for (
unsigned int c=0; c<spacedim*spacedim; c++)
391 data.shape_values[i][j*spacedim*spacedim+c] = fv_vec[c];
398 for (
unsigned int i = 0; i < fe_data.
n_points; i++)
399 for (
unsigned int j = 0; j < fe_data.
n_dofs; j++)
402 for (
unsigned int c=0; c<spacedim*spacedim; c++)
403 data.shape_gradients[i][j*spacedim*spacedim+c] = grads.col(c);
409 template<
unsigned int dim,
unsigned int spacedim>
416 ASSERT_DBG(fe_sys !=
nullptr).error(
"Mixed system must be represented by FESystem.");
417 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
430 unsigned int comp_offset = 0;
431 unsigned int shape_offset = 0;
432 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
435 for (
unsigned int i=0; i<fe_data.
n_points; i++)
436 for (
unsigned int n=0; n<fe_sys->
fe()[f]->n_dofs(); n++)
437 for (
unsigned int c=0; c<fe_sys->
fe()[f]->n_components(); c++)
440 comp_offset += fe_sys->
fe()[f]->n_components();
441 shape_offset += fe_sys->
fe()[f]->n_dofs()*fe_sys->
function_space_->n_components();
449 unsigned int comp_offset = 0;
450 unsigned int shape_offset = 0;
451 for (
unsigned int f=0; f<fe_sys->
fe().size(); f++)
454 for (
unsigned int i=0; i<fe_data.
n_points; i++)
455 for (
unsigned int n=0; n<fe_sys->
fe()[f]->n_dofs(); n++)
456 for (
unsigned int c=0; c<fe_sys->
fe()[f]->n_components(); c++)
459 comp_offset += fe_sys->
fe()[f]->n_components();
460 shape_offset += fe_sys->
fe()[f]->n_dofs()*fe_sys->
function_space_->n_components();
467 template<
unsigned int dim,
unsigned int spacedim>
490 ASSERT(
false).error(
"Not implemented.");
501 template<
unsigned int dim,
unsigned int spacedim>
508 this->
allocate(_mapping, _quadrature, _fe, _flags);
518 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
520 for (
auto fe_sub : fe->
fe())
527 template<
unsigned int dim,
unsigned int spacedim>
531 this->
data.present_cell = &cell;
534 this->
mapping->fill_fe_values(cell,
550 template<
unsigned int dim,
unsigned int spacedim>
559 this->
allocate(_mapping, *q, _fe, _flags);
561 for (
unsigned int sid = 0; sid < RefElement<dim>::n_sides; sid++)
563 for (
unsigned int pid = 0; pid < RefElement<dim>::n_side_permutations; pid++)
576 ASSERT_DBG(fe !=
nullptr).error(
"Mixed system must be represented by FESystem.");
578 for (
auto fe_sub : fe->
fe())
585 template<
unsigned int dim,
unsigned int spacedim>
588 for (
unsigned int sid=0; sid<RefElement<dim>::n_sides; sid++)
590 for (
unsigned int pid=0; pid<RefElement<dim>::n_side_permutations; pid++)
603 template<
unsigned int dim,
unsigned int spacedim>
609 this->
data.present_cell = &cell;
614 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.
std::vector< unsigned int > get_tensor_components() const
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...
void allocate(Mapping< dim, spacedim > &_mapping, Quadrature< dim > &_quadrature, FiniteElement< dim > &_fe, UpdateFlags flags)
Allocates space for computed data.
MappingInternalData * mapping_data
Precomputed mapping data.
void initialize(FEValuesBase &fv)
UpdateFlags update_each(UpdateFlags flags)
Determine quantities to be recomputed on each cell.
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...
std::vector< unsigned int > get_scalar_components() const
FEInternalData * fe_data
Precomputed finite element data.
FiniteElement< dim > * get_fe() const
Returns the finite element in use.
void allocate(unsigned int size, UpdateFlags flags, unsigned n_comp)
Resize the data arrays.
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.
Transformed quadrature points.
std::vector< unsigned int > get_vector_components() const
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)
Quadrature< dim > * quadrature
The quadrature rule used to calculate integrals.
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...
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.
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.
std::shared_ptr< FunctionSpace > function_space_
Function space defining the FE.
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.
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...
const arma::vec::fixed< dim > & point(const unsigned int i) const
Returns the ith quadrature point.
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.