19 #ifndef OP_FUNCTION_HH_
20 #define OP_FUNCTION_HH_
83 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
88 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, Domain::n_nodes(dim)}) {
89 this->
domain_ = Domain::domain();
94 uint n_elems = Domain::n_mesh_entities(
ppv);
98 for (
uint i_elm=0; i_elm<n_elems; ++i_elm)
99 for (
uint i_col=0; i_col<Domain::n_nodes(dim); ++i_col)
100 for (
uint i_row=0; i_row<spacedim; ++i_row) {
101 result(i_row, i_col)(i_elm) = ( *Domain::node(
ppv, i_elm, i_col) )(i_row);
108 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
113 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, Domain::n_nodes(dim)-1})
115 this->
domain_ = Domain::domain();
122 for (
unsigned int i=0; i<spacedim; i++)
123 for (
unsigned int j=0; j<Domain::n_nodes(dim)-1; j++)
124 jac_value(i,j) = coords_value(i,j+1) - coords_value(i,0);
129 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
134 :
PatchOp<spacedim>(dim, pfev, quad, {1})
136 this->
domain_ = Domain::domain();
153 :
PatchOp<3>(1, pfev, quad, {1})
160 uint n_sides =
ppv.n_mesh_items();
163 for (
uint i=0;i<n_sides; ++i) {
164 jac_det_value(0,0)(i) = 1.0;
173 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
178 :
PatchOp<spacedim>(dim, pfev, quad, {dim, spacedim})
180 this->
domain_ = Domain::domain();
187 inv_jac_value = eigen_arena_tools::inverse<spacedim, dim>(jac_value);
195 template<
unsigned int spacedim = 3>
214 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
219 :
PatchOp<spacedim>(dim, pfev, quad, {dim+1})
221 this->
domain_ = Domain::domain();
222 double n_points = quad->size();
224 for (
uint i_p=0; i_p<n_points; ++i_p) {
226 for (
uint i_c=0; i_c<dim+1; ++i_c)
227 this->
result_(i_c)(i_p) = ref_bar_coords[i_c];
235 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
240 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim})
242 this->
domain_ = Domain::domain();
253 uint n_elems =
ppv.n_mesh_items();
254 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> coords_vec(spacedim, Domain::n_nodes(dim));
255 Eigen::Matrix<ArenaOVec<double>, spacedim, Domain::n_nodes(dim)> coords_ovec;
256 for (
uint i=0; i<spacedim*Domain::n_nodes(dim); ++i) {
263 Eigen::Vector<ArenaOVec<double>, Eigen::Dynamic> ref_bary_ovec(Domain::n_nodes(dim));
264 for (
uint i=0; i<Domain::n_nodes(dim); ++i) {
269 Eigen::Vector<ArenaOVec<double>, Eigen::Dynamic> result_ovec = coords_ovec * ref_bary_ovec;
270 for (
uint i=0; i<spacedim; ++i) {
271 result_vec(i) = result_ovec(i).get_vec();
280 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
285 :
PatchOp<spacedim>(dim, pfev, quad, {1})
287 this->
domain_ = Domain::domain();
291 for (
uint i=0; i<point_weights_vec.size(); ++i)
292 this->
result_(0)(i) = point_weights_vec[i];
302 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
307 :
PatchOp<spacedim>(dim, pfev, quad, {1})
309 this->
domain_ = Domain::domain();
320 uint n_elems =
ppv.n_mesh_items();
331 template<
unsigned int dim,
unsigned int spacedim>
336 :
PatchOp<spacedim>(dim, pfev, quad, {1})
354 template<
unsigned int dim,
unsigned int spacedim = 3>
359 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim})
371 uint n_sides =
ppv.n_mesh_items();
372 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_value(dim, spacedim);
373 for (
uint i=0; i<dim*spacedim; ++i) {
376 for (
uint i_c=0; i_c<dim*spacedim; ++i_c) {
383 Eigen::VectorXd A(3);
384 for (
uint i=0; i<normal_value(0).data_size(); ++i) {
385 A(0) = normal_value(0)(i);
386 A(1) = normal_value(1)(i);
387 A(2) = normal_value(2)(i);
388 norm_vec(i) = A.norm();
391 for (
uint i=0; i<spacedim; ++i) {
392 normal_value(i) = normal_value(i) / norm_vec;
398 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
403 :
PatchOp<spacedim>(dim, pfev, quad, {1}, fe->n_dofs())
405 this->
domain_ = Domain::domain();
406 uint n_points = quad->size();
410 for (
unsigned int i_p = 0; i_p < n_points; i_p++)
411 for (
unsigned int i_dof = 0; i_dof < this->
n_dofs_; i_dof++)
412 ref_scalar_value(i_dof)(i_p) = fe->shape_value(i_dof, quad->point<dim>(i_p));
419 template<
unsigned int dim,
unsigned int spacedim>
424 :
PatchOp<spacedim>(dim, pfev, quad, {dim+1}, fe->n_dofs())
427 uint n_points = quad->size();
431 for (
unsigned int s=0; s<dim+1; ++s) {
433 for (
unsigned int i_p = 0; i_p < n_points; i_p++)
434 for (
unsigned int i_dof = 0; i_dof < this->
n_dofs_; i_dof++)
435 ref_scalar_value(s, i_dof)(i_p) = fe->shape_value(i_dof, side_q.
point<dim>(i_p));
443 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
448 :
PatchOp<spacedim>(dim, pfev, quad, {fe->n_components()}, fe->n_dofs())
450 this->
domain_ = Domain::domain();
451 uint n_points = quad->size();
457 for (
uint i_p=0; i_p<n_points; ++i_p)
458 for (
uint c=0; c<fe->n_components(); ++c)
459 ref_vector_value(c, i_dof)(i_p) = fe->shape_value(i_dof, quad->point<dim>(i_p), c);
466 template<
unsigned int dim,
unsigned int spacedim>
471 :
PatchOp<spacedim>(dim, pfev, quad, {dim+1, fe->n_components()}, fe->n_dofs())
474 uint n_points = quad->size();
479 for (
unsigned int s=0; s<dim+1; ++s) {
481 for (
unsigned int i_p = 0; i_p < n_points; i_p++)
482 for (
unsigned int i_dof = 0; i_dof < this->
n_dofs_; i_dof++)
483 for (
uint c=0; c<fe->n_components(); ++c)
484 ref_vector_value(s,fe->n_components()*i_dof+c)(i_p) = fe->shape_value(i_dof, side_q.
point<dim>(i_p), c);
492 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
497 :
PatchOp<spacedim>(dim, pfev, quad, {dim, 1}, fe->n_dofs())
499 this->
domain_ = Domain::domain();
500 uint n_points = quad->size();
504 for (
uint i_row=0; i_row<ref_scalar_value.rows(); ++i_row)
506 for (
uint i_p=0; i_p<n_points; ++i_p)
507 ref_scalar_value(i_row, i_dof)(i_p) = fe->shape_grad(i_dof, quad->point<dim>(i_p))[i_row];
514 template<
unsigned int dim,
unsigned int spacedim>
519 :
PatchOp<spacedim>(dim, pfev, quad, {dim+1, dim}, fe->n_dofs())
522 uint n_points = quad->size();
526 for (
unsigned int s=0; s<dim+1; ++s) {
529 for (
uint i_p=0; i_p<n_points; ++i_p)
530 for (
uint c=0; c<dim; ++c)
531 ref_scalar_value(s,dim*i_dof+c)(i_p) = fe->shape_grad(i_dof, side_q.
point<dim>(i_p))[c];
539 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
544 :
PatchOp<spacedim>(dim, pfev, quad, {dim, fe->n_components()}, fe->n_dofs())
546 this->
domain_ = Domain::domain();
547 uint n_points = quad->size();
553 for (
uint i_dim=0; i_dim<dim; ++i_dim)
555 for (
uint i_p=0; i_p<n_points; ++i_p)
556 ref_vector_value(i_dim,
n_comp*i_dof+i_c)(i_p) = fe->shape_grad(i_dof, quad->point<dim>(i_p), i_c)[i_dim];
564 template<
unsigned int dim,
unsigned int spacedim>
569 :
PatchOp<spacedim>(dim, pfev, quad, {(dim+1)*dim, fe->n_components()}, fe->n_dofs())
572 uint n_points = quad->size();
573 uint n_sides = dim+1;
578 for (
uint i_sd=0; i_sd<n_sides; ++i_sd) {
581 for (
uint i_dim=0; i_dim<dim; ++i_dim)
583 for (
uint i_p=0; i_p<n_points; ++i_p)
584 ref_vector_value(i_sd*dim+i_dim,
n_comp*i_dof+i_c)(i_p) = fe->shape_grad(i_dof, side_q.
point<dim>(i_p), i_c)[i_dim];
592 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
597 :
PatchOp<spacedim>(dim, pfev, quad, {1}, fe->n_dofs())
600 this->
domain_ = Domain::domain();
609 uint n_elem = this->
ppv().n_mesh_items();
612 for (
uint i=0; i<n_elem; ++i) {
617 Eigen::Vector<ArenaOVec<double>, Eigen::Dynamic> ref_ovec(
n_dofs);
622 Eigen::Vector<ArenaOVec<double>, Eigen::Dynamic> result_ovec = elem_ovec * ref_ovec;
624 result_vec(i) = result_ovec(i).
get_vec();
630 template<
unsigned int dim,
unsigned int spacedim>
635 :
PatchOp<spacedim>(dim, pfev, quad, {1}, fe->n_dofs())
645 uint n_sides =
ppv.n_mesh_items();
654 for (
uint i_pt=0; i_pt<n_patch_points; ++i_pt) {
655 result_vec(i_dof)(i_pt) = ref_vec(
ppv.int_table_(
ref_side_on_quads)(i_pt), i_dof)(i_pt / n_sides);
662 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
669 this->
domain_ = Domain::domain();
678 uint n_elem = this->
ppv().n_mesh_items();
681 for (
uint i=0; i<n_elem; ++i) {
686 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> ref_shape_ovec(spacedim,
n_dofs);
688 ref_shape_ovec(c) =
ArenaOVec(ref_shape_vec(c));
691 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> result_ovec = elem_ovec * ref_shape_ovec;
693 result_vec(c) = result_ovec(c).
get_vec();
701 template<
unsigned int dim,
unsigned int spacedim>
715 uint n_sides =
ppv.n_mesh_items();
727 for (
uint i_pt=0; i_pt<n_patch_points; ++i_pt)
728 for (
uint c=0; c<spacedim; c++)
729 result_vec(c,i_dof)(i_pt) = ref_shape_vec(
ppv.int_table_(
ref_side_on_quads)(i_pt),3*i_dof+c)(i_pt / n_sides);
738 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
745 this->
domain_ = Domain::domain();
761 uint n_elems =
ppv.n_mesh_items();
762 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> jac_vec(spacedim, dim);
763 for (
uint i=0; i<spacedim*dim; ++i) {
766 for (
uint i_c=0; i_c<spacedim*dim; ++i_c) {
772 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> jac_div_det_vec = jac_vec / jac_det_vec;
773 Eigen::Matrix<ArenaOVec<double>, spacedim, dim> jac_div_det_ovec;
774 for (
uint c=0; c<spacedim*dim; ++c) {
778 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> ref_shape_ovec(dim,
n_dofs);
783 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> result_ovec = jac_div_det_ovec * ref_shape_ovec;
785 result_vec(c) = result_ovec(c).get_vec();
794 template<
unsigned int dim,
unsigned int spacedim>
810 uint n_sides =
ppv.n_mesh_items();
821 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> jac_vec(spacedim, dim);
822 for (
uint i=0; i<spacedim*dim; ++i) {
825 for (
uint i_c=0; i_c<spacedim*dim; ++i_c) {
832 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> jac_div_det_vec = jac_vec / jac_det_vec;
833 Eigen::Matrix<ArenaOVec<double>, spacedim, dim> jac_div_det_ovec;
834 for (
uint c=0; c<spacedim*dim; ++c) {
835 jac_div_det_ovec(c) =
ArenaOVec(jac_div_det_vec(c));
840 for (
uint i=0; i<n_points_per_side; ++i) {
841 side_points_vec(i) = 1.0;
844 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> expand_jac_div_det_ovec = jac_div_det_ovec * side_points_ovec;
845 Eigen::Matrix<ArenaVec<double>, spacedim, dim> expand_jac_div_det_vec;
846 for (
uint c=0; c<spacedim*this->
dim_; ++c) {
847 expand_jac_div_det_vec(c) = expand_jac_div_det_ovec(c).
get_vec();
850 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> expand_ref_vec(dim,
n_dofs);
855 for (
uint i_pt=0; i_pt<n_points_per_side; ++i_pt) {
856 uint i_begin = i_pt * n_sides;
857 for (
uint i_sd=0; i_sd<n_sides; ++i_sd) {
858 for (
uint i_c=0; i_c<dim; ++i_c) {
859 expand_ref_vec(i_c, i_dof)(i_begin + i_sd) = ref_shape_vec(
ppv.int_table_(
ref_side_on_sides)(i_sd), i_dof*dim+i_c)(i_pt);
866 result_vec = expand_jac_div_det_vec * expand_ref_vec;
877 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
882 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim}, fe->n_dofs()),
in_op_(
nullptr)
884 this->
domain_ = Domain::domain();
885 switch (fe->fe_type()) {
893 ASSERT_PERMANENT(
false).error(
"Shape vector for FEVectorContravariant is not implemented yet!\n");
899 in_op_ =
new VectorPiolaShape<dim, Domain, spacedim>(pfev, quad, fe, *
this);
903 ASSERT(
false).error(
"Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
917 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
922 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, 1}, fe->n_dofs())
925 this->
domain_ = Domain::domain();
938 uint n_elems =
ppv.n_mesh_items();
939 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_vec(dim, spacedim);
940 for (
uint i=0; i<dim*spacedim; ++i) {
943 for (
uint i_c=0; i_c<dim*spacedim; ++i_c) {
947 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> ref_grads_ovec(this->
dim_, n_dofs);
949 ref_grads_ovec(i) =
ArenaOVec(ref_grads_vec(i));
952 Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
953 for (
uint i=0; i<this->
dim_*spacedim; ++i) {
954 inv_jac_ovec(i) =
ArenaOVec(inv_jac_vec(i));
958 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
960 result_vec(i) = result_ovec(i).get_vec();
966 template<
unsigned int dim,
unsigned int spacedim>
971 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, 1}, fe->n_dofs())
985 uint n_points = ref_shape_grads(0).data_size();
986 uint n_sides =
ppv.n_mesh_items();
993 Eigen::Matrix<ArenaVec<double>, dim, 3> inv_jac_expd_value;
994 for (
uint i=0; i<dim*3; ++i) {
996 for (
uint j=0; j<n_patch_points; ++j)
1001 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> ref_shape_grads_expd(dim,
n_dofs);
1006 for (
uint i_pt=0; i_pt<n_points; ++i_pt) {
1007 uint i_begin = i_pt * n_sides;
1008 for (
uint i_sd=0; i_sd<n_sides; ++i_sd) {
1009 for (
uint i_c=0; i_c<dim; ++i_c) {
1010 ref_shape_grads_expd(i_c, i_dof)(i_begin + i_sd) = ref_shape_grads(
ppv.int_table_(
ref_side_on_sides)(i_sd), i_dof*dim+i_c)(i_pt);
1017 grad_scalar_shape_value = inv_jac_expd_value.transpose() * ref_shape_grads_expd;
1022 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
1027 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs()),
dispatch_op_(dispatch_op)
1029 this->
domain_ = Domain::domain();
1043 uint n_elems =
ppv.n_mesh_items();
1044 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_vec(dim, spacedim);
1045 for (
uint i=0; i<dim*spacedim; ++i) {
1048 for (
uint i_c=0; i_c<dim*spacedim; ++i_c) {
1052 Eigen::Matrix<ArenaOVec<double>, dim, 3> inv_jac_ovec;
1053 for (
uint i=0; i<dim*spacedim; ++i) {
1054 inv_jac_ovec(i) =
ArenaOVec(inv_jac_vec(i));
1057 Eigen::Matrix<ArenaOVec<double>, dim, spacedim> ref_grads_ovec;
1059 for (
uint i=0; i<spacedim*dim; ++i) {
1060 ref_grads_ovec(i) =
ArenaOVec(ref_grads_vec(i_dof*3*dim + i));
1063 Eigen::Matrix<ArenaOVec<double>, spacedim, spacedim> result_ovec = inv_jac_ovec.transpose() * ref_grads_ovec;
1064 for (
uint i=0; i<spacedim; ++i) {
1065 for (
uint j=0; j<spacedim; ++j) {
1066 result_vec(j,i+spacedim*i_dof) = result_ovec(i,j).get_vec();
1077 template<
unsigned int dim,
unsigned int spacedim>
1082 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs()),
dispatch_op_(dispatch_op)
1095 uint n_points = ref_vector_grad(0).data_size();
1096 uint n_patch_sides =
ppv.n_mesh_items();
1100 Eigen::Matrix<ArenaVec<double>, dim, 3> inv_jac_expd_value;
1101 for (
uint i=0; i<dim*3; ++i) {
1103 for (
uint j=0; j<n_patch_points; ++j)
1108 Eigen::Matrix<ArenaVec<double>, dim, 3> ref_shape_grads_expd;
1109 for (
uint i=0; i<spacedim*dim; ++i) {
1114 for (
uint i_pt=0; i_pt<n_points; ++i_pt) {
1115 uint i_begin = i_pt * n_patch_sides;
1116 for (
uint i_sd=0; i_sd<n_patch_sides; ++i_sd) {
1117 for (
uint i_dim=0; i_dim<dim; ++i_dim) {
1118 for (
uint i_c=0; i_c<spacedim; ++i_c) {
1119 ref_shape_grads_expd(i_dim, i_c)(i_begin + i_sd) = ref_vector_grad(
ppv.int_table_(
ref_side_on_sides)(i_sd)*dim+i_dim, 3*i_dof+i_c)(i_pt);
1127 res_submat = (inv_jac_expd_value.transpose() * ref_shape_grads_expd).transpose();
1136 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
1141 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs()),
dispatch_op_(dispatch_op)
1143 this->
domain_ = Domain::domain();
1158 uint n_elems =
ppv.n_mesh_items();
1163 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_vec(dim, spacedim);
1164 for (
uint i=0; i<dim*spacedim; ++i) {
1167 for (
uint i_c=0; i_c<spacedim*dim; ++i_c) {
1173 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_div_det_vec = inv_jac_vec / jac_det_vec;
1174 Eigen::Matrix<ArenaOVec<double>, dim, spacedim> inv_jac_div_det_ovec;
1175 for (
uint c=0; c<dim*spacedim; ++c) {
1179 Eigen::Matrix<ArenaVec<double>, spacedim, dim> jac_vec;
1180 for (
uint i_c=0; i_c<spacedim*dim; ++i_c) {
1182 for (
uint i_el=0; i_el<n_elems; ++i_el)
1183 for (
uint i_pt=0; i_pt<n_ref_points; ++i_pt) {
1188 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> ref_grads_ovec(dim, dim);
1190 for (
uint i=0; i<dim*dim; ++i) {
1191 ref_grads_ovec(i) =
ArenaOVec(ref_grads_vec(i_dof*dim*dim + i));
1194 Eigen::Matrix<ArenaOVec<double>, spacedim, dim> inv_jac_multi_ref_ovec = inv_jac_div_det_ovec.transpose() * ref_grads_ovec;
1195 Eigen::Matrix<ArenaVec<double>, spacedim, dim> inv_jac_multi_ref_vec;
1196 for (
uint i=0; i<spacedim*dim; ++i) {
1197 inv_jac_multi_ref_vec(i) = inv_jac_multi_ref_ovec(i);
1200 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> sub_result_vec = inv_jac_multi_ref_vec * jac_vec.transpose();
1201 for (
uint i=0; i<spacedim; ++i) {
1202 for (
uint j=0; j<spacedim; ++j) {
1203 result_vec(j,i+spacedim*i_dof) = sub_result_vec(i,j);
1214 template<
unsigned int dim,
unsigned int spacedim>
1219 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs()),
dispatch_op_(dispatch_op)
1235 uint n_sides =
ppv.n_mesh_items();
1237 uint n_points = ref_grads_vec(0).data_size();
1242 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_vec(dim, spacedim);
1243 for (
uint i_c=0; i_c<dim*spacedim; ++i_c) {
1247 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> jac_vec(spacedim, dim);
1248 Eigen::Matrix<ArenaOVec<double>, spacedim, dim> jac_ovec;
1249 for (
uint i_c=0; i_c<spacedim*dim; ++i_c) {
1258 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> inv_jac_div_det_vec = inv_jac_vec / jac_det_vec;
1259 Eigen::Matrix<ArenaOVec<double>, dim, spacedim> inv_jac_div_det_ovec;
1260 for (
uint c=0; c<dim*spacedim; ++c) {
1266 for (
uint i=0; i<n_points_per_side; ++i) {
1267 side_points_vec(i) = 1.0;
1270 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> expand_inv_jac_div_det_ovec = inv_jac_div_det_ovec * side_points_ovec;
1271 Eigen::Matrix<ArenaVec<double>, dim, spacedim> expand_inv_jac_div_det_vec;
1272 for (
uint c=0; c<spacedim*dim; ++c) {
1273 expand_inv_jac_div_det_vec(c) = expand_inv_jac_div_det_ovec(c).
get_vec();
1275 Eigen::Matrix<ArenaOVec<double>, Eigen::Dynamic, Eigen::Dynamic> expand_jac_ovec = jac_ovec * side_points_ovec;
1276 Eigen::Matrix<ArenaVec<double>, spacedim, dim> expand_jac_vec;
1277 for (
uint c=0; c<spacedim*dim; ++c) {
1278 expand_jac_vec(c) = expand_jac_ovec(c).
get_vec();
1281 Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> expand_ref_vec(dim, dim);
1282 for (
uint c=0; c<dim*dim; c++)
1287 for (
uint i_pt=0; i_pt<n_points; ++i_pt) {
1288 uint i_begin = i_pt * n_sides;
1289 for (
uint i_sd=0; i_sd<n_sides; ++i_sd) {
1290 for (
uint i_dim=0; i_dim<dim; ++i_dim) {
1291 for (
uint i_c=0; i_c<dim; ++i_c) {
1292 expand_ref_vec(i_dim, i_c)(i_begin + i_sd) = ref_grads_vec(
ppv.int_table_(
ref_side_on_sides)(i_sd)*dim+i_dim, dim*i_dof+i_c)(i_pt);
1300 res_submat = expand_inv_jac_div_det_vec.transpose() * expand_ref_vec * expand_jac_vec.transpose();
1311 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
1316 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs()),
in_op_(
nullptr)
1318 this->
domain_ = Domain::domain();
1319 switch (fe->fe_type()) {
1327 ASSERT_PERMANENT(
false).error(
"Shape vector for FEVectorContravariant is not implemented yet!\n");
1333 in_op_ =
new GradVectorPiolaShape<dim, Domain, spacedim>(pfev, quad, fe, *
this);
1337 ASSERT(
false).error(
"Type of FiniteElement of grad_vector_shape accessor must be FEVector, FEVectorPiola or FEVectorContravariant!\n");
1351 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
1356 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs())
1358 this->
domain_ = Domain::domain();
1363 for (
uint i_dof=0; i_dof<this->
n_dofs(); ++i_dof) {
1364 Eigen::Map< Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> > grad_vector_dof = this->
input_ops(0)->
result_sub_matrix(i_dof);
1365 Eigen::Map< Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> > sym_grad_dof = this->
result_sub_matrix(i_dof);
1366 sym_grad_dof = 0.5 * (grad_vector_dof.transpose() + grad_vector_dof);
1372 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
1377 :
PatchOp<spacedim>(dim, pfev, quad, {1}, fe->n_dofs())
1379 this->
domain_ = Domain::domain();
1386 for (
uint i_dof=0; i_dof<this->
n_dofs(); ++i_dof) {
1387 Eigen::Map< Eigen::Matrix<ArenaVec<double>, Eigen::Dynamic, Eigen::Dynamic> > grad_vector_dof = this->
input_ops(0)->
result_sub_matrix(i_dof);
1388 divergence_value(i_dof) = grad_vector_dof(0,0) + grad_vector_dof(1,1) + grad_vector_dof(2,2);
1394 template<
unsigned int dim,
class Domain,
unsigned int spacedim = 3>
1399 :
PatchOp<spacedim>(dim, pfev, quad, {spacedim, spacedim}, fe->n_dofs())
1401 this->
domain_ = Domain::domain();
#define ASSERT_PERMANENT(expr)
Allow use shorter versions of macro names if these names is not used with external library.
#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.
Abstract class for the description of a general finite element on a reference simplex in dim dimensio...
Class used as template type for type resolution Bulk / Side.
static ElemDomain domain()
static constexpr uint n_nodes(uint dim)
static NodeAccessor< 3 > node(PatchPointValues< 3 > &ppv, unsigned int i_elm, unsigned int i_n)
Return i_n-th node of i_elm-th element stored in PatchPointValues::elem_dim_list_.
static uint n_mesh_entities(PatchPointValues< 3 > &ppv)
Return number of mesh entities (in this case elements) on patch.
void eval() override
Reinit function of operation. Implementation in descendants.
Coords(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
Dispatch class of vector values.
PatchOp< spacedim > * in_op_
void eval() override
Reinit function of operation. Implementation in descendants.
DispatchGradVectorShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Dispatch class of vector values.
void eval() override
Reinit function of operation. Implementation in descendants.
PatchOp< spacedim > * in_op_
DispatchVectorShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
FuncHelper()
Forbidden constructor.
static void fill_reduce_element_data_vec(PatchPointValues< spacedim > &ppv, ArenaVec< double > &source, ArenaVec< double > &target)
void eval() override
Reinit function of operation. Implementation in descendants.
GradScalarShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Evaluates gradient scalar values.
GradScalarShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
void eval() override
Reinit function of operation. Implementation in descendants.
GradVectorPiolaShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
PatchOp< spacedim > & dispatch_op_
Evaluates vector values (FEType == FEVectorPiola)
void eval() override
Reinit function of operation. Implementation in descendants.
GradVectorPiolaShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
PatchOp< spacedim > & dispatch_op_
void eval() override
Reinit function of operation. Implementation in descendants.
PatchOp< spacedim > & dispatch_op_
GradVectorShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
Evaluates vector values (FEType == FEVector)
GradVectorShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
PatchOp< spacedim > & dispatch_op_
void eval() override
Reinit function of operation. Implementation in descendants.
void eval() override
Reinit function of operation. Implementation in descendants.
InvJac(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
JacDet(PatchFEValues< 3 > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Evaluates Jacobian determinants on Bulk (Element) / Side.
JacDet(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Evaluates Jacobians on Bulk (Element) / Side.
Jac(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
JxW(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
JxW(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Evaluates normal vector on side quadrature points.
void eval() override
Reinit function of operation. Implementation in descendants.
NormalVec(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
Class represents zero operation of Join quantities.
void eval() override
Reinit function of operation. Implementation in descendants.
OpZero(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Evaluates coordinates of quadrature points - not implemented yet.
PtCoords(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
void eval() override
Reinit function of operation. Implementation in descendants.
RefBaryCoords(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
Template specialization of previous: Domain=SideDomain.
void eval() override
Reinit function of operation. Implementation in descendants.
RefGradScalar(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Fixed operation of gradient scalar reference values.
void eval() override
Reinit function of operation. Implementation in descendants.
RefGradScalar(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Template specialization of previous: Domain=SideDomain.
RefGradVector(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Fixed operation of gradient scalar reference values.
void eval() override
Reinit function of operation. Implementation in descendants.
RefGradVector(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Template specialization of previous: Domain=SideDomain.
RefScalar(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Fixed operation of scalar shape reference values.
RefScalar(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Template specialization of previous: Domain=SideDomain.
RefVector(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Fixed operation of vector shape reference values.
void eval() override
Reinit function of operation. Implementation in descendants.
RefVector(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
ScalarShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Evaluates scalar shape values.
void eval() override
Reinit function of operation. Implementation in descendants.
ScalarShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Class used as template type for type resolution Bulk / Side.
static constexpr uint n_nodes(uint dim)
static uint n_mesh_entities(PatchPointValues< 3 > &ppv)
Return number of mesh entities (in this case sides) on patch.
static ElemDomain domain()
static NodeAccessor< 3 > node(PatchPointValues< 3 > &ppv, unsigned int i_elm, unsigned int i_n)
Return i_n-th node of i_elm-th side stored in PatchPointValues::side_list_.
Evaluates vector vector divergence.
void eval() override
Reinit function of operation. Implementation in descendants.
VectorDivergence(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
VectorPiolaShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
PatchOp< spacedim > & dispatch_op_
Evaluates vector values (FEType == FEVectorPiola)
void eval() override
Reinit function of operation. Implementation in descendants.
VectorPiolaShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
PatchOp< spacedim > & dispatch_op_
PatchOp< spacedim > & dispatch_op_
VectorShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Evaluates vector values (FEType == FEVector)
PatchOp< spacedim > & dispatch_op_
void eval() override
Reinit function of operation. Implementation in descendants.
VectorShape(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe, PatchOp< spacedim > &dispatch_op)
Constructor.
Evaluates vector symmetric gradients.
void eval() override
Reinit function of operation. Implementation in descendants.
VectorSymGrad(PatchFEValues< spacedim > &pfev, const Quadrature *quad, std::shared_ptr< FiniteElement< dim >> fe)
Constructor.
Weights(PatchFEValues< spacedim > &pfev, const Quadrature *quad)
Constructor.
void eval() override
Reinit function of operation. Implementation in descendants.
Class represents element or FE operations.
Eigen::Vector< ArenaVec< double >, Eigen::Dynamic > result_
Result matrix of operation.
Eigen::Map< Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > > result_sub_matrix(uint i_dof)
Return map referenced result of DOF values as Eigen::Matrix.
const Quadrature * quad_
Pointer to Quadrature.
uint n_dofs() const
Getter for n_dofs_.
std::vector< PatchOp< spacedim > * > input_ops_
Indices of operations in PatchPointValues::operations_ vector on which PatchOp is depended.
PatchOp< spacedim > * input_ops(uint i_op) const
Return pointer to operation of i_op index in input operation vector.
Eigen::Map< Eigen::Matrix< ArenaVec< double >, Eigen::Dynamic, Eigen::Dynamic > > result_matrix()
Return map referenced result as Eigen::Vector.
void allocate_const_result(ArenaVec< double > &value_vec)
ElemDomain domain_
Flag: BulkOp = 0, SideOp = 1.
void allocate_result(size_t data_size, PatchArena &arena)
PatchPointValues< spacedim > & ppv() const
Return reference of PatchPointValues.
virtual void eval()=0
Reinit function of operation. Implementation in descendants.
uint n_dofs_
Number of DOFs of FE operations (or 1 in case of element operations)
PatchFEValues< spacedim > * patch_fe_
Pointer to PatchFEValues object.
uint n_mesh_items() const
Getter for n_mesh_items__.
ElemDimList< spacedim > * elem_dim_list_
Number and list of elements on patch.
std::vector< Side > side_list_
List of sides on patch.
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.
static BaryPoint local_to_bary(const LocalPoint &lp)
Converts from local to barycentric coordinates.
static Eigen::Matrix< ArenaVec< double >, dim, 1 > normal_vector_array(ArenaVec< uint > loc_side_idx_vec)
Class FEValues calculates finite element data on the actual cells such as shape function values,...
Base class of FE operations.
ElemDomain
Distinguishes bulk and side domain.
@ patch_elem_on_domain
Index of patch element for each bulk element or side.
@ ref_side_on_sides
Ref index of side in element for each side in patch.
@ ref_side_on_quads
Ref index of side in element for each quadrature point in patch.