21 #ifndef PATCH_POINT_VALUES_HH_
22 #define PATCH_POINT_VALUES_HH_
24 #include <Eigen/Dense>
35 template<
unsigned int spacedim>
class ElOp;
133 template<
unsigned int spacedim = 3>
197 elOp.allocate_result(sizes[elOp.size_type()], *
patch_arena_);
217 for (
uint i_col=0; i_col<coords.n_cols; ++i_col)
218 for (
uint i_row=0; i_row<coords.n_rows; ++i_row) {
219 coords_mat(i_row, i_col)(i_elem) = coords(i_row, i_col);
238 for (
uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
239 for (
uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
240 coords_mat(i_row, i_col)(i_elem) = elm_coords(i_row, i_col);
248 for (
uint i_col=0; i_col<side_coords.n_cols; ++i_col)
249 for (
uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
250 coords_mat(i_row, i_col)(i_elem) = side_coords(i_row, i_col);
268 uint point_pos = i_point_on_elem *
n_elems_ + elem_table_row;
287 uint point_pos = i_point_on_side *
n_elems_ + elem_table_row;
333 ElOp<spacedim> op_accessor(this->
dim_, shape, reinit_f, size_type, input_ops_vec, n_dofs);
371 const auto &op_matrix =
operations_[op_idx].result_matrix();
373 for (
uint i=0; i<3; ++i)
374 val(i) = op_matrix(i)(op_matrix_idx);
387 const auto &op_matrix =
operations_[op_idx].result_matrix();
389 for (
uint i=0; i<3; ++i)
390 for (
uint j=0; j<3; ++j)
391 val(i,j) = op_matrix(i,j)(op_matrix_idx);
414 stream <<
"Point vals: " << std::endl;
416 const auto &
mat = op.result_matrix();
417 for (
uint i_mat=0; i_mat<
mat.rows()*
mat.cols(); ++i_mat) {
418 if (
mat(i_mat).data_size()==0) stream <<
"<empty>";
420 const double *vals =
mat(i_mat).data_ptr();
421 for (
size_t i_val=0; i_val<
mat(i_mat).data_size(); ++i_val)
422 stream << vals[i_val] <<
" ";
426 stream <<
" --- end of operation ---" << std::endl;
432 if (
int_table_(i_row).data_size()==0) stream <<
"<empty>";
435 for (
size_t i_val=0; i_val<
int_table_(i_row).data_size(); ++i_val)
436 stream << vals[i_val] <<
" ";
453 {
"el_coords",
"jacobian",
"inv_jac",
"jac_det",
"pt_coords",
"weights",
"JxW",
"",
"",
"",
"",
"" },
454 {
"el_coords",
"el_jac",
"el_inv_jac",
"side_coords",
"side_jac",
"side_jac_det",
"pt_coords",
"weights",
"JxW",
"normal_vec",
"",
"",
"",
"",
"" }
456 stream << std::setfill(
' ') <<
" Operation" << setw(12) <<
"" <<
"Shape" << setw(2) <<
""
457 <<
"n DOFs" << setw(2) <<
"" <<
"Input operations" << endl;
459 stream <<
" " << std::left << setw(20) << op_names[bulk_side][i] <<
"" <<
" " << setw(6) <<
operations_[i].format_shape() <<
"" <<
" "
460 << setw(7) <<
operations_[i].n_dofs() <<
"" <<
" ";
462 for (
auto i_o : input_ops) stream << op_names[bulk_side][i_o] <<
" ";
505 friend class ElOp<spacedim>;
506 template<
unsigned int dim>
508 template<
unsigned int dim>
510 template<
unsigned int dim>
517 template<
unsigned int spacedim = 3>
532 if (shape_vec.size() == 1) shape_vec.push_back(1);
588 result_ = Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>(
shape_[0],
shape_[1]);
594 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> &
result_matrix() {
599 const Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> &
result_matrix()
const {
607 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>
result_;
625 auto &op = operations[vector_sym_grad_op_idx];
626 auto &grad_vector_value = operations[ op.input_ops()[0] ].result_matrix();
627 auto &sym_grad_value = op.result_matrix();
628 sym_grad_value = 0.5 * (grad_vector_value.transpose() + grad_vector_value);
632 auto &op = operations[vector_divergence_op_idx];
633 auto &grad_vector_value = operations[ op.input_ops()[0] ].result_matrix();
634 auto &divergence_value = op.result_matrix();
635 divergence_value(0,0) = grad_vector_value(0,0) + grad_vector_value(1,1) + grad_vector_value(2,2);
642 template<
unsigned int dim>
646 auto &jac_value = op.result_matrix();
647 const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
648 for (
unsigned int i=0; i<3; i++)
649 for (
unsigned int j=0; j<dim; j++)
650 jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
652 template<
unsigned int dim>
656 auto &inv_jac_value = op.result_matrix();
657 const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
658 inv_jac_value = eigen_arena_tools::inverse<3, dim>(jac_value);
660 template<
unsigned int dim>
664 auto &jac_det_value = op.result_matrix();
665 const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
666 jac_det_value(0,0) = eigen_arena_tools::determinant<3, dim>(jac_value).abs();
675 op.allocate_result(point_weights.size(), *arena);
676 auto &weights_value = op.result_matrix();
677 for (
uint i=0; i<point_weights.size(); ++i)
678 weights_value(0,0)(i) = point_weights[i];
682 auto &weights_value = operations[ op.input_ops()[0] ].result_matrix();
683 auto &jac_det_value = operations[ op.input_ops()[1] ].result_matrix();
687 auto &jxw_value = op.result_matrix();
688 jxw_value(0,0) = jxw_ovec.
get_vec();
692 auto &op = operations[scalar_shape_op_idx];
693 uint n_dofs = shape_values.size();
694 uint n_points = shape_values[0].size();
695 uint n_elem = op.result_matrix()(0).data_size() / n_points;
697 auto &shape_matrix = op.result_matrix();
700 for (
uint i=0; i<n_elem; ++i) {
705 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
706 for (
uint i_p=0; i_p<n_points; ++i_p)
707 ref_vec(i_dof * n_points + i_p) = shape_values[i_dof][i_p];
709 shape_matrix(0) = shape_ovec.
get_vec();
712 auto &op = operations[vector_shape_op_idx];
713 uint n_dofs = shape_values.size();
714 uint n_points = shape_values[0].size();
715 uint n_elem = op.result_matrix()(0).data_size() / n_points;
717 auto &shape_matrix = op.result_matrix();
718 Eigen::Matrix<ArenaVec<double>, 3, 1> ref_shape_vec;
719 Eigen::Matrix<ArenaOVec<double>, 3, 1> ref_shape_ovec;
721 for (
uint i=0; i<n_elem; ++i) {
724 for (
uint c=0; c<3; ++c) {
725 ref_shape_vec(c) =
ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
726 ref_shape_ovec(c) =
ArenaOVec(ref_shape_vec(c));
729 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
730 for (
uint i_p=0; i_p<n_points; ++i_p)
731 for (
uint c=0; c<3; ++c)
732 ref_shape_vec(c,0)(i_dof * n_points + i_p) = shape_values[i_dof][i_p](c);
733 Eigen::Matrix<ArenaOVec<double>, 1, 3> shape_omatrix = elem_ovec * ref_shape_ovec.transpose();
734 for (
uint c=0; c<3; ++c)
735 shape_matrix(c) = shape_omatrix(0,c).
get_vec();
751 template<
unsigned int dim>
754 auto &op = operations[scalar_shape_grads_op_idx];
755 auto &inv_jac_vec = operations[ op.input_ops()[0] ].result_matrix();
756 uint n_points = ref_shape_grads.size();
757 uint n_dofs = ref_shape_grads[0].size();
759 Eigen::Matrix<ArenaVec<double>, dim, 1> ref_grads_vec;
760 Eigen::Matrix<ArenaOVec<double>, dim, 1> ref_grads_ovec;
761 for (
uint i=0; i<ref_grads_vec.rows(); ++i) {
762 ref_grads_vec(i) =
ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
763 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
764 for (
uint i_p=0; i_p<n_points; ++i_p)
765 ref_grads_vec(i)(i_dof * n_points + i_p) = ref_shape_grads[i_p][i_dof](i);
766 ref_grads_ovec(i) =
ArenaOVec(ref_grads_vec(i));
769 Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
770 for (
uint i=0; i<dim*3; ++i) {
771 inv_jac_ovec(i) =
ArenaOVec(inv_jac_vec(i));
774 auto &result_vec = op.result_matrix();
775 Eigen::Matrix<ArenaOVec<double>, 3, 1> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
776 for (
uint i=0; i<3; ++i) {
777 result_vec(i) = result_ovec(i).get_vec();
780 template<
unsigned int dim>
783 auto &op = operations[vector_shape_grads_op_idx];
784 auto &inv_jac_vec = operations[ op.input_ops()[0] ].result_matrix();
785 uint n_points = ref_shape_grads.size();
786 uint n_dofs = ref_shape_grads[0].size();
788 Eigen::Matrix<ArenaVec<double>, dim, 3> ref_grads_vec;
789 Eigen::Matrix<ArenaOVec<double>, dim, 3> ref_grads_ovec;
790 for (
uint i=0; i<ref_grads_vec.rows()*ref_grads_vec.cols(); ++i) {
791 ref_grads_vec(i) =
ArenaVec<double>(n_points * n_dofs, op.result_matrix()(0).arena());
792 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
793 for (
uint i_p=0; i_p<n_points; ++i_p)
794 ref_grads_vec(i)(i_dof * n_points + i_p) = ref_shape_grads[i_p][i_dof](i);
795 ref_grads_ovec(i) =
ArenaOVec(ref_grads_vec(i));
798 Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
799 for (
uint i=0; i<dim*3; ++i) {
800 inv_jac_ovec(i) =
ArenaOVec(inv_jac_vec(i));
803 auto &result_vec = op.result_matrix();
804 Eigen::Matrix<ArenaOVec<double>, 3, 3> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
805 for (
uint i=0; i<9; ++i) {
806 result_vec(i) = result_ovec(i).get_vec();
809 template<
unsigned int dim>
838 template<
unsigned int dim>
874 template<
unsigned int dim>
878 auto &jac_value = op.result_matrix();
879 const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
880 for (
unsigned int i=0; i<3; i++)
881 for (
unsigned int j=0; j<dim; j++)
882 jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
884 template<
unsigned int dim>
887 auto &inv_jac_value = op.result_matrix();
888 const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
889 inv_jac_value = eigen_arena_tools::inverse<3, dim>(jac_value);
891 template<
unsigned int dim>
895 auto &jac_value = op.result_matrix();
896 const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
897 for (
unsigned int i=0; i<3; i++)
898 for (
unsigned int j=0; j<dim-1; j++)
899 jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
901 template<
unsigned int dim>
905 auto &jac_det_value = op.result_matrix();
906 const auto &jac_value = operations[ op.input_ops()[0] ].result_matrix();
916 op.allocate_result(point_weights.size(), *arena);
917 auto &weights_value = op.result_matrix();
918 for (
uint i=0; i<point_weights.size(); ++i)
919 weights_value(0,0)(i) = point_weights[i];
923 auto &weights_value = operations[ op.input_ops()[0] ].result_matrix();
924 auto &jac_det_value = operations[ op.input_ops()[1] ].result_matrix();
928 auto &jxw_value = op.result_matrix();
929 jxw_value(0,0) = jxw_ovec.
get_vec();
931 template<
unsigned int dim>
934 auto &normal_value = op.result_matrix();
935 auto &inv_jac_mat_value = operations[ op.input_ops()[0] ].result_matrix();
938 ArenaVec<double> norm_vec( normal_value(0).data_size(), normal_value(0).arena() );
939 Eigen::VectorXd A(3);
940 for (
uint i=0; i<normal_value(0).data_size(); ++i) {
941 A(0) = normal_value(0)(i);
942 A(1) = normal_value(1)(i);
943 A(2) = normal_value(2)(i);
944 norm_vec(i) = A.norm();
947 size_t points_per_side = el_table(4).
data_size() / el_table(3).data_size();
948 size_t n_points = el_table(3).data_size();
949 for (
uint i=0; i<3; ++i) {
950 normal_value(i) = normal_value(i) / norm_vec;
951 ArenaVec<double> expand_vec( normal_value(i).data_size() * points_per_side, normal_value(i).arena() );
953 expand_vec(j) = normal_value(i)(j % n_points);
955 normal_value(i) = expand_vec;
960 uint n_dofs = shape_values[0][0].size();
961 uint n_sides = el_table(3).data_size();
962 uint n_patch_points = el_table(4).data_size();
964 auto &op = operations[scalar_shape_op_idx];
965 auto &scalar_shape_value = op.result_matrix();
966 scalar_shape_value(0) =
ArenaVec<double>(n_dofs*n_patch_points, scalar_shape_value(0).arena());
968 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
969 uint dof_shift = i_dof * n_patch_points;
970 for (
uint i_pt=0; i_pt<n_patch_points; ++i_pt)
971 scalar_shape_value(0)(i_pt + dof_shift) = shape_values[el_table(4)(i_pt)][i_pt / n_sides][i_dof];
976 uint n_dofs = shape_values[0][0].size();
977 uint n_sides = el_table(3).data_size();
978 uint n_patch_points = el_table(4).data_size();
980 auto &op = operations[vector_shape_op_idx];
981 auto &vector_shape_value = op.result_matrix();
982 for (
uint c=0; c<3; c++)
983 vector_shape_value(c) =
ArenaVec<double>(n_dofs*n_patch_points, vector_shape_value(c).arena());
985 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
986 uint dof_shift = i_dof * n_patch_points;
987 for (
uint i_pt=0; i_pt<n_patch_points; ++i_pt)
988 for (
uint c=0; c<3; c++)
989 vector_shape_value(c)(i_pt + dof_shift) = shape_values[el_table(4)(i_pt)][i_pt / n_sides][i_dof](c);
996 template<
unsigned int dim>
999 uint n_points = ref_shape_grads[0].size();
1000 uint n_dofs = ref_shape_grads[0][0].size();
1001 uint n_sides = el_table(3).data_size();
1002 uint n_patch_points = el_table(4).data_size();
1005 auto &op = operations[scalar_shape_grads_op_idx];
1006 auto &grad_scalar_shape_value = op.result_matrix();
1009 auto &inv_jac_value = operations[ op.input_ops()[0] ].result_matrix();
1010 Eigen::Matrix<ArenaVec<double>, dim, 3> inv_jac_expd_value;
1011 for (
uint i=0; i<dim*3; ++i) {
1012 inv_jac_expd_value(i) =
ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(i).arena() );
1013 for (
uint j=0; j<n_dofs*n_patch_points; ++j)
1014 inv_jac_expd_value(i)(j) = inv_jac_value(i)(j%n_sides);
1018 Eigen::Matrix<ArenaVec<double>, dim, 1> ref_shape_grads_expd;
1019 for (
uint i=0; i<dim; ++i) {
1020 ref_shape_grads_expd(i) =
ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(0).arena() );
1022 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
1023 for (
uint i_pt=0; i_pt<n_points; ++i_pt) {
1024 uint i_begin = (i_dof * n_points + i_pt) * n_sides;
1025 for (
uint i_sd=0; i_sd<n_sides; ++i_sd) {
1026 for (
uint i_c=0; i_c<dim; ++i_c) {
1027 ref_shape_grads_expd(i_c)(i_begin + i_sd) = ref_shape_grads[el_table(3)(i_sd)][i_pt][i_dof][i_c];
1034 grad_scalar_shape_value = inv_jac_expd_value.transpose() * ref_shape_grads_expd;
1036 template<
unsigned int dim>
1039 uint n_points = ref_shape_grads[0].size();
1040 uint n_dofs = ref_shape_grads[0][0].size();
1041 uint n_sides = el_table(3).data_size();
1042 uint n_patch_points = el_table(4).data_size();
1045 auto &op = operations[vector_shape_grads_op_idx];
1046 auto &grad_scalar_shape_value = op.result_matrix();
1049 auto &inv_jac_value = operations[ op.input_ops()[0] ].result_matrix();
1050 Eigen::Matrix<ArenaVec<double>, dim, 3> inv_jac_expd_value;
1051 for (
uint i=0; i<dim*3; ++i) {
1052 inv_jac_expd_value(i) =
ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(i).arena() );
1053 for (
uint j=0; j<n_dofs*n_patch_points; ++j)
1054 inv_jac_expd_value(i)(j) = inv_jac_value(i)(j%n_sides);
1058 Eigen::Matrix<ArenaVec<double>, dim, 3> ref_shape_grads_expd;
1059 for (
uint i=0; i<3*dim; ++i) {
1060 ref_shape_grads_expd(i) =
ArenaVec<double>( n_dofs*n_patch_points, inv_jac_value(0).arena() );
1062 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
1063 for (
uint i_pt=0; i_pt<n_points; ++i_pt) {
1064 uint i_begin = (i_dof * n_points + i_pt) * n_sides;
1065 for (
uint i_sd=0; i_sd<n_sides; ++i_sd) {
1066 for (
uint i_c=0; i_c<3*dim; ++i_c) {
1067 ref_shape_grads_expd(i_c)(i_begin + i_sd) = ref_shape_grads[el_table(3)(i_sd)][i_pt][i_dof][i_c];
1074 grad_scalar_shape_value = inv_jac_expd_value.transpose() * ref_shape_grads_expd;
1088 template<
unsigned int dim>
1091 template<
unsigned int dim>
1105 auto &result_vec = op.result_matrix();
1106 for (
uint i=0;i<result_vec(0,0).data_size(); ++i) {
1107 result_vec(0,0)(i) = 1.0;
1116 template<
unsigned int spacedim = 3>
1139 template<
unsigned int dim>
1171 template<
unsigned int spacedim = 3>
1194 template<
unsigned int dim>
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
ArenaVec< T > get_vec() const
Convert ArenaOVec to ArenaVec and its.
size_t data_size() const
Getter for data_size_.
Class represents element or FE operations.
ElOp(uint dim, std::initializer_list< uint > shape, ReinitFunction reinit_f, OpSizeType size_type, std::vector< uint > input_ops={}, uint n_dofs=1)
Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > & result_matrix()
Return map referenced result as Eigen::Matrix.
void allocate_result(size_t data_size, PatchArena &arena)
std::string format_shape() const
std::vector< uint > input_ops_
Indices of operations in PatchPointValues::operations_ vector on which ElOp is depended.
std::vector< uint > set_shape_vec(std::initializer_list< uint > shape) const
Aligns shape_vec to 2 items (equal to matrix number of dimensions)
Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > result_
Result matrix of operation.
std::vector< uint > shape_
Shape of stored data (size of vector or number of rows and cols of matrix)
const std::vector< uint > & shape() const
Getter for shape_.
const std::vector< uint > & input_ops() const
Getter for input_ops_.
OpSizeType size_type_
Type of operation by size of vector (element, point or fixed size)
const Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > & result_matrix() const
Same as previous but return const reference.
uint n_dofs() const
Getter for n_dofs_.
ReinitFunction reinit_func
Pointer to patch reinit function of element data table specialized by operation.
void reinit_function(std::vector< ElOp< spacedim >> &operations, IntTableArena &int_table)
Call reinit function on element table if function is defined.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter for dimension.
OpSizeType size_type() const
Getter for size_type_.
Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum.
PatchPointValues(uint dim, uint quad_order, AssemblyArena &asm_arena)
Constructor.
void init()
Initialize operations vector.
Bulk Side specialization, order of item in operations_ vector corresponds to the SideOps enum.
PatchPointValues(uint dim, uint quad_order, AssemblyArena &asm_arena)
Constructor.
void init()
Initialize operations vector.
void print_data_tables(ostream &stream, bool points, bool ints) const
PatchPointValues(uint dim, AssemblyArena &asm_arena)
AssemblyArena & asm_arena_
Reference to global assembly arena of PatchFeValues.
std::vector< OpSizeType > row_sizes_
hold sizes of rows by type of operation
std::vector< uint > points_map_
Map of point patch indices to ElOp::result_ and int_table_ tables.
std::vector< unsigned int > fe_ops_indices_
Indices of FE operations in operations_ vector.
uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint i_point_on_elem)
uint register_element(arma::mat coords, uint element_patch_idx)
void resize_tables(uint n_elems, uint n_points)
Resize data tables. Method is called before reinit of patch.
unsigned int get_fe_op(FEOps fe_op) const
Return index of FE operation in operations_ vector.
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)
void set_fe_op(FEOps fe_op, unsigned int op_vec_idx)
Set index of FE operation in operations_ vector.
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Getter for quadrature.
uint n_elems() const
Getter for n_elems_.
uint i_elem_
Index of registered element in table, helper value used during patch creating.
void reset()
Reset number of columns (points and elements)
uint n_points() const
Getter for n_points_.
ElOp< spacedim > & make_fixed_op(std::initializer_list< uint > shape, ReinitFunction reinit_f)
std::vector< uint > elements_map_
Map of element patch indices to ElOp::result_ and int_table_ tables.
uint dim() const
Getter for dim_.
Tensor tensor_value(uint op_idx, uint point_idx, uint i_dof=0) const
Scalar scalar_value(uint op_idx, uint point_idx, uint i_dof=0) const
void init_finalize(PatchArena *patch_arena)
uint register_side(arma::mat elm_coords, arma::mat side_coords, uint side_idx)
uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx, uint i_point_on_side)
uint n_points_
Number of points in patch.
void print_operations(ostream &stream, uint bulk_side) const
Quadrature * quad_
Quadrature of given dimension and order passed in constructor.
ElOp< spacedim > & make_new_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, OpSizeType size_type=pointOp)
std::vector< OpSizeType > int_sizes_
Set size and type of rows of int_table_, value is set implicitly in constructor of descendants.
Vector vector_value(uint op_idx, uint point_idx, uint i_dof=0) const
PatchArena * patch_arena_
Pointer to global patch arena of PatchFeValues.
uint n_elems_
Number of elements in patch.
Symmetric Gauss-Legendre quadrature formulae on simplices.
Base class for quadrature rules on simplices in arbitrary dimensions.
const std::vector< double > & get_weights() const
Return a reference to the whole array of weights.
static Eigen::Matrix< ArenaVec< double >, dim, 1 > normal_vector_array(ArenaVec< uint > loc_side_idx_vec)
Class ElementValues calculates data related to transformation of reference cell to actual cell (Jacob...
ArmaMat< double, N, M > mat
@ opCoords
operations evaluated on quadrature points
@ opInvJac
inverse Jacobian
@ opWeights
weight of quadrature point
@ opJac
Jacobian of element.
@ opElCoords
operations evaluated on elements
@ opJxW
JxW value of quadrature point.
@ opNormalVec
normal vector of quadrature point
@ opSideCoords
operations evaluated on sides
@ opJxW
JxW value of quadrature point.
@ opCoords
operations evaluated on quadrature points
@ opElCoords
operations evaluated on elements
@ opElJac
Jacobian of element.
@ opWeights
weight of quadrature point
@ opSideJac
Jacobian of element.
@ OpNItems
Holds number of valid FE operations and value of invalid FE operation.
@ opVectorShape
Vector shape operation.
@ opVectorDivergence
Vector divergence.
@ opVectorSymGrad
Vector symmetric gradient.
@ opGradScalarShape
Scalar shape gradient.
@ opGradVectorShape
Vector shape gradient.
@ opScalarShape
Scalar shape operation.
std::function< void(std::vector< ElOp< 3 > > &, IntTableArena &)> ReinitFunction
Type for conciseness.
OpSizeType
Distinguishes operations by type and size of output rows.
@ pointOp
operation is evaluated on quadrature points
@ elemOp
operation is evaluated on elements or sides
@ fixedSizeOp
operation has fixed size and it is filled during initialization
Definitions of particular quadrature rules on simplices.
Defines reinit operations on bulk points.
static void ptop_vector_contravariant_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::mat > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)
static void elop_inv_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_shape(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::vec3 > > shape_values, uint vector_shape_op_idx)
static void ptop_vector_contravariant_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::vec3 > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::mat > > ref_shape_grads, uint scalar_shape_grads_op_idx)
static void ptop_vector_piola_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::vec3 > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void elop_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_piola_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED std::vector< std::vector< arma::mat > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)
static void elop_jac_det(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, PatchArena *arena, const std::vector< double > &point_weights)
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_vector_shape_grads(std::vector< ElOp< 3 >> &operations, std::vector< std::vector< arma::mat > > ref_shape_grads, uint vector_shape_grads_op_idx)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
Defines common functionality of reinit operations.
static void ptop_vector_divergence(std::vector< ElOp< 3 >> &operations, uint vector_divergence_op_idx)
Common reinit function of vector divergence on bulk and side points.
static void ptop_vector_sym_grad(std::vector< ElOp< 3 >> &operations, uint vector_sym_grad_op_idx)
Common reinit function of vector symmetric gradient on bulk and side points.
static void op_base(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
Defines reinit operations on side points.
static void ptop_vector_shape(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, uint vector_shape_op_idx)
static void ptop_vector_shape_grads(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, uint vector_shape_grads_op_idx)
static void elop_el_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void elop_sd_jac_det(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void elop_sd_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_normal_vec(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, PatchArena *arena, const std::vector< double > &point_weights)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, IntTableArena &el_table, std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, uint scalar_shape_grads_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)
static void ptop_vector_contravariant_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)
static void ptop_vector_contravariant_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void elop_el_inv_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_piola_shape(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::vec3 > > > shape_values, FMT_UNUSED uint vector_shape_op_idx)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table)
static void ptop_vector_piola_shape_grads(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED IntTableArena &el_table, FMT_UNUSED std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, FMT_UNUSED uint vector_shape_grads_op_idx)