21 #ifndef PATCH_FE_VALUES_HH_
22 #define PATCH_FE_VALUES_HH_
52 template <
class ValueType>
72 template <
class ValueType>
97 template <
class ValueType>
116 unsigned int n_dofs_side,
unsigned int op_idx_bulk,
unsigned int op_idx_side,
unsigned int join_idx)
176 template<
unsigned int dim>
185 template<
unsigned int FE_dim>
188 if (fe_sys !=
nullptr) {
189 return fe_sys->
fe()[component_idx];
191 ASSERT_EQ(component_idx, 0).warning(
"Non-zero component_idx can only be used for FESystem.");
203 template<
unsigned int FE_dim>
207 arma::mat shape_values(fe->n_dofs(), fe->n_components());
208 for (
unsigned int i=0; i<q->
size(); i++)
210 for (
unsigned int j=0; j<fe->n_dofs(); j++)
212 for (
unsigned int c=0; c<fe->n_components(); c++)
213 shape_values(j,c) = fe->shape_value(j, q->
point<FE_dim>(i), c);
215 ref_shape_vals[i][j] = trans(shape_values.row(j));
219 return ref_shape_vals;
230 template<
unsigned int FE_dim>
234 arma::mat shape_values(fe->n_dofs(), fe->n_components());
236 for (
unsigned int sid=0; sid<FE_dim+1; sid++) {
238 for (
unsigned int i=0; i<quad.size(); i++)
240 for (
unsigned int j=0; j<fe->n_dofs(); j++)
242 for (
unsigned int c=0; c<fe->n_components(); c++) {
243 shape_values(j,c) = fe->shape_value(j, quad.template point<FE_dim>(i), c);
246 ref_shape_vals[sid][i][j] = trans(shape_values.row(j));
251 return ref_shape_vals;
256 template<
unsigned int dim>
301 auto fe_component = this->
fe_comp(
fe_, component_idx);
302 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
308 for (
unsigned int j = 0; j < fe_component->n_dofs(); j++) {
309 shape_values[j][i] = ref_shape_vals[i][j][0];
335 auto fe_component = this->
fe_comp(
fe_, component_idx);
336 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
342 bulk_reinit::ptop_scalar_shape_grads<dim>(operations, ref_shape_grads, scalar_shape_grads_op_idx);
366 for (
unsigned int i_pt=0; i_pt<q->
size(); i_pt++)
368 for (
unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
371 for (
unsigned int c=0; c<fe->n_components(); c++)
372 grad.col(c) += fe->shape_grad(i_dof, q->
point<dim>(i_pt), c);
374 ref_shape_grads[i_pt][i_dof] = grad;
378 return ref_shape_grads;
382 std::shared_ptr< FiniteElement<dim> >
fe_;
386 template<
unsigned int dim>
428 auto fe_component = this->
fe_comp(
fe_, component_idx);
429 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
437 for (
unsigned int s=0; s<dim+1; ++s)
439 for (
unsigned int j = 0; j < fe_component->n_dofs(); j++) {
440 shape_values[s][i][j] = ref_shape_vals[s][i][j][0];
455 auto fe_component = this->
fe_comp(
fe_, component_idx);
456 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
462 side_reinit::ptop_scalar_shape_grads<dim>(operations, el_table, ref_shape_grads, scalar_shape_grads_op_idx);
484 for (
unsigned int sid=0; sid<dim+1; sid++) {
486 for (
unsigned int i_pt=0; i_pt<quad.size(); i_pt++)
488 for (
unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
491 for (
unsigned int c=0; c<fe->n_components(); c++)
492 grad.col(c) += fe->shape_grad(i_dof, quad.template point<dim>(i_pt), c);
494 ref_shape_grads[sid][i_pt][i_dof] = grad;
499 return ref_shape_grads;
503 std::shared_ptr< FiniteElement<dim> >
fe_;
507 template<
unsigned int dim>
524 ASSERT_EQ(fe_component_low->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
529 for (
unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
530 shape_values_bulk[i][j] = ref_shape_vals_bulk[i][j][0];
541 ASSERT_EQ(fe_component_high->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
548 for (
unsigned int s=0; s<dim+1; ++s)
550 for (
unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
551 shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j][0];
561 fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, 0) );
562 unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
564 fe_component_low->n_dofs(), fe_component_high->n_dofs(), op_idx_bulk, op_idx_side, end_idx) );
593 template<
unsigned int spacedim = 3>
640 used_quads_[0] =
false; used_quads_[1] =
false;
644 : asm_arena_(1024 * 1024, 256),
645 patch_arena_(nullptr),
654 used_quads_[0] =
false; used_quads_[1] =
false;
661 if (patch_arena_!=
nullptr)
667 if (is_bulk)
return patch_point_vals_bulk_[dim-1].get_quadrature();
668 else return patch_point_vals_side_[dim-1].get_quadrature();
678 template<
unsigned int DIM>
681 if ( _quadrature.
dim() == DIM ) {
682 used_quads_[0] =
true;
683 patch_point_vals_bulk_[DIM-1].initialize();
685 used_quads_[1] =
true;
686 patch_point_vals_side_[DIM-1].initialize();
692 patch_arena_ = asm_arena_.get_child_arena();
693 for (
unsigned int i=0; i<3; ++i) {
694 if (used_quads_[0]) patch_point_vals_bulk_[i].init_finalize(patch_arena_);
695 if (used_quads_[1]) patch_point_vals_side_[i].init_finalize(patch_arena_);
702 for (
unsigned int i=0; i<3; ++i) {
703 if (used_quads_[0]) patch_point_vals_bulk_[i].reset();
704 if (used_quads_[1]) patch_point_vals_side_[i].reset();
706 patch_arena_->reset();
712 for (
unsigned int i=0; i<3; ++i) {
713 if (used_quads_[0]) patch_point_vals_bulk_[i].reinit_patch();
714 if (used_quads_[1]) patch_point_vals_side_[i].reinit_patch();
721 template<
unsigned int dim>
723 ASSERT((dim>=0) && (dim<=3))(dim).error(
"Dimension must be 0, 1, 2 or 3.");
728 template<
unsigned int dim>
730 ASSERT((dim>0) && (dim<=3))(dim).error(
"Dimension must be 1, 2 or 3.");
735 template<
unsigned int dim>
737 ASSERT((dim>0) && (dim<=3))(dim).error(
"Dimension must be 1, 2 or 3.");
742 template<
unsigned int dim>
745 return JoinValues<dim>(&patch_point_vals_bulk_[dim-2], &patch_point_vals_side_[dim-1], fe_);
752 for (
uint i=0; i<3; ++i) {
753 if (used_quads_[0]) patch_point_vals_bulk_[i].resize_tables(table_sizes.
elem_sizes_[0][i], table_sizes.
point_sizes_[0][i]);
754 if (used_quads_[1]) patch_point_vals_side_[i].resize_tables(table_sizes.
elem_sizes_[1][i], table_sizes.
point_sizes_[1][i]);
761 switch (cell.
dim()) {
764 return patch_point_vals_bulk_[0].register_element(coords, element_patch_idx);
768 return patch_point_vals_bulk_[1].register_element(coords, element_patch_idx);
772 return patch_point_vals_bulk_[2].register_element(coords, element_patch_idx);
784 for (
unsigned int n=0; n<cell_side.
dim(); n++)
785 for (
unsigned int c=0; c<spacedim; c++)
786 side_coords(c,n) = (*cell_side.
side().
node(n))[c];
790 switch (cell.
dim()) {
793 return patch_point_vals_side_[0].register_side(elm_coords, side_coords, cell_side.
side_idx());
797 return patch_point_vals_side_[1].register_side(elm_coords, side_coords, cell_side.
side_idx());
801 return patch_point_vals_side_[2].register_side(elm_coords, side_coords, cell_side.
side_idx());
812 return patch_point_vals_bulk_[cell.
dim()-1].register_bulk_point(elem_table_row, value_patch_idx, cell.
elm_idx(), i_point_on_elem);
817 return patch_point_vals_side_[cell_side.
dim()-1].register_side_point(elem_table_row, value_patch_idx, cell_side.
elem_idx(),
818 cell_side.
side_idx(), i_point_on_side);
823 stream << endl <<
"Table of patch FE data:" << endl;
824 for (
uint i=0; i<3; ++i) {
825 stream << std::setfill(
'-') << setw(100) <<
"" << endl;
826 stream <<
"Bulk, dimension " << (i+1) << endl;
827 patch_point_vals_bulk_[i].print_data_tables(stream, points, ints);
830 for (
uint i=0; i<3; ++i) {
831 stream << std::setfill(
'-') << setw(100) <<
"" << endl;
832 stream <<
"Side, dimension " << (i+1) << endl;
833 patch_point_vals_side_[i].print_data_tables(stream, points, ints);
835 stream << std::setfill(
'=') << setw(100) <<
"" << endl;
840 stream << endl <<
"Table of patch FE operations:" << endl;
841 for (
uint i=0; i<3; ++i) {
842 stream << std::setfill(
'-') << setw(100) <<
"" << endl;
843 stream <<
"Bulk, dimension " << (i+1) << endl;
844 patch_point_vals_bulk_[i].print_operations(stream, 0);
846 for (
uint i=0; i<3; ++i) {
847 stream << std::setfill(
'-') << setw(100) <<
"" << endl;
848 stream <<
"Side, dimension " << (i+1) << endl;
849 patch_point_vals_side_[i].print_operations(stream, 1);
851 stream << std::setfill(
'=') << setw(100) <<
"" << endl;
863 template <
class ValueType>
865 template <
class ValueType>
867 template <
class ValueType>
872 template <
class ValueType>
875 return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
881 return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
887 return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
890 template <
class ValueType>
893 return patch_point_vals_.scalar_value(op_idx_, value_cache_idx);
899 return patch_point_vals_.vector_value(op_idx_, value_cache_idx);
905 return patch_point_vals_.tensor_value(op_idx_, value_cache_idx);
908 template <
class ValueType>
911 return patch_point_vals_.scalar_value(op_idx_, value_cache_idx, shape_idx);
917 return patch_point_vals_.vector_value(op_idx_, value_cache_idx, shape_idx);
923 return patch_point_vals_.tensor_value(op_idx_, value_cache_idx, shape_idx);
926 template <
class ValueType>
929 return patch_point_vals_.scalar_value(op_idx_, value_cache_idx, shape_idx);
935 return patch_point_vals_.vector_value(op_idx_, value_cache_idx, shape_idx);
940 unsigned int value_cache_idx = point.elm_cache_map()->element_eval_point(point.elem_patch_idx(), point.eval_point_idx());
941 return patch_point_vals_.tensor_value(op_idx_, value_cache_idx, shape_idx);
945 template <
class ValueType>
957 Vector vect; vect.zeros();
963 Tensor tens; tens.zeros();
967 template <
class ValueType>
979 Vector vect; vect.zeros();
985 Tensor tens; tens.zeros();
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
std::vector< std::vector< std::vector< arma::vec > > > ref_shape_values_side(Quadrature *q, std::shared_ptr< FiniteElement< FE_dim >> fe)
Precomputed values of basis functions at the side quadrature points.
std::shared_ptr< FiniteElement< FE_dim > > fe_comp(std::shared_ptr< FiniteElement< FE_dim > > fe, uint component_idx)
Return FiniteElement of component_idx for FESystem or fe for other types.
std::vector< std::vector< arma::vec > > ref_shape_values_bulk(Quadrature *q, std::shared_ptr< FiniteElement< FE_dim >> fe)
Precomputed values of basis functions at the bulk quadrature points.
Base point accessor class.
const ElementCacheMap * elm_cache_map() const
unsigned int elem_patch_idx() const
unsigned int eval_point_idx() const
Return index in EvalPoints object.
ElQ< Scalar > JxW()
Register the product of Jacobian determinant and the quadrature weight at bulk quadrature points.
PatchPointValues< 3 > & patch_point_vals_
ElQ< Vector > coords()
Create bulk accessor of coords entity.
std::vector< std::vector< arma::mat > > ref_shape_gradients(std::shared_ptr< FiniteElement< dim >> fe)
Precomputed gradients of basis functions at the quadrature points.
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
std::shared_ptr< FiniteElement< dim > > fe_
FeQ< Vector > grad_scalar_shape(uint component_idx=0)
Return the value of the function_no-th gradient shape function at the p bulk quadrature point.
FeQ< Scalar > scalar_shape(uint component_idx=0)
Return the value of the function_no-th shape function at the p bulk quadrature point.
BulkValues(PatchPointValues< 3 > &patch_point_vals, MixedPtr< FiniteElement > fe)
Constructor.
Cell accessor allow iterate over DOF handler cells.
unsigned int dim() const
Return dimension of element appropriate to cell.
unsigned int elm_idx() const
Return serial idx to element of loc_ele_idx_.
ElementAccessor< 3 > elm() const
Return ElementAccessor to element of loc_ele_idx_.
Side accessor allows to iterate over sides of DOF handler cell.
unsigned int elem_idx() const
Side side() const
Return Side of given cell and side_idx.
const DHCellAccessor & cell() const
Return DHCellAccessor appropriate to the side.
unsigned int dim() const
Return dimension of element appropriate to the side.
unsigned int side_idx() const
ElQ(PatchPointValues< 3 > &patch_point_vals, unsigned int op_idx)
Constructor.
unsigned int op_idx_
Index of operation in patch_point_vals_.operations vector.
ElQ()=delete
Forbidden default constructor.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
ValueType operator()(FMT_UNUSED const SidePoint &point)
ValueType operator()(FMT_UNUSED const BulkPoint &point)
int element_eval_point(unsigned int i_elem_in_cache, unsigned int i_eval_point) const
Compound finite element on dim dimensional simplex.
const std::vector< std::shared_ptr< FiniteElement< dim > > > & fe() const
Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
FeQ(PatchPointValues< 3 > &patch_point_vals, unsigned int op_idx, unsigned int n_dofs)
unsigned int n_dofs_
Number of DOFs.
ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point)
unsigned int op_idx_
Index of operation in patch_point_vals_.operations vector.
FeQ()=delete
Forbidden default constructor.
ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const SidePoint &point)
Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
unsigned int n_dofs_both() const
unsigned int local_idx() const
Return local index of DOF (on low / high-dim) - should be private method.
bool operator==(const JoinShapeAccessor< ValueType > &other) const
Comparison of accessors.
JoinShapeAccessor(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int op_idx_bulk, unsigned int op_idx_side, unsigned int join_idx)
PatchPointValues< 3 > * patch_point_vals_bulk_
Pointer to bulk PatchPointValues.
unsigned int op_idx_side_
Index of operation in patch_point_vals_side_.operations vector.
unsigned int op_idx_bulk_
Index of operation in patch_point_vals_bulk_.operations vector.
unsigned int n_dofs_high() const
JoinShapeAccessor()
Default constructor.
void inc()
Iterates to next item.
unsigned int join_idx_
Index of processed DOF.
unsigned int n_dofs_high_
Number of DOFs on high-dim element.
unsigned int join_idx() const
Return global index of DOF.
unsigned int n_dofs_low() const
ValueType operator()(const BulkPoint &point)
unsigned int n_dofs_low_
Number of DOFs on low-dim element.
PatchPointValues< 3 > * patch_point_vals_side_
Pointer to side PatchPointValues.
JoinValues(FMT_UNUSED PatchPointValues< 3 > *patch_point_vals_bulk, FMT_UNUSED PatchPointValues< 3 > *patch_point_vals_side, FMT_UNUSED MixedPtr< FiniteElement > fe)
Constructor.
Range< JoinShapeAccessor< Scalar > > scalar_join_shape(FMT_UNUSED uint component_idx=0)
std::shared_ptr< FiniteElement< dim > > fe_high_dim_
JoinValues(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, MixedPtr< FiniteElement > fe)
Constructor.
PatchPointValues< 3 > * patch_point_vals_bulk_
Range< JoinShapeAccessor< Scalar > > scalar_join_shape(uint component_idx=0)
std::shared_ptr< FiniteElement< dim-1 > > fe_low_dim_
PatchPointValues< 3 > * patch_point_vals_side_
static ElementMap element_map(ElementAccessor< 3 > elm)
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
PatchArena * patch_arena_
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx, uint i_point_on_elem)
Register bulk point to patch_point_vals_ table by dimension of element.
~PatchFEValues()
Destructor.
std::array< FeBulk::PatchPointValues< spacedim >, 3 > patch_point_vals_bulk_
Sub objects of bulk data of dimensions 1,2,3.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
void init_finalize()
Finalize initialization, creates child (patch) arena and passes it to PatchPointValue objects.
void print_data_tables(ostream &stream, bool points, bool ints, bool only_bulk=true) const
Temporary development method.
JoinValues< dim > join_values()
Return JoinValue object of dimension given by template parameter.
Quadrature * get_quadrature(uint dim, bool is_bulk) const
Return bulk or side quadrature of given dimension.
BulkValues< dim > bulk_values()
Return BulkValue object of dimension given by template parameter.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx, uint i_point_on_side)
Register side point to patch_point_vals_ table by dimension of side.
void initialize(Quadrature &_quadrature)
Initialize structures and calculates cell-independent data.
void print_operations(ostream &stream) const
Temporary development method.
std::array< FeSide::PatchPointValues< spacedim >, 3 > patch_point_vals_side_
Sub objects of side data of dimensions 1,2,3.
MixedPtr< FiniteElement > fe_
Mixed of shared pointers of FiniteElement object.
SideValues< dim > side_values()
Return SideValue object of dimension given by template parameter.
void resize_tables(TableSizes table_sizes)
Resize tables of patch_point_vals_.
PatchFEValues(unsigned int quad_order, MixedPtr< FiniteElement > fe)
void reinit_patch()
Reinit data.
unsigned int n_dofs() const
Returns the number of shape functions.
void reset()
Reset PatchpointValues structures.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs, OpSizeType size_type=pointOp)
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Getter for quadrature.
uint dim() const
Getter for dim_.
Scalar scalar_value(uint op_idx, uint point_idx, uint i_dof=0) const
Base class for quadrature rules on simplices in arbitrary dimensions.
Quadrature make_from_side(unsigned int sid) const
unsigned int size() const
Returns number of quadrature points.
Armor::ArmaVec< double, point_dim > point(unsigned int i) const
Returns the ith quadrature point.
General point a+ side_begin_ + ccessor allow iterate over quadrature points of given side defined in ...
unsigned int eval_point_idx() const
Return index in EvalPoints object.
SideValues(PatchPointValues< 3 > &patch_point_vals, MixedPtr< FiniteElement > fe)
Constructor.
ElQ< Vector > normal_vector()
Register the normal vector to a side at side quadrature points.
std::shared_ptr< FiniteElement< dim > > fe_
ElQ< Scalar > JxW()
Same as BulkValues::JxW but register at side quadrature points.
ElQ< Vector > coords()
Create side accessor of coords entity.
FeQ< Scalar > scalar_shape(uint component_idx=0)
Same as BulkValues::scalar_shape but register at side quadrature points.
FeQ< Vector > grad_scalar_shape(uint component_idx=0)
Same as BulkValues::grad_scalar_shape but register at side quadrature points.
std::vector< std::vector< std::vector< arma::mat > > > ref_shape_gradients(std::shared_ptr< FiniteElement< dim >> fe)
Precomputed gradients of basis functions at the quadrature points.
PatchPointValues< 3 > & patch_point_vals_
ElQ< Scalar > determinant()
Create bulk accessor of jac determinant entity.
NodeAccessor< 3 > node(unsigned int i) const
Returns node for given local index i on the side.
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
Class FESystem for compound finite elements.
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Iter< Object > make_iter(Object obj)
Class MappingP1 implements the affine transformation of the unit cell onto the actual cell.
ArmaMat< double, N, M > mat
@ opCoords
operations evaluated on quadrature points
@ opInvJac
inverse Jacobian
@ opJxW
JxW value of quadrature point.
@ opNormalVec
normal vector of quadrature point
std::vector< std::array< uint, 3 > > DimPointTable
Holds triplet (dim; bulk/side; idx of point in subtable)
Store finite element data on the actual patch such as shape function values, gradients,...
Definitions of particular quadrature rules on simplices.
Class RefElement defines numbering of vertices, sides, calculation of normal vectors etc.
Struct for pre-computing number of elements, sides, bulk points and side points on each dimension.
std::vector< std::vector< uint > > point_sizes_
void reset()
Set all values to zero.
std::vector< std::vector< uint > > elem_sizes_
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< double > > shape_values, uint scalar_shape_op_idx)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< double > > > shape_values, uint scalar_shape_op_idx)
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.