21 #ifndef PATCH_POINT_VALUES_HH_
22 #define PATCH_POINT_VALUES_HH_
24 #include <Eigen/Dense>
33 template<
unsigned int spacedim>
class ElOp;
119 template<
unsigned int spacedim = 3>
197 for (
uint i_col=0; i_col<coords.n_cols; ++i_col)
198 for (
uint i_row=0; i_row<coords.n_rows; ++i_row) {
215 for (
uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
216 for (
uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
222 for (
uint i_col=0; i_col<side_coords.n_cols; ++i_col)
223 for (
uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
354 for (
uint i=0; i<3; ++i)
367 for (
uint i=0; i<3; ++i)
368 for (
uint j=0; j<3; ++j)
393 for (
uint i_col=0; i_col<3; ++i_col)
394 stream <<
int_vals_(i_col)(i_row) <<
" ";
410 {
"el_coords",
"jacobian",
"inv_jac",
"jac_det",
"pt_coords",
"weights",
"JxW",
"",
"",
"",
"",
"" },
411 {
"el_coords",
"el_jac",
"el_inv_jac",
"side_coords",
"side_jac",
"side_jac_det",
"exp_el_coords",
"exp_el_jac",
"exp_el_inv_jac",
412 "exp_side_coords",
"exp_side_jac",
"exp_side_jac_det",
"pt_coords",
"weights",
"JxW",
"normal_vec",
"",
"",
"",
"",
"" }
414 stream <<
"Rows sizes:" << endl;
418 stream << std::setfill(
'=') << setw(100) <<
"" << endl;
419 stream << std::setfill(
' ') <<
" Operation" << setw(12) <<
"" <<
"Shape" << setw(2) <<
""
420 <<
"Result row" << setw(2) <<
"" <<
"n DOFs" << setw(2) <<
"" <<
"Input operations" << endl;
422 stream <<
" " << std::left << setw(20) << op_names[bulk_side][i] <<
"" <<
" " << setw(6) <<
operations_[i].format_shape() <<
"" <<
" "
423 << setw(11) <<
operations_[i].result_row() <<
"" <<
" " << setw(7) <<
operations_[i].n_dofs() <<
"" <<
" ";
425 for (
auto i_o : input_ops) stream << op_names[bulk_side][i_o] <<
" ";
463 friend class ElOp<spacedim>;
464 template<
unsigned int dim>
466 template<
unsigned int dim>
468 template<
unsigned int dim>
475 template<
unsigned int spacedim = 3>
540 template<
unsigned int dim1,
unsigned int dim2>
541 Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>
value(
TableDbl &op_results,
uint i_dof = 0)
const {
542 return Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>(op_results.data() +
result_row_ + i_dof *
n_comp(), dim1, dim2);
547 return Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>>(op_results(
result_row_).data(), dim1, dim2);
552 return Eigen::Map<Eigen::Vector<double, Eigen::Dynamic>>(op_results(
result_row_).data(), op_results(
result_row_).rows());
578 uint size = op_results(row_begin).rows();
579 for (
int i_pt=size-1; i_pt>=0; --i_pt) {
580 uint el_table_idx = el_table(1)(i_pt);
581 for (
uint i_q=row_begin; i_q<row_end; ++i_q)
582 op_results(i_q)(i_pt) = op_results(i_q)(el_table_idx);
590 template<
unsigned int dim>
594 auto jac_value = op.value<3, dim>(op_results);
595 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
596 jac_value = eigen_tools::jacobian<3,dim>(coords_value);
598 template<
unsigned int dim>
602 auto inv_jac_value = op.value<dim, 3>(op_results);
603 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
604 inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
606 template<
unsigned int dim>
610 auto &jac_det_value = op_results(op.result_row());
611 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
612 jac_det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim>>(jac_value).
array().abs();
621 ArrayDbl &result_row = op_results( op.result_row() );
622 result_row.resize(point_weights.rows());
623 result_row << point_weights;
628 Eigen::Map<Eigen::Vector<double, Eigen::Dynamic>> jac_det_value = operations[op.input_ops()[1]].vector_value(op_results);
629 Eigen::Map<Eigen::Vector<double, Eigen::Dynamic>> weights_value = operations[op.input_ops()[0]].vector_value(op_results);
630 Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>> result_value = op.matrix_value(op_results, weights_value.rows(), jac_det_value.rows());
631 result_value = weights_value * jac_det_value.transpose();
640 auto &op = operations[scalar_shape_op_idx];
641 uint n_points = shape_values.size();
643 for (
uint i_row=0; i_row<shape_values[0].size(); ++i_row) {
644 ArrayDbl &result_row = op_results( op.result_row()+i_row );
645 for (
uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
646 result_row(i_pt) = shape_values[i_pt % n_points][i_row];
649 template<
unsigned int dim>
652 auto &op = operations[scalar_shape_grads_op_idx];
653 uint n_points = ref_shape_grads.size();
654 uint n_dofs = ref_shape_grads[0].size();
656 Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
657 ref_shape_grads_expd.resize(dim*n_dofs);
658 for (
uint i=0; i<ref_shape_grads_expd.rows(); ++i)
659 ref_shape_grads_expd(i).resize(op_results(0).
rows());
661 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
662 for (
uint i_c=0; i_c<dim; ++i_c) {
663 ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
664 for (
uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
665 shape_grad_row(i_pt) = ref_shape_grads[i_pt % n_points][i_dof][i_c];
668 auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
669 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
670 auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
671 Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
672 shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
682 template<
unsigned int dim>
686 auto jac_value = op.value<3, dim>(op_results);
687 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
688 jac_value = eigen_tools::jacobian<3,dim>(coords_value);
690 template<
unsigned int dim>
693 auto inv_jac_value = op.value<dim, 3>(op_results);
694 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
695 inv_jac_value = eigen_tools::inverse<3,dim>(jac_value);
697 template<
unsigned int dim>
701 auto jac_value = op.value<3, dim-1>(op_results);
702 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
705 template<
unsigned int dim>
709 ArrayDbl &det_value = op_results( op.result_row() );
710 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim-1>(op_results);
732 template<
unsigned int dim>
748 ArrayDbl &result_row = op_results( op.result_row() );
749 auto n_points = point_weights.size();
750 for (
uint i=0; i<result_row.rows(); ++i)
751 result_row(i) = point_weights[i%n_points];
755 ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
756 ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
757 ArrayDbl &result_row = op_results( op.result_row() );
758 result_row = jac_det_row * weights_row;
760 template<
unsigned int dim>
763 auto normal_value = op.value<3, 1>(op_results);
764 auto inv_jac_mat_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
768 norm_vec.resize(normal_value(0).
rows());
769 Eigen::VectorXd A(3);
771 for (
uint i=0; i<normal_value(0).rows(); ++i) {
772 A(0) = normal_value(0)(i);
773 A(1) = normal_value(1)(i);
774 A(2) = normal_value(2)(i);
775 norm_vec(i) = A.norm();
777 for (
uint i=0; i<3; ++i) {
778 normal_value(i) /= norm_vec;
783 auto &op = operations[scalar_shape_op_idx];
784 uint n_points = shape_values[0].size();
786 for (
uint i_row=0; i_row<shape_values[0][0].size(); ++i_row) {
787 ArrayDbl &result_row = op_results( op.result_row()+i_row );
788 for (
uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
789 result_row(i_pt) = shape_values[el_table(3)(i_pt)][i_pt % n_points][i_row];
792 template<
unsigned int dim>
795 auto &op = operations[scalar_shape_grads_op_idx];
796 uint n_points = ref_shape_grads[0].size();
797 uint n_dofs = ref_shape_grads[0][0].size();
799 Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
800 ref_shape_grads_expd.resize(dim*n_dofs);
801 for (
uint i=0; i<ref_shape_grads_expd.rows(); ++i)
802 ref_shape_grads_expd(i).resize(op_results(0).
rows());
804 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
805 for (
uint i_c=0; i_c<dim; ++i_c) {
806 ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
807 for (
uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
808 shape_grad_row(i_pt) = ref_shape_grads[el_table(3)(i_pt)][i_pt % n_points][i_dof][i_c];
811 auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
812 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
813 auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
814 Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
815 shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
829 ArrayDbl &result_vec = op_results( op.result_row() );
830 for (
uint i=0;i<result_vec.size(); ++i) {
844 template<
unsigned int spacedim = 3>
866 template<
unsigned int dim>
882 ArrayDbl point_weights(point_weights_vec.size());
883 for (
uint i=0; i<point_weights_vec.size(); ++i)
884 point_weights(i) = point_weights_vec[i];
901 template<
unsigned int spacedim = 3>
923 template<
unsigned int dim>
947 this->
make_expansion( sd_jac, {spacedim, this->
dim_-1}, &side_reinit::expd_sd_jac<dim> );
Class represents element or FE operations.
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 > shape_
Shape of stored data (size of vector or number of rows and cols of matrix)
uint result_row() const
Getter for result_row_.
const std::vector< uint > & shape() const
Getter for shape_.
const std::vector< uint > & input_ops() const
Getter for input_ops_.
uint result_row_
First row to scalar, vector or matrix result.
uint n_dofs() const
Getter for n_dofs_.
ReinitFunction reinit_func
Pointer to patch reinit function of element data table specialized by operation.
ElOp(uint dim, std::initializer_list< uint > shape, uint result_row, ReinitFunction reinit_f, std::vector< uint > input_ops={}, uint n_dofs=1)
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter for dimension.
Eigen::Map< Eigen::Matrix< ArrayDbl, dim1, dim2 > > value(TableDbl &op_results, uint i_dof=0) const
Return map referenced Eigen::Matrix of given dimension.
Eigen::Map< Eigen::Vector< double, Eigen::Dynamic > > vector_value(TableDbl &op_results) const
Return map referenced Eigen::Matrix of given dimensions.
Eigen::Map< Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > > matrix_value(TableDbl &op_results, uint dim1, uint dim2) const
Return map referenced Eigen::Matrix of given dimensions.
void reinit_function(std::vector< ElOp< spacedim >> &operations, TableDbl &data_table, TableInt &int_table)
Call reinit function on element table if function is defined.
Bulk data specialization, order of item in operations_ vector corresponds to the BulkOps enum.
PatchPointValues(uint dim, uint quad_order)
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)
Constructor.
void init()
Initialize operations vector.
void print_data_tables(ostream &stream, bool points, bool ints) const
ElOp< spacedim > & make_expansion(ElOp< spacedim > &el_op, std::initializer_list< uint > shape, ReinitFunction reinit_f)
std::vector< OpSizeType > row_sizes_
hold sizes of rows by type of operation
void initialize(uint int_cols)
uint n_rows() const
Getter for n_rows_.
std::vector< uint > points_map_
Map of point patch indices to point_vals_ and int_vals_ tables.
uint register_side_point(uint elem_table_row, uint value_patch_idx, uint elem_idx, uint side_idx)
uint n_rows_
Number of columns of point_vals table.
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.
Tensor tensor_val(uint result_row, uint point_idx) const
uint register_side(arma::mat elm_coords, arma::mat side_coords)
PatchPointValues(uint dim)
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 n_elems() const
Getter for n_elems_.
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 el_vals_ table.
uint dim() const
Getter for dim_.
Vector vector_val(uint result_row, uint point_idx) const
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.
Scalar scalar_val(uint result_row, uint point_idx) const
ElOp< spacedim > & make_new_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec, OpSizeType size_type=pointOp)
uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx)
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::Vector< ArrayDbl, dim > normal_vector_array(Eigen::Array< uint, Eigen::Dynamic, 1 > loc_side_idx_array)
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.
@ opExpansionElCoords
operation executed expansion to quadrature point (value of element / side > values on quadrature poin...
@ opExpansionElJac
expands Jacobian on element
@ opNormalVec
normal vector of quadrature point
@ opSideCoords
operations evaluated on sides
@ opExpansionSideCoords
expands coordinates on side
@ opJxW
JxW value of quadrature point.
@ opCoords
operations evaluated on quadrature points
@ opExpansionElInvJac
expands inverse Jacobian on element
@ opExpansionSideJac
expands Jacobian on side
@ opElCoords
operations evaluated on elements
@ opElJac
Jacobian of element.
@ opWeights
weight of quadrature point
@ opSideJac
Jacobian of element.
std::function< void(std::vector< ElOp< 3 > > &, TableDbl &, TableInt &)> 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_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< std::vector< arma::mat > > ref_shape_grads, uint scalar_shape_grads_op_idx)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
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 elop_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, ArrayDbl point_weights)
static void elop_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
Defines common functionality of reinit operations.
static void expand_data(ElOp< 3 > &op, TableDbl &op_results, TableInt &el_table)
static void op_base(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
Defines reinit operations on side points.
static void ptop_normal_vec(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void elop_el_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_weights(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, std::vector< double > point_weights)
static void elop_sd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_sd_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
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)
static void expd_el_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_sd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_JxW(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_el_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_el_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void elop_el_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void expd_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void ptop_coords(FMT_UNUSED std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void elop_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, FMT_UNUSED TableInt &el_table)
static void ptop_scalar_shape_grads(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table, std::vector< std::vector< std::vector< arma::mat > > > ref_shape_grads, uint scalar_shape_grads_op_idx)