21 #ifndef PATCH_POINT_VALUES_HH_
22 #define PATCH_POINT_VALUES_HH_
24 #include <Eigen/Dense>
34 template<
unsigned int spacedim>
class ElOp;
120 template<
unsigned int spacedim = 3>
181 elOp.allocate_result(sizes[elOp.size_type()],
patch_arena_);
206 for (
uint i_col=0; i_col<coords.n_cols; ++i_col)
207 for (
uint i_row=0; i_row<coords.n_rows; ++i_row) {
209 coords_mat(i_row, i_col)(i_elem) = coords(i_row, i_col);
225 for (
uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
226 for (
uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
232 for (
uint i_col=0; i_col<side_coords.n_cols; ++i_col)
233 for (
uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
364 for (
uint i=0; i<3; ++i)
377 for (
uint i=0; i<3; ++i)
378 for (
uint j=0; j<3; ++j)
394 const auto &
mat = op.result_matrix();
395 for (
uint i_mat=0; i_mat<
mat.rows()*
mat.cols(); ++i_mat) {
396 if (
mat(i_mat).data_size()==0) stream <<
"<empty>";
398 const double *vals =
mat(i_mat).data_ptr();
399 for (
size_t i_val=0; i_val<
mat(i_mat).data_size(); ++i_val)
400 stream << vals[i_val] <<
" ";
410 for (
uint i_col=0; i_col<3; ++i_col)
411 stream <<
int_vals_(i_col)(i_row) <<
" ";
427 {
"el_coords",
"jacobian",
"inv_jac",
"jac_det",
"pt_coords",
"weights",
"JxW",
"",
"",
"",
"",
"" },
428 {
"el_coords",
"el_jac",
"el_inv_jac",
"side_coords",
"side_jac",
"side_jac_det",
"exp_el_coords",
"exp_el_jac",
"exp_el_inv_jac",
429 "exp_side_coords",
"exp_side_jac",
"exp_side_jac_det",
"pt_coords",
"weights",
"JxW",
"normal_vec",
"",
"",
"",
"",
"" }
431 stream << std::setfill(
' ') <<
" Operation" << setw(12) <<
"" <<
"Shape" << setw(2) <<
""
432 <<
"Result row" << setw(2) <<
"" <<
"n DOFs" << setw(2) <<
"" <<
"Input operations" << endl;
434 stream <<
" " << std::left << setw(20) << op_names[bulk_side][i] <<
"" <<
" " << setw(6) <<
operations_[i].format_shape() <<
"" <<
" "
435 << setw(11) <<
operations_[i].result_row() <<
"" <<
" " << setw(7) <<
operations_[i].n_dofs() <<
"" <<
" ";
437 for (
auto i_o : input_ops) stream << op_names[bulk_side][i_o] <<
" ";
476 friend class ElOp<spacedim>;
477 template<
unsigned int dim>
479 template<
unsigned int dim>
481 template<
unsigned int dim>
488 template<
unsigned int spacedim = 3>
503 if (shape_vec.size() == 1) shape_vec.push_back(1);
564 result_ = Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>(
shape_[0],
shape_[1]);
570 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> &
result_matrix() {
578 const Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> &
result_matrix()
const {
586 template<
unsigned int dim1,
unsigned int dim2>
587 Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>
value(
TableDbl &op_results,
uint i_dof = 0)
const {
588 return Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>(op_results.data() +
result_row_ + i_dof *
n_comp(), dim1, dim2);
593 return Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>>(op_results(
result_row_).data(), dim1, dim2);
598 return Eigen::Map<Eigen::Vector<double, Eigen::Dynamic>>(op_results(
result_row_).data(), op_results(
result_row_).rows());
606 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic>
result_;
627 uint size = op_results(row_begin).rows();
628 for (
int i_pt=size-1; i_pt>=0; --i_pt) {
629 uint el_table_idx = el_table(1)(i_pt);
630 for (
uint i_q=row_begin; i_q<row_end; ++i_q)
631 op_results(i_q)(i_pt) = op_results(i_q)(el_table_idx);
639 template<
unsigned int dim>
643 auto &jac_value = op.result_matrix();
644 const auto &coords_value = operations[ op.input_ops()[0] ].result_matrix();
645 for (
unsigned int i=0; i<3; i++)
646 for (
unsigned int j=0; j<dim; j++)
647 jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
649 template<
unsigned int dim>
653 auto inv_jac_value = op.value<dim, 3>(op_results);
654 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
655 inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
657 template<
unsigned int dim>
661 auto &jac_det_value = op_results(op.result_row());
662 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
663 jac_det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim>>(jac_value).
array().abs();
672 ArrayDbl &result_row = op_results( op.result_row() );
673 result_row.resize(point_weights.rows());
674 result_row << point_weights;
679 Eigen::Map<Eigen::Vector<double, Eigen::Dynamic>> jac_det_value = operations[op.input_ops()[1]].vector_value(op_results);
680 Eigen::Map<Eigen::Vector<double, Eigen::Dynamic>> weights_value = operations[op.input_ops()[0]].vector_value(op_results);
681 Eigen::Map<Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>> result_value = op.matrix_value(op_results, weights_value.rows(), jac_det_value.rows());
682 result_value = weights_value * jac_det_value.transpose();
691 auto &op = operations[scalar_shape_op_idx];
692 uint n_points = shape_values.size();
694 for (
uint i_row=0; i_row<shape_values[0].size(); ++i_row) {
695 ArrayDbl &result_row = op_results( op.result_row()+i_row );
696 for (
uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
697 result_row(i_pt) = shape_values[i_pt % n_points][i_row];
700 template<
unsigned int dim>
703 auto &op = operations[scalar_shape_grads_op_idx];
704 uint n_points = ref_shape_grads.size();
705 uint n_dofs = ref_shape_grads[0].size();
707 Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
708 ref_shape_grads_expd.resize(dim*n_dofs);
709 for (
uint i=0; i<ref_shape_grads_expd.rows(); ++i)
710 ref_shape_grads_expd(i).resize(op_results(0).
rows());
712 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
713 for (
uint i_c=0; i_c<dim; ++i_c) {
714 ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
715 for (
uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
716 shape_grad_row(i_pt) = ref_shape_grads[i_pt % n_points][i_dof][i_c];
719 auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
720 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
721 auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
722 Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
723 shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
733 template<
unsigned int dim>
737 auto jac_value = op.value<3, dim>(op_results);
738 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
739 jac_value = eigen_tools::jacobian<3,dim>(coords_value);
741 template<
unsigned int dim>
744 auto inv_jac_value = op.value<dim, 3>(op_results);
745 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
746 inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
748 template<
unsigned int dim>
752 auto jac_value = op.value<3, dim-1>(op_results);
753 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
756 template<
unsigned int dim>
760 ArrayDbl &det_value = op_results( op.result_row() );
761 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim-1>(op_results);
783 template<
unsigned int dim>
799 ArrayDbl &result_row = op_results( op.result_row() );
800 auto n_points = point_weights.size();
801 for (
uint i=0; i<result_row.rows(); ++i)
802 result_row(i) = point_weights[i%n_points];
806 ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
807 ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
808 ArrayDbl &result_row = op_results( op.result_row() );
809 result_row = jac_det_row * weights_row;
811 template<
unsigned int dim>
814 auto normal_value = op.value<3, 1>(op_results);
815 auto inv_jac_mat_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
819 norm_vec.resize(normal_value(0).
rows());
820 Eigen::VectorXd A(3);
822 for (
uint i=0; i<normal_value(0).rows(); ++i) {
823 A(0) = normal_value(0)(i);
824 A(1) = normal_value(1)(i);
825 A(2) = normal_value(2)(i);
826 norm_vec(i) = A.norm();
828 for (
uint i=0; i<3; ++i) {
829 normal_value(i) /= norm_vec;
834 auto &op = operations[scalar_shape_op_idx];
835 uint n_points = shape_values[0].size();
837 for (
uint i_row=0; i_row<shape_values[0][0].size(); ++i_row) {
838 ArrayDbl &result_row = op_results( op.result_row()+i_row );
839 for (
uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
840 result_row(i_pt) = shape_values[el_table(3)(i_pt)][i_pt % n_points][i_row];
843 template<
unsigned int dim>
846 auto &op = operations[scalar_shape_grads_op_idx];
847 uint n_points = ref_shape_grads[0].size();
848 uint n_dofs = ref_shape_grads[0][0].size();
850 Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
851 ref_shape_grads_expd.resize(dim*n_dofs);
852 for (
uint i=0; i<ref_shape_grads_expd.rows(); ++i)
853 ref_shape_grads_expd(i).resize(op_results(0).
rows());
855 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
856 for (
uint i_c=0; i_c<dim; ++i_c) {
857 ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
858 for (
uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
859 shape_grad_row(i_pt) = ref_shape_grads[el_table(3)(i_pt)][i_pt % n_points][i_dof][i_c];
862 auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
863 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
864 auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
865 Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
866 shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
880 ArrayDbl &result_vec = op_results( op.result_row() );
881 for (
uint i=0;i<result_vec.size(); ++i) {
895 template<
unsigned int spacedim = 3>
917 template<
unsigned int dim>
956 template<
unsigned int spacedim = 3>
978 template<
unsigned int dim>
1002 this->
make_expansion( sd_jac, {spacedim, this->
dim_-1}, &side_reinit::expd_sd_jac<dim> );
#define ASSERT_EQ(a, b)
Definition of comparative assert macro (EQual) only for debug mode.
Class represents element or FE operations.
Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > & result_matrix()
Return map referenced result as Eigen::Matrix.
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)
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_.
OpSizeType size_type_
Type of operation by size of vector (element, point or fixed size)
uint result_row_
First row to scalar, vector or matrix result TODO replace.
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.
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.
void allocate_result(size_t data_size, AssemblyArena &arena)
ElOp(uint dim, std::initializer_list< uint > shape, uint result_row, ReinitFunction reinit_f, OpSizeType size_type, std::vector< uint > input_ops={}, uint n_dofs=1)
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.
OpSizeType size_type() const
Getter for size_type_.
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, AssemblyArena &patch_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 &patch_arena)
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)
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.
PatchPointValues(uint dim, AssemblyArena &patch_arena)
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.
AssemblyArena & patch_arena_
Reference to global Arena of PatchFeValues.
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 elop_jac(std::vector< ElOp< 3 >> &operations, FMT_UNUSED TableDbl &op_results, FMT_UNUSED TableInt &el_table)
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 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)