21 #ifndef PATCH_POINT_VALUES_HH_
22 #define PATCH_POINT_VALUES_HH_
24 #include <Eigen/Dense>
33 template<
unsigned int spacedim>
class ElOp;
102 template<
unsigned int spacedim = 3>
175 for (
uint i_col=0; i_col<coords.n_cols; ++i_col)
176 for (
uint i_row=0; i_row<coords.n_rows; ++i_row) {
188 for (
uint i_col=0; i_col<elm_coords.n_cols; ++i_col)
189 for (
uint i_row=0; i_row<elm_coords.n_rows; ++i_row) {
195 for (
uint i_col=0; i_col<side_coords.n_cols; ++i_col)
196 for (
uint i_row=0; i_row<side_coords.n_rows; ++i_row) {
283 for (
uint i=0; i<3; ++i)
290 for (
uint i=0; i<3; ++i)
291 for (
uint j=0; j<3; ++j)
297 void print(
bool points,
bool ints)
const {
298 std::cout <<
"** Dimension: " <<
dim_ << std::endl;
304 std::cout << std::endl;
306 std::cout << std::endl;
309 std::cout <<
"Int vals: " <<
int_vals_.rows() <<
" - " <<
int_vals_.cols() << std::endl;
311 for (
uint i_col=0; i_col<3; ++i_col)
312 std::cout <<
int_vals_(i_col)(i_row) <<
" ";
313 std::cout << std::endl;
315 std::cout << std::endl;
317 std::cout <<
"*****************" << std::endl;
352 friend class ElOp<spacedim>;
353 template<
unsigned int dim>
355 template<
unsigned int dim>
357 template<
unsigned int dim>
364 template<
unsigned int spacedim = 3>
404 template<
unsigned int dim1,
unsigned int dim2>
405 Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>
value(
TableDbl &op_results,
uint i_dof = 0)
const {
406 return Eigen::Map<Eigen::Matrix<ArrayDbl, dim1, dim2>>(op_results.data() +
result_row_ + i_dof *
n_comp(), dim1, dim2);
432 uint size = op_results(row_begin).rows();
433 for (
int i_pt=size-1; i_pt>=0; --i_pt) {
434 uint el_table_idx = el_table(1)(i_pt);
435 for (
uint i_q=row_begin; i_q<row_end; ++i_q)
436 op_results(i_q)(i_pt) = op_results(i_q)(el_table_idx);
444 template<
unsigned int dim>
448 auto jac_value = op.value<3, dim>(op_results);
449 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
450 jac_value = eigen_tools::jacobian<3,dim>(coords_value);
452 template<
unsigned int dim>
456 auto inv_jac_value = op.value<dim, 3>(op_results);
457 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
458 inv_jac_value = eigen_tools::inverse<3, dim>(jac_value);
460 template<
unsigned int dim>
464 auto &jac_det_value = op_results(op.result_row());
465 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
466 jac_det_value = eigen_tools::determinant<Eigen::Matrix<ArrayDbl, 3, dim>>(jac_value);
493 ArrayDbl &result_row = op_results( op.result_row() );
494 auto n_points = point_weights.size();
495 for (
uint i=0; i<result_row.rows(); ++i)
496 result_row(i) = point_weights[i%n_points];
500 ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
501 ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
502 ArrayDbl &result_row = op_results( op.result_row() );
503 result_row = jac_det_row * weights_row;
507 auto &op = operations[scalar_shape_op_idx];
508 uint n_points = shape_values.size();
510 for (
uint i_row=0; i_row<shape_values[0].size(); ++i_row) {
511 ArrayDbl &result_row = op_results( op.result_row()+i_row );
512 for (
uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
513 result_row(i_pt) = shape_values[i_pt % n_points][i_row];
516 template<
unsigned int dim>
519 auto &op = operations[scalar_shape_grads_op_idx];
520 uint n_points = ref_shape_grads.size();
521 uint n_dofs = ref_shape_grads[0].size();
523 Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
524 ref_shape_grads_expd.resize(dim*n_dofs);
525 for (
uint i=0; i<ref_shape_grads_expd.rows(); ++i)
526 ref_shape_grads_expd(i).resize(op_results(0).
rows());
528 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
529 for (
uint i_c=0; i_c<dim; ++i_c) {
530 ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
531 for (
uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
532 shape_grad_row(i_pt) = ref_shape_grads[i_pt % n_points][i_dof][i_c];
535 auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
536 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
537 auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
538 Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
539 shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
549 template<
unsigned int dim>
553 auto jac_value = op.value<3, dim>(op_results);
554 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim+1>(op_results);
555 jac_value = eigen_tools::jacobian<3,dim>(coords_value);
557 template<
unsigned int dim>
560 auto inv_jac_value = op.value<dim, 3>(op_results);
561 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
562 inv_jac_value = eigen_tools::inverse<3,dim>(jac_value);
564 template<
unsigned int dim>
568 auto jac_value = op.value<3, dim-1>(op_results);
569 auto coords_value = operations[ op.input_ops()[0] ].value<3, dim>(op_results);
572 template<
unsigned int dim>
576 ArrayDbl &det_value = op_results( op.result_row() );
577 auto jac_value = operations[ op.input_ops()[0] ].value<3, dim-1>(op_results);
614 ArrayDbl &result_row = op_results( op.result_row() );
615 auto n_points = point_weights.size();
616 for (
uint i=0; i<result_row.rows(); ++i)
617 result_row(i) = point_weights[i%n_points];
621 ArrayDbl &weights_row = op_results( operations[op.input_ops()[0]].result_row() );
622 ArrayDbl &jac_det_row = op_results( operations[op.input_ops()[1]].result_row() );
623 ArrayDbl &result_row = op_results( op.result_row() );
624 result_row = jac_det_row * weights_row;
626 template<
unsigned int dim>
629 auto normal_value = op.value<3, 1>(op_results);
630 auto inv_jac_mat_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
635 auto &op = operations[scalar_shape_op_idx];
636 uint n_points = shape_values[0].size();
638 for (
uint i_row=0; i_row<shape_values[0][0].size(); ++i_row) {
639 ArrayDbl &result_row = op_results( op.result_row()+i_row );
640 for (
uint i_pt=0; i_pt<result_row.rows(); ++i_pt)
641 result_row(i_pt) = shape_values[el_table(3)(i_pt)][i_pt % n_points][i_row];
644 template<
unsigned int dim>
647 auto &op = operations[scalar_shape_grads_op_idx];
648 uint n_points = ref_shape_grads[0].size();
649 uint n_dofs = ref_shape_grads[0][0].size();
651 Eigen::Vector<ArrayDbl, Eigen::Dynamic> ref_shape_grads_expd;
652 ref_shape_grads_expd.resize(dim*n_dofs);
653 for (
uint i=0; i<ref_shape_grads_expd.rows(); ++i)
654 ref_shape_grads_expd(i).resize(op_results(0).
rows());
656 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof)
657 for (
uint i_c=0; i_c<dim; ++i_c) {
658 ArrayDbl &shape_grad_row = ref_shape_grads_expd(i_dof*dim+i_c);
659 for (
uint i_pt=0; i_pt<shape_grad_row.rows(); ++i_pt)
660 shape_grad_row(i_pt) = ref_shape_grads[el_table(3)(i_pt)][i_pt % n_points][i_dof][i_c];
663 auto inv_jac_value = operations[ op.input_ops()[0] ].value<dim, 3>(op_results);
664 for (
uint i_dof=0; i_dof<n_dofs; ++i_dof) {
665 auto shape_grad_value = op.value<3, 1>(op_results, i_dof);
666 Eigen::Map<Eigen::Matrix<ArrayDbl, dim, 1>> ref_shape_grads_dof_value(ref_shape_grads_expd.data() + dim*i_dof, dim, 1);
667 shape_grad_value = inv_jac_value.transpose() * ref_shape_grads_dof_value;
677 template<
unsigned int spacedim = 3>
699 template<
unsigned int dim>
740 template<
unsigned int spacedim = 3>
762 template<
unsigned int dim>
Class represents FE operations.
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 of result_row_.
const std::vector< uint > & shape() const
Getter of shape_.
const std::vector< uint > & input_ops() const
Getter of input_ops_.
uint result_row_
First row to scalar, vector or matrix result.
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)
Constructor.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
uint dim() const
Getter of 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.
uint n_comp() const
Number of components computed from shape_ vector.
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.
Class for storing FE data of quadrature points.
ElOp< spacedim > & make_expansion(ElOp< spacedim > &el_op, std::initializer_list< uint > shape, ReinitFunction reinit_f)
Add accessor to operations_ vector.
uint add_rows(uint n_added)
Temporary method allow acces to old structure PatchFEValues::DimPatchFEValues.
void initialize(uint int_cols)
void resize_tables(uint n_points)
Resize data tables.
uint n_rows() const
Getter of 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)
Register side point, add to int_vals_ table.
uint n_rows_
Number of columns of point_vals table.
uint register_element(arma::mat coords, uint element_patch_idx)
Register element, add to point_vals_ table.
Tensor tensor_val(uint result_row, uint point_idx) const
uint register_side(arma::mat elm_coords, arma::mat side_coords)
Register side, add to point_vals_ table.
PatchPointValues(uint dim)
std::vector< ElOp< spacedim > > operations_
Vector of all defined operations.
Quadrature * get_quadrature() const
Return quadrature.
uint n_elems() const
Getter of n_elems_.
void reset()
Reset number of rows (points and elements)
uint n_points() const
Getter of n_points_.
std::vector< uint > elements_map_
Map of element patch indices to el_vals_ table.
uint dim() const
Getter of dim_.
void print(bool points, bool ints) const
Temporary development method.
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.
Vector vector_val(uint result_row, uint point_idx) const
void reinit_patch()
Reinit data.
uint n_points_
Number of points in patch.
Quadrature * quad_
Quadrature of given dimension and order passed in constructor.
Scalar scalar_val(uint result_row, uint point_idx) const
uint register_bulk_point(uint elem_table_row, uint value_patch_idx, uint elem_idx)
Register bulk point, add to int_vals_ table.
uint n_elems_
Number of elements in patch.
ElOp< spacedim > & make_new_op(std::initializer_list< uint > shape, ReinitFunction reinit_f, std::vector< uint > input_ops_vec)
Add accessor to operations_ vector.
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
BulkOps
enum of bulk operation
@ opExpansionJac
expands Jacobian
@ opCoords
operations evaluated on quadrature points
@ opInvJac
inverse Jacobian
@ opWeights
weight of quadrature point
@ opJac
Jacobian of element.
@ opExpansionInvJac
expands inverse Jacobian
@ opExpansionCoords
operation executed expansion to quadrature point (value of element > values on quadrature points)
@ opElCoords
operations evaluated on elements
@ opJxW
JxW value of quadrature point.
SideOps
enum of side operation
@ 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.
Definitions of particular quadrature rules on simplices.
Defines reinit operations on bulk points.
static void expd_jac_det(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, 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 elop_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 expd_jac(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_coords(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
static void expd_inv_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, 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 expd_sd_jac(std::vector< ElOp< 3 >> &operations, TableDbl &op_results, TableInt &el_table)
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 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)