21 #ifndef PATCH_FE_VALUES_HH_
22 #define PATCH_FE_VALUES_HH_
50 template <
class ValueType>
70 template <
class ValueType>
95 template <
class ValueType>
114 unsigned int begin,
unsigned int begin_side,
unsigned int n_dofs_bulk,
unsigned int n_dofs_side,
unsigned int join_idx)
174 template<
unsigned int dim>
183 template<
unsigned int FE_dim>
186 if (fe_sys !=
nullptr) {
187 return fe_sys->
fe()[component_idx];
189 ASSERT_EQ(component_idx, 0).warning(
"Non-zero component_idx can only be used for FESystem.");
201 template<
unsigned int FE_dim>
205 arma::mat shape_values(fe->n_dofs(), fe->n_components());
206 for (
unsigned int i=0; i<q->
size(); i++)
208 for (
unsigned int j=0; j<fe->n_dofs(); j++)
210 for (
unsigned int c=0; c<fe->n_components(); c++)
211 shape_values(j,c) = fe->shape_value(j, q->
point<FE_dim>(i), c);
213 ref_shape_vals[i][j] = trans(shape_values.row(j));
217 return ref_shape_vals;
228 template<
unsigned int FE_dim>
232 arma::mat shape_values(fe->n_dofs(), fe->n_components());
234 for (
unsigned int sid=0; sid<FE_dim+1; sid++) {
236 for (
unsigned int i=0; i<quad.size(); i++)
238 for (
unsigned int j=0; j<fe->n_dofs(); j++)
240 for (
unsigned int c=0; c<fe->n_components(); c++) {
241 shape_values(j,c) = fe->shape_value(j, quad.template point<FE_dim>(i), c);
244 ref_shape_vals[sid][i][j] = trans(shape_values.row(j));
249 return ref_shape_vals;
254 template<
unsigned int dim>
302 auto fe_component = this->
fe_comp(
fe_, component_idx);
303 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
309 for (
unsigned int j = 0; j < fe_component->n_dofs(); j++) {
310 shape_values[i][j] = ref_shape_vals[i][j][0];
317 uint begin = scalar_shape_bulk_op.result_row();
336 auto fe_component = this->
fe_comp(
fe_, component_idx);
337 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
343 bulk_reinit::ptop_scalar_shape_grads<dim>(operations, op_results, ref_shape_grads, scalar_shape_grads_op_idx);
346 uint begin = grad_scalar_shape_bulk_op.result_row();
367 for (
unsigned int i_pt=0; i_pt<q->
size(); i_pt++)
369 for (
unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
372 for (
unsigned int c=0; c<fe->n_components(); c++)
373 grad.col(c) += fe->shape_grad(i_dof, q->
point<dim>(i_pt), c);
375 ref_shape_grads[i_pt][i_dof] = grad;
379 return ref_shape_grads;
383 std::shared_ptr< FiniteElement<dim> >
fe_;
387 template<
unsigned int dim>
433 auto fe_component = this->
fe_comp(
fe_, component_idx);
434 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
442 for (
unsigned int s=0; s<dim+1; ++s)
444 for (
unsigned int j = 0; j < fe_component->n_dofs(); j++) {
445 shape_values[s][i][j] = ref_shape_vals[s][i][j][0];
452 uint begin = scalar_shape_bulk_op.result_row();
460 auto fe_component = this->
fe_comp(
fe_, component_idx);
461 ASSERT_EQ(fe_component->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of grad_scalar_shape accessor must be FEScalar!\n");
467 side_reinit::ptop_scalar_shape_grads<dim>(operations, op_results, el_table, ref_shape_grads, scalar_shape_grads_op_idx);
470 uint begin = grad_scalar_shape_side_op.result_row();
489 for (
unsigned int sid=0; sid<dim+1; sid++) {
491 for (
unsigned int i_pt=0; i_pt<quad.size(); i_pt++)
493 for (
unsigned int i_dof=0; i_dof<fe->n_dofs(); i_dof++)
496 for (
unsigned int c=0; c<fe->n_components(); c++)
497 grad.col(c) += fe->shape_grad(i_dof, quad.template point<dim>(i_pt), c);
499 ref_shape_grads[sid][i_pt][i_dof] = grad;
504 return ref_shape_grads;
508 std::shared_ptr< FiniteElement<dim> >
fe_;
512 template<
unsigned int dim>
529 ASSERT_EQ(fe_component_low->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
534 for (
unsigned int j = 0; j < fe_component_low->n_dofs(); j++) {
535 shape_values_bulk[i][j] = ref_shape_vals_bulk[i][j][0];
542 uint begin_bulk = grad_scalar_shape_bulk_op.result_row();
546 ASSERT_EQ(fe_component_high->fe_type(),
FEType::FEScalar).error(
"Type of FiniteElement of scalar_shape accessor must be FEScalar!\n");
553 for (
unsigned int s=0; s<dim+1; ++s)
555 for (
unsigned int j = 0; j < fe_component_high->n_dofs(); j++) {
556 shape_values_side[s][i][j] = ref_shape_vals_side[s][i][j][0];
563 uint begin_side = grad_scalar_shape_side_op.result_row();
566 begin_bulk, begin_side, fe_component_low->n_dofs(), fe_component_high->n_dofs(), 0) );
567 unsigned int end_idx = fe_component_low->n_dofs() + fe_component_high->n_dofs();
569 begin_bulk, begin_side, fe_component_low->n_dofs(), fe_component_high->n_dofs(), end_idx) );
598 template<
unsigned int spacedim = 3>
649 template<
unsigned int DIM>
661 template<
unsigned int DIM>
667 template<
unsigned int DIM>
683 template<
class MapType>
690 inline void set_shape_value(
unsigned int i_point,
unsigned int i_func_comp,
double val)
699 inline void set_shape_gradient(
unsigned int i_point,
unsigned int i_func_comp, arma::vec::fixed<spacedim> val)
712 inline double shape_value(
const unsigned int function_no,
const unsigned int point_no)
const
727 inline arma::vec::fixed<spacedim>
shape_grad(
const unsigned int function_no,
const unsigned int point_no)
const
893 dim_fe_side_vals_({DimPatchFEValues(0), DimPatchFEValues(0), DimPatchFEValues(0)}),
901 dim_fe_side_vals_({DimPatchFEValues(n_quad_points), DimPatchFEValues(n_quad_points), DimPatchFEValues(n_quad_points)}),
930 template<
unsigned int DIM>
934 if ( _quadrature.
dim() == DIM ) {
950 for (
unsigned int i=0; i<3; ++i) {
957 void reinit(std::array<PatchElementsList, 4> patch_elements)
959 for (
unsigned int i=0; i<3; ++i) {
968 for (
unsigned int i=0; i<3; ++i) {
977 template<
unsigned int dim>
979 ASSERT((dim>=0) && (dim<=3))(dim).error(
"Dimension must be 0, 1, 2 or 3.");
984 template<
unsigned int dim>
986 ASSERT((dim>0) && (dim<=3))(dim).error(
"Dimension must be 1, 2 or 3.");
991 template<
unsigned int dim>
993 ASSERT((dim>0) && (dim<=3))(dim).error(
"Dimension must be 1, 2 or 3.");
998 template<
unsigned int dim>
1000 ASSERT((dim>1) && (dim<=3))(dim).error(
"Dimension must be 2 or 3.");
1018 for (
uint i=0; i<3; ++i) {
1027 switch (cell.
dim()) {
1050 for (
unsigned int n=0; n<cell_side.
dim(); n++)
1051 for (
unsigned int c=0; c<spacedim; c++)
1052 side_coords(c,n) = (*cell_side.
side().
node(n))[c];
1056 switch (cell.
dim()) {
1087 void print(
bool points,
bool ints,
bool only_bulk=
true)
const {
1088 for (
uint i=0; i<3; ++i)
1091 for (
uint i=0; i<3; ++i)
1097 stream << endl <<
"Table of patch FE operations:" << endl;
1098 for (
uint i=0; i<3; ++i) {
1099 stream << std::setfill(
'-') << setw(100) <<
"" << endl;
1103 for (
uint i=0; i<3; ++i) {
1104 stream << std::setfill(
'-') << setw(100) <<
"" << endl;
1108 stream << std::setfill(
'=') << setw(100) <<
"" << endl;
1124 template <
class ValueType>
1126 template <
class ValueType>
1128 template <
class ValueType>
1133 template <
class ValueType>
1136 return patch_point_vals_.scalar_val(
begin_, value_cache_idx);
1142 return patch_point_vals_.vector_val(
begin_, value_cache_idx);
1148 return patch_point_vals_.tensor_val(
begin_, value_cache_idx);
1151 template <
class ValueType>
1154 return patch_point_vals_.scalar_val(
begin_, value_cache_idx);
1160 return patch_point_vals_.vector_val(
begin_, value_cache_idx);
1166 return patch_point_vals_.tensor_val(
begin_, value_cache_idx);
1169 template <
class ValueType>
1172 return patch_point_vals_.scalar_val(
begin_+shape_idx, value_cache_idx);
1178 return patch_point_vals_.vector_val(
begin_+3*shape_idx, value_cache_idx);
1183 Tensor tens; tens.zeros();
1187 template <
class ValueType>
1190 return patch_point_vals_.scalar_val(
begin_+shape_idx, value_cache_idx);
1196 return patch_point_vals_.vector_val(
begin_+3*shape_idx, value_cache_idx);
1201 Tensor tens; tens.zeros();
1206 template <
class ValueType>
1218 Vector vect; vect.zeros();
1224 Tensor tens; tens.zeros();
1228 template <
class ValueType>
1240 Vector vect; vect.zeros();
1246 Tensor tens; tens.zeros();
#define ASSERT_LT(a, b)
Definition of comparative assert macro (Less Than) only for debug mode.
#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()=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)
ElQ(PatchPointValues< 3 > &patch_point_vals, unsigned int begin)
Constructor.
int element_eval_point(unsigned int i_elem_in_cache, unsigned int i_eval_point) const
Class for computation of data on cell and side.
Structure for storing the precomputed finite element data.
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.
unsigned int begin_
Index of the first component of the Quantity. Size is given by ValueType.
PatchPointValues< 3 > & patch_point_vals_
Reference to PatchPointValues.
unsigned int n_dofs_
Number of DOFs.
ValueType operator()(FMT_UNUSED unsigned int shape_idx, FMT_UNUSED const BulkPoint &point)
FeQ(PatchPointValues< 3 > &patch_point_vals, unsigned int begin, unsigned int n_dofs)
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 begin_side_
Index of the first component of the side Quantity. Size is given by ValueType.
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.
unsigned int begin_
Index of the first component of the bulk Quantity. Size is given by ValueType.
PatchPointValues< 3 > * patch_point_vals_bulk_
Pointer to bulk PatchPointValues.
JoinShapeAccessor(PatchPointValues< 3 > *patch_point_vals_bulk, PatchPointValues< 3 > *patch_point_vals_side, unsigned int begin, unsigned int begin_side, unsigned int n_dofs_bulk, unsigned int n_dofs_side, unsigned int join_idx)
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)
Subobject holds FE data of one dimension (0,1,2,3)
arma::vec::fixed< spacedim > normal_vector(const SidePoint &p)
Returns the normal vector to a side at given quadrature point.
void fill_data_specialized(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....
unsigned int max_n_elem_
Maximal number of elements on patch.
unsigned int n_dofs_
Number of finite element dofs.
std::shared_ptr< FEInternalData > init_fe_data(const FiniteElement< DIM > &fe, const Quadrature &q)
Precompute finite element data on reference element.
unsigned int used_size_
Number of elements / sides on patch. Must be less or equal to size of element_data vector.
double shape_value(const unsigned int function_no, const BulkPoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
double shape_value(const unsigned int function_no, const unsigned int point_no) const
Return the value of the function_no-th shape function at the point_no-th quadrature point.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const BulkPoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
void resize(unsigned int max_size)
Set size of ElementFEData. Important: Use only during the initialization of FESystem !
std::vector< ElementFEData > element_data_
Data of elements / sides on patch.
double JxW(const BulkPoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
std::map< unsigned int, unsigned int > element_patch_map_
Map of element patch indexes to element_data_.
DimPatchFEValues(unsigned int max_size=0)
Constructor.
std::vector< unsigned int > fe_sys_n_components_
Numbers of components of FESystem sub-elements in reference space.
void allocate(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags flags)
Allocates space for computed data.
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.
void set_shape_gradient(unsigned int i_point, unsigned int i_func_comp, arma::vec::fixed< spacedim > val)
void reinit(PatchElementsList patch_elements)
std::vector< unsigned int > fe_sys_n_space_components_
Numbers of components of FESystem sub-elements in real space.
void set_shape_value(unsigned int i_point, unsigned int i_func_comp, double val)
void initialize(Quadrature &_quadrature, FiniteElement< DIM > &_fe, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
unsigned int n_components_
Number of components of the FE.
MeshObjectType object_type_
Distinguishes using of PatchFEValues for storing data of elements or sides.
double JxW(const SidePoint &p)
Return the product of Jacobian determinant and the quadrature weight at given quadrature point.
unsigned int n_points_
Number of integration points.
UpdateFlags update_flags
Flags that indicate which finite element quantities are to be computed.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const unsigned int point_no) const
Return the gradient of the function_no-th shape function at the point_no-th quadrature point.
unsigned int max_size() const
std::vector< shared_ptr< FEInternalData > > side_fe_data_
Precomputed FE data (shape functions on reference element) for all side quadrature points.
std::shared_ptr< FEInternalData > fe_data_
Precomputed finite element data.
unsigned int used_size() const
FEType fe_type_
Type of finite element (scalar, vector, tensor).
unsigned int patch_data_idx_
Patch index of processed element / side.
unsigned int dim_
Dimension of reference space.
std::vector< std::vector< unsigned int > > fe_sys_dofs_
Dof indices of FESystem sub-elements.
std::vector< PatchFEValues< spacedim >::DimPatchFEValues > fe_values_vec
Vector of FEValues for sub-elements of FESystem.
double shape_value(const unsigned int function_no, const SidePoint &p) const
Return the value of the function_no-th shape function at the p quadrature point.
arma::vec::fixed< spacedim > shape_grad(const unsigned int function_no, const SidePoint &p) const
Return the gradient of the function_no-th shape function at the p quadrature point.
std::vector< std::vector< arma::vec::fixed< spacedim > > > shape_gradients_
std::shared_ptr< ElementValues< spacedim > > elm_values_
Auxiliary object for calculation of element-dependent data.
std::vector< std::vector< double > > shape_values_
Shape functions evaluated at the quadrature points.
std::array< DimPatchFEValues, 3 > dim_fe_vals_
Sub objects of dimensions 1,2,3.
uint register_element(DHCellAccessor cell, uint element_patch_idx)
Register element to patch_point_vals_ table by dimension of element.
std::array< DimPatchFEValues, 3 > dim_fe_side_vals_
void print(bool points, bool ints, bool only_bulk=true) const
Temporary development method.
~PatchFEValues()
Destructor.
std::array< FeBulk::PatchPointValues< spacedim >, 3 > patch_point_vals_bulk_
uint register_bulk_point(DHCellAccessor cell, uint elem_table_row, uint value_patch_idx)
Register bulk point to patch_point_vals_ table by dimension of element.
uint register_side(DHCellSide cell_side)
Register side to patch_point_vals_ table by dimension of side.
DimPointTable dim_point_table_
void initialize(Quadrature &_quadrature, UpdateFlags _flags)
Initialize structures and calculates cell-independent data.
JoinValues< dim > join_values()
Return JoinValue object of dimension given by template parameter.
void reinit(std::array< PatchElementsList, 4 > patch_elements)
Reinit data - old data storing, temporary.
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.
void print_operations(ostream &stream) const
Temporary development method.
void resize_tables(std::vector< std::vector< uint > > dim_sizes)
Resize tables of patch_point_vals_.
std::array< FeSide::PatchPointValues< spacedim >, 3 > patch_point_vals_side_
DimPointTable & dim_point_table(unsigned int new_size)
Resize dim_point_table_ if actual size is less than new_size and return reference.
MixedPtr< FiniteElement > fe_
Mixed of shared pointers of FiniteElement object.
SideValues< dim > side_values()
Return SideValue object of dimension given by template parameter.
void reinit_patch()
Reinit data.
PatchFEValues(unsigned int n_quad_points, unsigned int quad_order, MixedPtr< FiniteElement > fe)
unsigned int n_dofs() const
Returns the number of shape functions.
uint register_side_point(DHCellSide cell_side, uint elem_table_row, uint value_patch_idx)
Register side point to patch_point_vals_ table by dimension of side.
void reset()
Reset PatchpointValues structures.
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Return quadrature.
uint dim() const
Getter of dim_.
ElOp< spacedim > & make_fe_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, uint n_dofs)
Add accessor to operations_ vector.
Scalar scalar_val(uint result_row, uint point_idx) 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 ...
virtual unsigned int side_idx() const =0
Intermediate step in implementation of PatcFEValues.
unsigned int eval_point_idx() const
Return index in EvalPoints object.
unsigned int local_point_idx() const
Return local index in quadrature. Temporary method - intermediate step in implementation of PatcFEVal...
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.
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< std::vector< double > > shape_values, uint scalar_shape_op_idx)
static void ptop_scalar_shape(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &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.
UpdateFlags
Enum type UpdateFlags indicates which quantities are to be recomputed on each finite element cell.