7#ifndef __SURFACE_SLIDING_CONSTRAINS_HPP__
8#define __SURFACE_SLIDING_CONSTRAINS_HPP__
25 boost::shared_ptr<VectorDouble> &active_variables_ptr)
29 MoFEMErrorCode
doWork(
int side, EntityType type,
30 EntitiesFieldData::EntData &data) {
32 if (type == MBVERTEX) {
33 for (
unsigned int dd = 0; dd != data.getFieldData().size(); ++dd)
40 template <
int SizeLambda>
41 struct OpGetActiveDofsPositions
46 boost::shared_ptr<VectorDouble> &active_variables_ptr)
50 MoFEMErrorCode
doWork(
int side, EntityType type,
51 EntitiesFieldData::EntData &data) {
53 if (type == MBVERTEX) {
54 for (
unsigned int dd = 0; dd != data.getFieldData().size(); ++dd)
61 template <
int SizeLambda,
int SizePositions>
67 boost::shared_ptr<VectorDouble> &results_ptr)
72 MoFEMErrorCode
doWork(
int side, EntityType type,
73 EntitiesFieldData::EntData &data) {
77 VectorInt &indices = data.getIndices();
81 else if (indices.size() == SizeLambda)
83 else if (indices.size() == SizePositions)
87 "Element %s: Data inconsistency nb of indices %ld",
90 CHKERR VecSetOption(
getSNESf(), VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
96 template <
int SizeLambda,
int SizePositions>
102 const std::string field_name_col,
103 boost::shared_ptr<MatrixDouble> &jacobian_ptr)
110 MoFEMErrorCode
doWork(
int row_side,
int col_side, EntityType row_type,
112 EntitiesFieldData::EntData &row_data,
113 EntitiesFieldData::EntData &col_data) {
115 if (row_type != MBVERTEX)
117 if (col_type != MBVERTEX)
119 VectorInt &row_indices = row_data.getIndices();
120 VectorInt &col_indices = col_data.getIndices();
121 if (row_indices.empty() || col_indices.empty())
125 if (row_indices.size() == SizeLambda)
127 else if (row_indices.size() == SizePositions)
128 shift_row = SizeLambda;
131 "Data inconsistency");
134 if (col_indices.size() == SizeLambda)
136 else if (col_indices.size() == SizePositions)
137 shift_col = SizeLambda;
140 "Data inconsistency");
142 MatrixDouble jac(row_indices.size(), col_indices.size());
143 for (
unsigned int rr = 0; rr != row_indices.size(); ++rr) {
144 for (
unsigned int cc = 0; cc != col_indices.size(); ++cc) {
145 jac(rr, cc) = (*jacobianPtr)(shift_row + rr, shift_col + cc);
146 if (jac(rr, cc) != jac(rr, cc))
148 "Jacobian assemble not a number jac(rr,cc) = %3.4f",
152 CHKERR MatSetValues(
getSNESB(), row_indices.size(), &row_indices[0],
153 col_indices.size(), &col_indices[0], &jac(0, 0),
334 struct DriverElementOrientation {
370 if (
B != PETSC_NULLPTR) {
374 if (
F != PETSC_NULLPTR) {
394 PetscOptionsBegin(PETSC_COMM_WORLD,
"",
395 "Get surface sliding constrains element scaling",
397 CHKERR PetscOptionsScalar(
"-surface_sliding_alpha",
"scaling parameter",
"",
410 CHKERRABORT(PETSC_COMM_WORLD,
ierr);
428 boost::shared_ptr<VectorDouble> &active_variables_ptr,
429 boost::shared_ptr<VectorDouble> &results_ptr,
430 boost::shared_ptr<MatrixDouble> &jacobian_ptr,
440 MoFEMErrorCode
doWork(
int side, EntityType type,
441 EntitiesFieldData::EntData &data) {
443 if (type != MBVERTEX)
458 int nb_gauss_pts = data.getN().size1();
459 int nb_base_functions = data.getN().size2();
460 auto t_base1 = data.getFTensor0N();
461 auto t_base2 = data.getFTensor0N();
462 auto t_diff_base = data.getFTensor1DiffN<2>();
470 ublas::vector<adouble> c_vec(3);
471 ublas::vector<adouble> f_vec(9);
477 ublas::vector<adouble> lambda_dofs(3);
478 for (
int dd = 0; dd != 3; ++dd) {
479 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
481 ublas::vector<adouble> position_dofs(9);
482 for (
int dd = 0; dd != 9; ++dd) {
483 position_dofs[dd] <<= (*activeVariablesPtr)[3 + dd];
488 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
491 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
496 t_position_ksi(
i) = 0;
497 t_position_eta(
i) = 0;
500 for (
int bb = 0; bb != nb_base_functions; ++bb) {
501 t_position(
i) += t_base1 * t_position_dofs(
i);
502 t_position_ksi(
i) += t_diff_base(N0) * t_position_dofs(
i);
503 t_position_eta(
i) += t_diff_base(N1) * t_position_dofs(
i);
504 lambda += t_base1 * t_lambda_dof;
511 t_delta(
i) = t_position(
i) - t_coord_ref(
i);
516 &f_vec[0], &f_vec[1], &f_vec[2]);
518 double w = t_w * (0.5 *
aLpha * eo);
520 for (
int bb = 0; bb != nb_base_functions; ++bb) {
521 double alpha = w * t_base2;
522 t_c += alpha * t_normal(
i) * t_delta(
i);
523 t_f(
i) += alpha * (
lambda * t_normal(
i));
533 for (
int rr = 0; rr != 3; ++rr)
536 for (
int rr = 0; rr != 9; ++rr)
537 f_vec(rr) >>= (*resultsPtr)[3 + rr];
542 double *jac_ptr[3 + 9];
543 for (
int rr = 0; rr != 3 + 9; ++rr)
551 "ADOL-C function evaluation with error");
559 const std::string lagrange_multipliers_field_name,
560 const std::string material_field_name,
561 const double *alpha =
nullptr) {
564 if (alpha !=
nullptr) {
568 boost::shared_ptr<VectorDouble> active_variables_ptr(
569 new VectorDouble(3 + 9));
570 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(3 + 9));
571 boost::shared_ptr<MatrixDouble> jacobian_ptr(
572 new MatrixDouble(3 + 9, 3 + 9));
576 lagrange_multipliers_field_name, active_variables_ptr));
578 material_field_name, active_variables_ptr));
580 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
590 lagrange_multipliers_field_name, active_variables_ptr));
592 material_field_name, active_variables_ptr));
594 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
597 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
599 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
601 material_field_name, material_field_name, jacobian_ptr));
608 const std::string lagrange_multipliers_field_name,
609 const std::string material_field_name) {
612 boost::shared_ptr<VectorDouble> active_variables_ptr(
613 new VectorDouble(3 + 9));
614 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(3 + 9));
615 boost::shared_ptr<MatrixDouble> jacobian_ptr(
616 new MatrixDouble(3 + 9, 3 + 9));
621 lagrange_multipliers_field_name, active_variables_ptr));
623 material_field_name, active_variables_ptr));
625 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
628 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
636 struct CalculateEdgeBase {
641 double def_val[] = {0, 0, 0};
642 CHKERR moab.tag_get_handle(
"EDGE_BASE0", 3, MB_TYPE_DOUBLE, th0,
643 MB_TAG_CREAT | MB_TAG_SPARSE, def_val);
644 CHKERR moab.tag_get_handle(
"EDGE_BASE1", 3, MB_TYPE_DOUBLE, th1,
645 MB_TAG_CREAT | MB_TAG_SPARSE, def_val);
646 int def_numb_val[] = {-1};
647 CHKERR moab.tag_get_handle(
"PATCH_NUMBER", 1, MB_TYPE_INTEGER, th2,
648 MB_TAG_CREAT | MB_TAG_SPARSE, def_numb_val);
649 int def_orientation_val[] = {1};
650 CHKERR moab.tag_get_handle(
"PATCH_ORIENTATION", 1, MB_TYPE_INTEGER, th3,
651 MB_TAG_CREAT | MB_TAG_SPARSE,
652 def_orientation_val);
661 auto get_edges = [&](
const Range &ents) {
663 CHKERR moab.get_adjacencies(ents, 1,
false, edges,
664 moab::Interface::UNION);
668 auto get_face_adj = [&, get_edges](
const Range &faces) {
670 CHKERR moab.get_adjacencies(subtract(get_edges(faces), edges), 2,
false,
671 adj_faces, moab::Interface::UNION);
672 return intersect(adj_faces, tris);
675 auto get_patch = [&, get_face_adj](
const EntityHandle face) {
677 patch_ents.insert(face);
680 nb0 = patch_ents.size();
681 patch_ents.merge(get_face_adj(patch_ents));
682 }
while (nb0 != patch_ents.size());
686 auto get_patches = [&]() {
687 std::vector<Range> patches;
688 while (!tris.empty()) {
689 patches.push_back(get_patch(tris[0]));
690 tris = subtract(tris, patches.back());
695 Tag th0, th1, th2, th3;
698 auto patches = get_patches();
700 for (
auto patch : patches) {
701 std::vector<int> tags_vals(patch.size(), pp);
702 CHKERR moab.tag_set_data(th2, patch, &*tags_vals.begin());
710 moab::Interface &moab,
Range edges,
Range tris,
711 bool number_pathes =
true,
712 boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
713 surface_orientation =
nullptr,
716 Tag th0, th1, th2, th3;
721 for (
auto edge : edges) {
723 CHKERR moab.get_adjacencies(&edge, 1, 2,
false, adj_faces);
724 adj_faces = intersect(adj_faces, tris);
725 if (adj_faces.size() != 2) {
727 "Should be 2 faces adjacent to edge but is %ld",
730 VectorDouble3
v[2] = {VectorDouble3(3), VectorDouble3(3)};
731 auto get_tensor_from_vec = [](VectorDouble3 &
v) {
733 &
v[0], &
v[1], &
v[2]);
740 std::array<double, 3> areas;
741 auto calculate_normals = [&, get_tensor_from_vec]() {
743 for (
auto face : adj_faces) {
744 double &x = (
v[ff])[0];
745 double &y = (
v[ff])[1];
746 double &z = (
v[ff])[2];
747 moab::Util::normal(&moab, face, x, y, z);
748 auto t_n = get_tensor_from_vec(
v[ff]);
749 areas[ff] = sqrt(t_n(
i) * t_n(
i));
753 CHKERR moab.tag_get_data(th3, &face, 1, &orientation);
754 if (orientation == -1) {
757 if (surface_orientation && m_field_ptr) {
758 CHKERR surface_orientation->getElementOrientation(*m_field_ptr,
760 int eo = surface_orientation->elementOrientation;
770 auto get_patch_number = [&]() {
771 std::vector<int> p = {0, 0};
772 CHKERR moab.tag_get_data(th2, adj_faces, &*p.begin());
776 std::array<EntityHandle, 2> faces_handles = {adj_faces[0],
778 auto order_normals = [&, get_patch_number]() {
779 auto p = get_patch_number();
782 faces_handles[0] = adj_faces[1];
783 faces_handles[1] = adj_faces[0];
789 auto t_n0 = get_tensor_from_vec(
v[0]);
790 auto t_n1 = get_tensor_from_vec(
v[1]);
793 auto get_base_for_coplanar_faces = [&]() {
796 std::vector<EntityHandle> face_conn;
797 CHKERR moab.get_connectivity(&faces_handles[1], 1, face_conn,
true);
798 std::vector<EntityHandle> edge_conn;
799 CHKERR moab.get_connectivity(&edge, 1, edge_conn,
true);
800 std::sort(face_conn.begin(), face_conn.end());
801 std::sort(edge_conn.begin(), edge_conn.end());
804 if (face_conn[
n] != edge_conn[
n])
806 VectorDouble3 coords_face(3);
807 CHKERR moab.get_coords(&face_conn[
n], 1,
808 &*coords_face.data().begin());
809 VectorDouble6 coords_edge(6);
810 CHKERR moab.get_coords(&edge_conn[0], 2,
811 &*coords_edge.data().begin());
813 coords_edge[0], coords_edge[1], coords_edge[2]);
815 coords_edge[3], coords_edge[4], coords_edge[5]);
816 auto t_face_n = get_tensor_from_vec(coords_face);
817 t_face_n(
i) -= t_edge_n0(
i);
819 t_delta(
i) = t_edge_n1(
i) - t_edge_n0(
i);
820 t_delta(
i) /= sqrt(t_delta(
i) * t_delta(
i));
822 if (t_n1(
i) * t_face_n(
i) < 0)
828 auto get_base_for_not_planar_faces = [&]() {
831 t_n1(
i) /= sqrt(t_n1(
i) * t_n1(
i));
836 constexpr double tol = 1e-6;
837 if ((t_cross(
i) * t_cross(
i)) <
tol)
838 CHKERR get_base_for_coplanar_faces();
840 CHKERR get_base_for_not_planar_faces();
842 VectorDouble3 &v0 =
v[0];
843 VectorDouble3 &v1 =
v[1];
844 CHKERR moab.tag_set_data(th0, &edge, 1, &v0[0]);
845 CHKERR moab.tag_set_data(th1, &edge, 1, &v1[0]);
850 static MoFEMErrorCode
saveEdges(moab::Interface &moab, std::string name,
856 CHKERR moab.create_meshset(MESHSET_SET, meshset);
857 CHKERR moab.add_entities(meshset, edges);
858 if (faces !=
nullptr) {
859 CHKERR moab.add_entities(meshset, *faces);
861 CHKERR moab.write_file(name.c_str(),
"VTK",
"", &meshset, 1, ths, 4);
862 CHKERR moab.delete_entities(&meshset, 1);
884 if (
B != PETSC_NULLPTR) {
888 if (
F != PETSC_NULLPTR) {
907 PetscOptionsBegin(PETSC_COMM_WORLD,
"",
908 "Get edge sliding constrains element scaling",
910 CHKERR PetscOptionsScalar(
"-edge_sliding_alpha",
"scaling parameter",
"",
921 CHKERRABORT(PETSC_COMM_WORLD,
ierr);
936 boost::shared_ptr<VectorDouble> &active_variables_ptr,
937 boost::shared_ptr<VectorDouble> &results_ptr,
938 boost::shared_ptr<MatrixDouble> &jacobian_ptr,
939 bool evaluate_jacobian,
const double &alpha)
946 MoFEMErrorCode
doWork(
int side, EntityType type,
947 EntitiesFieldData::EntData &data) {
949 if (type != MBVERTEX)
957 Tag th0, th1, th2, th3;
967 VectorInt &indices = data.getIndices();
971 ublas::vector<adouble> lambda_dofs(4);
972 for (
int dd = 0; dd != 4; ++dd) {
973 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
975 ublas::vector<adouble> position_dofs(6);
976 for (
int dd = 0; dd != 6; ++dd) {
977 position_dofs[dd] <<= (*activeVariablesPtr)[4 + dd];
981 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
983 &position_dofs[3], &position_dofs[4], &position_dofs[5]);
986 t_tangent(
i) = t_node1(
i) - t_node0(
i);
987 adouble l = sqrt(t_tangent(
i) * t_tangent(
i));
991 t_dot0 = t_edge_base0(
i) * t_tangent(
i);
992 t_dot1 = t_edge_base1(
i) * t_tangent(
i);
995 t_base0(
i) = t_edge_base0(
i) - t_dot0 * t_tangent(
i);
996 t_base1(
i) = t_edge_base1(
i) - t_dot1 * t_tangent(
i);
997 t_base0(
i) /= sqrt(t_base0(
i) * t_base0(
i));
998 t_base1(
i) /= sqrt(t_base1(
i) * t_base1(
i));
1000 auto t_base_fun1 = data.getFTensor0N();
1001 auto t_base_fun2 = data.getFTensor0N();
1007 ublas::vector<adouble> c_vec(4);
1008 ublas::vector<adouble> f_vec(6);
1013 int nb_gauss_pts = data.getN().size1();
1014 int nb_base_functions = data.getN().size2();
1015 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
1018 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
1020 &lambda_dofs[0], &lambda_dofs[1]);
1024 for (
int bb = 0; bb != nb_base_functions; ++bb) {
1025 t_position(
i) += t_base_fun1 * t_position_dofs(
i);
1026 t_lambda(
j) += t_base_fun1 * t_lambda_dof(
j);
1032 t_delta(
i) = t_position(
i) - t_coord_ref(
i);
1036 double w = t_w *
aLpha;
1041 &f_vec[0], &f_vec[1], &f_vec[2]);
1042 for (
int bb = 0; bb != nb_base_functions; ++bb) {
1043 if (indices[2 * bb] != -1) {
1044 val = w * t_base_fun2 *
l;
1045 t_c(N0) += val * dot0;
1046 t_c(N1) += val * dot1;
1047 val1 = val * t_lambda(N0);
1048 val2 = val * t_lambda(N1);
1049 t_f(
i) += val1 * t_base0(
i) + val2 * t_base1(
i);
1060 for (
int rr = 0; rr != 4; ++rr) {
1061 c_vec[rr] >>= (*resultsPtr)[rr];
1063 for (
int rr = 0; rr != 6; ++rr) {
1064 f_vec(rr) >>= (*resultsPtr)[4 + rr];
1070 double *jac_ptr[4 + 6];
1071 for (
int rr = 0; rr != 4 + 6; ++rr) {
1072 jac_ptr[rr] = &(*jacobianPtr)(rr, 0);
1079 "ADOL-C function evaluation with error");
1088 const std::string lagrange_multipliers_field_name,
1089 const std::string material_field_name) {
1094 material_field_name);
1099 const std::string lagrange_multipliers_field_name,
1100 const std::string material_field_name,
1101 const double *alpha =
nullptr) {
1108 boost::shared_ptr<VectorDouble> active_variables_ptr(
1109 new VectorDouble(4 + 6));
1110 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(4 + 6));
1111 boost::shared_ptr<MatrixDouble> jacobian_ptr(
1112 new MatrixDouble(4 + 6, 4 + 6));
1116 lagrange_multipliers_field_name, active_variables_ptr));
1118 material_field_name, active_variables_ptr));
1120 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1121 jacobian_ptr,
false,
aLpha));
1130 lagrange_multipliers_field_name, active_variables_ptr));
1132 material_field_name, active_variables_ptr));
1134 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1135 jacobian_ptr,
true,
aLpha));
1137 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
1139 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
1141 material_field_name, material_field_name, jacobian_ptr));
1148 const std::string lagrange_multipliers_field_name,
1149 const std::string material_field_name) {
1155 material_field_name);
1161 const std::string lagrange_multipliers_field_name,
1162 const std::string material_field_name) {
1165 boost::shared_ptr<VectorDouble> active_variables_ptr(
1166 new VectorDouble(4 + 6));
1167 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(4 + 6));
1168 boost::shared_ptr<MatrixDouble> jacobian_ptr(
1169 new MatrixDouble(4 + 6, 4 + 6));
1174 lagrange_multipliers_field_name, active_variables_ptr));
1176 material_field_name, active_variables_ptr));
1178 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1179 jacobian_ptr,
true,
aLpha));
1181 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static PetscErrorCode ierr
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_OPERATION_UNSUCCESSFUL
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
auto cross(const Tensor1_Expr< A, T, 3, i > &a, const Tensor1_Expr< B, U, 3, j > &b, const Index< k, 3 > &)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
constexpr auto field_name
static MoFEMErrorCode numberSurfaces(moab::Interface &moab, Range edges, Range tris)
static MoFEMErrorCode createTag(moab::Interface &moab, Tag &th0, Tag &th1, Tag &th2, Tag &th3)
static MoFEMErrorCode setTags(moab::Interface &moab, Range edges, Range tris, bool number_pathes=true, boost::shared_ptr< SurfaceSlidingConstrains::DriverElementOrientation > surface_orientation=nullptr, MoFEM::Interface *m_field_ptr=nullptr)
static MoFEMErrorCode saveEdges(moab::Interface &moab, std::string name, Range edges, Range *faces=nullptr)
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
MyEdgeFE(MoFEM::Interface &m_field)
boost::shared_ptr< VectorDouble > activeVariablesPtr
boost::shared_ptr< VectorDouble > resultsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpJacobian(int tag, const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr, boost::shared_ptr< VectorDouble > &results_ptr, boost::shared_ptr< MatrixDouble > &jacobian_ptr, bool evaluate_jacobian, const double &alpha)
boost::shared_ptr< MatrixDouble > jacobianPtr
MoFEMErrorCode setOperators(int tag, Range edges, Range faces, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
EdgeSlidingConstrains(MoFEM::Interface &m_field)
boost::shared_ptr< MyEdgeFE > feRhsPtr
MoFEM::Interface & mField
MyEdgeFE & getLoopFeRhs()
MyEdgeFE & getLoopFeLhs()
MoFEMErrorCode setOperators(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name, const double *alpha=nullptr)
MoFEMErrorCode getOptions()
boost::shared_ptr< MyEdgeFE > feLhsPtr
MoFEMErrorCode setOperatorsConstrainOnly(int tag, Range edges, Range faces, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
MoFEMErrorCode setOperatorsConstrainOnly(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
boost::shared_ptr< MatrixDouble > jacobianPtr
OpAssembleLhs(const std::string field_name_row, const std::string field_name_col, boost::shared_ptr< MatrixDouble > &jacobian_ptr)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpAssembleRhs(const std::string field_name, boost::shared_ptr< VectorDouble > &results_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< VectorDouble > resultsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpGetActiveDofsLambda(const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr)
boost::shared_ptr< VectorDouble > activeVariablesPtr
OpGetActiveDofsPositions(const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< VectorDouble > activeVariablesPtr
Implementation of surface sliding constrains.
virtual ~GenericSliding()=default
virtual moab::Interface & get_moab()=0
bool sYmm
If true assume that matrix is symmetric structure.
Deprecated interface functions.
default operator for EDGE element
const EdgeElementForcesAndSourcesCore * getEdgeFE()
get pointer to this finite element
default operator for TRI element
FaceElementForcesAndSourcesCore * getFaceFE()
return pointer to Generic Triangle Finite Element object
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
std::string getFEName() const
Get name of the element.
@ OPCOL
operator doWork function is executed on FE columns
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
virtual MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Vec & snes_f
Reference to residual vector.
Mat & snes_B
Reference to preconditioner of Jacobian matrix.
Class implemented by user to detect face orientation.
virtual MoFEMErrorCode getElementOrientation(MoFEM::Interface &m_field, const EntityHandle face)
virtual MoFEMErrorCode getElementOrientation(MoFEM::Interface &m_field, const FEMethod *fe_method_ptr)
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
MyTriangleFE(MoFEM::Interface &m_field)
OpJacobian(int tag, const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr, boost::shared_ptr< VectorDouble > &results_ptr, boost::shared_ptr< MatrixDouble > &jacobian_ptr, DriverElementOrientation &orientation, bool evaluate_jacobian, const double &alpha)
boost::shared_ptr< VectorDouble > resultsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
DriverElementOrientation & oRientation
boost::shared_ptr< VectorDouble > activeVariablesPtr
boost::shared_ptr< MatrixDouble > jacobianPtr
Shape preserving constrains, i.e. nodes sliding on body surface.
boost::shared_ptr< MyTriangleFE > feRhsPtr
virtual ~SurfaceSlidingConstrains()=default
DriverElementOrientation & crackFrontOrientation
boost::shared_ptr< MyTriangleFE > feLhsPtr
MyTriangleFE & getLoopFeLhs()
MoFEM::Interface & mField
SurfaceSlidingConstrains(MoFEM::Interface &m_field, DriverElementOrientation &orientation)
MoFEMErrorCode getOptions()
MoFEMErrorCode setOperators(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name, const double *alpha=nullptr)
MyTriangleFE & getLoopFeRhs()
MoFEMErrorCode setOperatorsConstrainOnly(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name)