5#ifndef __SURFACE_SLIDING_CONSTRAINS_HPP__
6#define __SURFACE_SLIDING_CONSTRAINS_HPP__
23 boost::shared_ptr<VectorDouble> &active_variables_ptr)
27 MoFEMErrorCode
doWork(
int side, EntityType type,
28 EntitiesFieldData::EntData &data) {
30 if (type == MBVERTEX) {
31 for (
unsigned int dd = 0; dd != data.getFieldData().size(); ++dd)
38 template <
int SizeLambda>
44 boost::shared_ptr<VectorDouble> &active_variables_ptr)
48 MoFEMErrorCode
doWork(
int side, EntityType type,
49 EntitiesFieldData::EntData &data) {
51 if (type == MBVERTEX) {
52 for (
unsigned int dd = 0; dd != data.getFieldData().size(); ++dd)
59 template <
int SizeLambda,
int SizePositions>
65 boost::shared_ptr<VectorDouble> &results_ptr)
70 MoFEMErrorCode
doWork(
int side, EntityType type,
71 EntitiesFieldData::EntData &data) {
75 VectorInt &indices = data.getIndices();
79 else if (indices.size() == SizeLambda)
81 else if (indices.size() == SizePositions)
85 "Element %s: Data inconsistency nb of indices %d",
88 CHKERR VecSetOption(
getSNESf(), VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
94 template <
int SizeLambda,
int SizePositions>
100 const std::string field_name_col,
101 boost::shared_ptr<MatrixDouble> &jacobian_ptr)
108 MoFEMErrorCode
doWork(
int row_side,
int col_side, EntityType row_type,
110 EntitiesFieldData::EntData &row_data,
111 EntitiesFieldData::EntData &col_data) {
113 if (row_type != MBVERTEX)
115 if (col_type != MBVERTEX)
117 VectorInt &row_indices = row_data.getIndices();
118 VectorInt &col_indices = col_data.getIndices();
119 if (row_indices.empty() || col_indices.empty())
123 if (row_indices.size() == SizeLambda)
125 else if (row_indices.size() == SizePositions)
126 shift_row = SizeLambda;
129 "Data inconsistency");
132 if (col_indices.size() == SizeLambda)
134 else if (col_indices.size() == SizePositions)
135 shift_col = SizeLambda;
138 "Data inconsistency");
140 MatrixDouble jac(row_indices.size(), col_indices.size());
141 for (
unsigned int rr = 0; rr != row_indices.size(); ++rr) {
142 for (
unsigned int cc = 0; cc != col_indices.size(); ++cc) {
143 jac(rr, cc) = (*jacobianPtr)(shift_row + rr, shift_col + cc);
144 if (jac(rr, cc) != jac(rr, cc))
146 "Jacobian assemble not a number jac(rr,cc) = %3.4f",
150 CHKERR MatSetValues(
getSNESB(), row_indices.size(), &row_indices[0],
151 col_indices.size(), &col_indices[0], &jac(0, 0),
336 virtual MoFEMErrorCode
368 if (
B != PETSC_NULLPTR) {
372 if (
F != PETSC_NULLPTR) {
392 PetscOptionsBegin(PETSC_COMM_WORLD,
"",
393 "Get surface sliding constrains element scaling",
395 CHKERR PetscOptionsScalar(
"-surface_sliding_alpha",
"scaling parameter",
"",
408 CHKERRABORT(PETSC_COMM_WORLD,
ierr);
426 boost::shared_ptr<VectorDouble> &active_variables_ptr,
427 boost::shared_ptr<VectorDouble> &results_ptr,
428 boost::shared_ptr<MatrixDouble> &jacobian_ptr,
438 MoFEMErrorCode
doWork(
int side, EntityType type,
439 EntitiesFieldData::EntData &data) {
441 if (type != MBVERTEX)
454 int nb_gauss_pts = data.getN().size1();
455 int nb_base_functions = data.getN().size2();
456 auto t_base1 = data.getFTensor0N();
457 auto t_base2 = data.getFTensor0N();
458 auto t_diff_base = data.getFTensor1DiffN<2>();
466 ublas::vector<adouble> c_vec(3);
467 ublas::vector<adouble> f_vec(9);
473 ublas::vector<adouble> lambda_dofs(3);
474 for (
int dd = 0; dd != 3; ++dd) {
475 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
477 ublas::vector<adouble> position_dofs(9);
478 for (
int dd = 0; dd != 9; ++dd) {
479 position_dofs[dd] <<= (*activeVariablesPtr)[3 + dd];
484 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
487 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
492 t_position_ksi(
i) = 0;
493 t_position_eta(
i) = 0;
496 for (
int bb = 0; bb != nb_base_functions; ++bb) {
497 t_position(
i) += t_base1 * t_position_dofs(
i);
498 t_position_ksi(
i) += t_diff_base(N0) * t_position_dofs(
i);
499 t_position_eta(
i) += t_diff_base(N1) * t_position_dofs(
i);
500 lambda += t_base1 * t_lambda_dof;
507 t_delta(
i) = t_position(
i) - t_coord_ref(
i);
512 &f_vec[0], &f_vec[1], &f_vec[2]);
514 double w = t_w * (0.5 *
aLpha * eo);
516 for (
int bb = 0; bb != nb_base_functions; ++bb) {
517 double alpha = w * t_base2;
518 t_c += alpha * t_normal(
i) * t_delta(
i);
519 t_f(
i) += alpha * (
lambda * t_normal(
i));
529 for (
int rr = 0; rr != 3; ++rr)
532 for (
int rr = 0; rr != 9; ++rr)
533 f_vec(rr) >>= (*resultsPtr)[3 + rr];
538 double *jac_ptr[3 + 9];
539 for (
int rr = 0; rr != 3 + 9; ++rr)
547 "ADOL-C function evaluation with error");
555 const std::string lagrange_multipliers_field_name,
556 const std::string material_field_name,
557 const double *alpha =
nullptr) {
560 if (alpha !=
nullptr) {
564 boost::shared_ptr<VectorDouble> active_variables_ptr(
565 new VectorDouble(3 + 9));
566 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(3 + 9));
567 boost::shared_ptr<MatrixDouble> jacobian_ptr(
568 new MatrixDouble(3 + 9, 3 + 9));
572 lagrange_multipliers_field_name, active_variables_ptr));
574 material_field_name, active_variables_ptr));
576 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
586 lagrange_multipliers_field_name, active_variables_ptr));
588 material_field_name, active_variables_ptr));
590 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
593 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
595 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
597 material_field_name, material_field_name, jacobian_ptr));
604 const std::string lagrange_multipliers_field_name,
605 const std::string material_field_name) {
608 boost::shared_ptr<VectorDouble> active_variables_ptr(
609 new VectorDouble(3 + 9));
610 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(3 + 9));
611 boost::shared_ptr<MatrixDouble> jacobian_ptr(
612 new MatrixDouble(3 + 9, 3 + 9));
617 lagrange_multipliers_field_name, active_variables_ptr));
619 material_field_name, active_variables_ptr));
621 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
624 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
637 double def_val[] = {0, 0, 0};
638 CHKERR moab.tag_get_handle(
"EDGE_BASE0", 3, MB_TYPE_DOUBLE, th0,
639 MB_TAG_CREAT | MB_TAG_SPARSE, def_val);
640 CHKERR moab.tag_get_handle(
"EDGE_BASE1", 3, MB_TYPE_DOUBLE, th1,
641 MB_TAG_CREAT | MB_TAG_SPARSE, def_val);
642 int def_numb_val[] = {-1};
643 CHKERR moab.tag_get_handle(
"PATCH_NUMBER", 1, MB_TYPE_INTEGER, th2,
644 MB_TAG_CREAT | MB_TAG_SPARSE, def_numb_val);
645 int def_orientation_val[] = {1};
646 CHKERR moab.tag_get_handle(
"PATCH_ORIENTATION", 1, MB_TYPE_INTEGER, th3,
647 MB_TAG_CREAT | MB_TAG_SPARSE,
648 def_orientation_val);
657 auto get_edges = [&](
const Range &ents) {
659 CHKERR moab.get_adjacencies(ents, 1,
false, edges,
660 moab::Interface::UNION);
664 auto get_face_adj = [&, get_edges](
const Range &faces) {
666 CHKERR moab.get_adjacencies(subtract(get_edges(faces), edges), 2,
false,
667 adj_faces, moab::Interface::UNION);
668 return intersect(adj_faces, tris);
671 auto get_patch = [&, get_face_adj](
const EntityHandle face) {
673 patch_ents.insert(face);
676 nb0 = patch_ents.size();
677 patch_ents.merge(get_face_adj(patch_ents));
678 }
while (nb0 != patch_ents.size());
682 auto get_patches = [&]() {
683 std::vector<Range> patches;
684 while (!tris.empty()) {
685 patches.push_back(get_patch(tris[0]));
686 tris = subtract(tris, patches.back());
691 Tag th0, th1, th2, th3;
694 auto patches = get_patches();
696 for (
auto patch : patches) {
697 std::vector<int> tags_vals(patch.size(), pp);
698 CHKERR moab.tag_set_data(th2, patch, &*tags_vals.begin());
706 moab::Interface &moab,
Range edges,
Range tris,
707 bool number_pathes =
true,
708 boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
709 surface_orientation =
nullptr,
712 Tag th0, th1, th2, th3;
717 for (
auto edge : edges) {
719 CHKERR moab.get_adjacencies(&edge, 1, 2,
false, adj_faces);
720 adj_faces = intersect(adj_faces, tris);
721 if (adj_faces.size() != 2) {
723 "Should be 2 faces adjacent to edge but is %ld",
726 VectorDouble3
v[2] = {VectorDouble3(3), VectorDouble3(3)};
727 auto get_tensor_from_vec = [](VectorDouble3 &
v) {
729 &
v[0], &
v[1], &
v[2]);
736 std::array<double, 3> areas;
737 auto calculate_normals = [&, get_tensor_from_vec]() {
739 for (
auto face : adj_faces) {
740 double &x = (
v[ff])[0];
741 double &y = (
v[ff])[1];
742 double &z = (
v[ff])[2];
743 moab::Util::normal(&moab, face, x, y, z);
744 auto t_n = get_tensor_from_vec(
v[ff]);
745 areas[ff] = sqrt(t_n(
i) * t_n(
i));
749 CHKERR moab.tag_get_data(th3, &face, 1, &orientation);
750 if (orientation == -1) {
753 if (surface_orientation && m_field_ptr) {
754 CHKERR surface_orientation->getElementOrientation(*m_field_ptr,
756 int eo = surface_orientation->elementOrientation;
766 auto get_patch_number = [&]() {
767 std::vector<int> p = {0, 0};
768 CHKERR moab.tag_get_data(th2, adj_faces, &*p.begin());
772 std::array<EntityHandle, 2> faces_handles = {adj_faces[0],
774 auto order_normals = [&, get_patch_number]() {
775 auto p = get_patch_number();
778 faces_handles[0] = adj_faces[1];
779 faces_handles[1] = adj_faces[0];
785 auto t_n0 = get_tensor_from_vec(
v[0]);
786 auto t_n1 = get_tensor_from_vec(
v[1]);
789 auto get_base_for_coplanar_faces = [&]() {
792 std::vector<EntityHandle> face_conn;
793 CHKERR moab.get_connectivity(&faces_handles[1], 1, face_conn,
true);
794 std::vector<EntityHandle> edge_conn;
795 CHKERR moab.get_connectivity(&edge, 1, edge_conn,
true);
796 std::sort(face_conn.begin(), face_conn.end());
797 std::sort(edge_conn.begin(), edge_conn.end());
800 if (face_conn[
n] != edge_conn[
n])
802 VectorDouble3 coords_face(3);
803 CHKERR moab.get_coords(&face_conn[
n], 1,
804 &*coords_face.data().begin());
805 VectorDouble6 coords_edge(6);
806 CHKERR moab.get_coords(&edge_conn[0], 2,
807 &*coords_edge.data().begin());
809 coords_edge[0], coords_edge[1], coords_edge[2]);
811 coords_edge[3], coords_edge[4], coords_edge[5]);
812 auto t_face_n = get_tensor_from_vec(coords_face);
813 t_face_n(
i) -= t_edge_n0(
i);
815 t_delta(
i) = t_edge_n1(
i) - t_edge_n0(
i);
816 t_delta(
i) /= sqrt(t_delta(
i) * t_delta(
i));
818 if (t_n1(
i) * t_face_n(
i) < 0)
824 auto get_base_for_not_planar_faces = [&]() {
827 t_n1(
i) /= sqrt(t_n1(
i) * t_n1(
i));
832 constexpr double tol = 1e-6;
833 if ((t_cross(
i) * t_cross(
i)) <
tol)
834 CHKERR get_base_for_coplanar_faces();
836 CHKERR get_base_for_not_planar_faces();
838 VectorDouble3 &v0 =
v[0];
839 VectorDouble3 &v1 =
v[1];
840 CHKERR moab.tag_set_data(th0, &edge, 1, &v0[0]);
841 CHKERR moab.tag_set_data(th1, &edge, 1, &v1[0]);
846 static MoFEMErrorCode
saveEdges(moab::Interface &moab, std::string name,
852 CHKERR moab.create_meshset(MESHSET_SET, meshset);
853 CHKERR moab.add_entities(meshset, edges);
854 if (faces !=
nullptr) {
855 CHKERR moab.add_entities(meshset, *faces);
857 CHKERR moab.write_file(name.c_str(),
"VTK",
"", &meshset, 1, ths, 4);
858 CHKERR moab.delete_entities(&meshset, 1);
880 if (
B != PETSC_NULLPTR) {
884 if (
F != PETSC_NULLPTR) {
903 PetscOptionsBegin(PETSC_COMM_WORLD,
"",
904 "Get edge sliding constrains element scaling",
"none");
905 CHKERR PetscOptionsScalar(
"-edge_sliding_alpha",
"scaling parameter",
"",
916 CHKERRABORT(PETSC_COMM_WORLD,
ierr);
931 boost::shared_ptr<VectorDouble> &active_variables_ptr,
932 boost::shared_ptr<VectorDouble> &results_ptr,
933 boost::shared_ptr<MatrixDouble> &jacobian_ptr,
934 bool evaluate_jacobian,
const double &alpha)
941 MoFEMErrorCode
doWork(
int side, EntityType type,
942 EntitiesFieldData::EntData &data) {
944 if (type != MBVERTEX)
952 Tag th0, th1, th2, th3;
962 VectorInt &indices = data.getIndices();
966 ublas::vector<adouble> lambda_dofs(4);
967 for (
int dd = 0; dd != 4; ++dd) {
968 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
970 ublas::vector<adouble> position_dofs(6);
971 for (
int dd = 0; dd != 6; ++dd) {
972 position_dofs[dd] <<= (*activeVariablesPtr)[4 + dd];
976 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
978 &position_dofs[3], &position_dofs[4], &position_dofs[5]);
981 t_tangent(
i) = t_node1(
i) - t_node0(
i);
982 adouble l = sqrt(t_tangent(
i) * t_tangent(
i));
986 t_dot0 = t_edge_base0(
i) * t_tangent(
i);
987 t_dot1 = t_edge_base1(
i) * t_tangent(
i);
990 t_base0(
i) = t_edge_base0(
i) - t_dot0 * t_tangent(
i);
991 t_base1(
i) = t_edge_base1(
i) - t_dot1 * t_tangent(
i);
992 t_base0(
i) /= sqrt(t_base0(
i) * t_base0(
i));
993 t_base1(
i) /= sqrt(t_base1(
i) * t_base1(
i));
995 auto t_base_fun1 = data.getFTensor0N();
996 auto t_base_fun2 = data.getFTensor0N();
1002 ublas::vector<adouble> c_vec(4);
1003 ublas::vector<adouble> f_vec(6);
1008 int nb_gauss_pts = data.getN().size1();
1009 int nb_base_functions = data.getN().size2();
1010 for (
int gg = 0; gg != nb_gauss_pts; ++gg) {
1013 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
1015 &lambda_dofs[0], &lambda_dofs[1]);
1019 for (
int bb = 0; bb != nb_base_functions; ++bb) {
1020 t_position(
i) += t_base_fun1 * t_position_dofs(
i);
1021 t_lambda(
j) += t_base_fun1 * t_lambda_dof(
j);
1027 t_delta(
i) = t_position(
i) - t_coord_ref(
i);
1031 double w = t_w *
aLpha;
1036 &f_vec[0], &f_vec[1], &f_vec[2]);
1037 for (
int bb = 0; bb != nb_base_functions; ++bb) {
1038 if (indices[2 * bb] != -1) {
1039 val = w * t_base_fun2 *
l;
1040 t_c(N0) += val * dot0;
1041 t_c(N1) += val * dot1;
1042 val1 = val * t_lambda(N0);
1043 val2 = val * t_lambda(N1);
1044 t_f(
i) += val1 * t_base0(
i) + val2 * t_base1(
i);
1055 for (
int rr = 0; rr != 4; ++rr) {
1056 c_vec[rr] >>= (*resultsPtr)[rr];
1058 for (
int rr = 0; rr != 6; ++rr) {
1059 f_vec(rr) >>= (*resultsPtr)[4 + rr];
1065 double *jac_ptr[4 + 6];
1066 for (
int rr = 0; rr != 4 + 6; ++rr) {
1067 jac_ptr[rr] = &(*jacobianPtr)(rr, 0);
1074 "ADOL-C function evaluation with error");
1083 const std::string lagrange_multipliers_field_name,
1084 const std::string material_field_name) {
1089 material_field_name);
1094 const std::string lagrange_multipliers_field_name,
1095 const std::string material_field_name,
1096 const double *alpha =
nullptr) {
1103 boost::shared_ptr<VectorDouble> active_variables_ptr(
1104 new VectorDouble(4 + 6));
1105 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(4 + 6));
1106 boost::shared_ptr<MatrixDouble> jacobian_ptr(
1107 new MatrixDouble(4 + 6, 4 + 6));
1111 lagrange_multipliers_field_name, active_variables_ptr));
1113 material_field_name, active_variables_ptr));
1115 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1116 jacobian_ptr,
false,
aLpha));
1125 lagrange_multipliers_field_name, active_variables_ptr));
1127 material_field_name, active_variables_ptr));
1129 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1130 jacobian_ptr,
true,
aLpha));
1132 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
1134 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
1136 material_field_name, material_field_name, jacobian_ptr));
1143 const std::string lagrange_multipliers_field_name,
1144 const std::string material_field_name) {
1150 material_field_name);
1156 const std::string lagrange_multipliers_field_name,
1157 const std::string material_field_name) {
1160 boost::shared_ptr<VectorDouble> active_variables_ptr(
1161 new VectorDouble(4 + 6));
1162 boost::shared_ptr<VectorDouble> results_ptr(
new VectorDouble(4 + 6));
1163 boost::shared_ptr<MatrixDouble> jacobian_ptr(
1164 new MatrixDouble(4 + 6, 4 + 6));
1169 lagrange_multipliers_field_name, active_variables_ptr));
1171 material_field_name, active_variables_ptr));
1173 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1174 jacobian_ptr,
true,
aLpha));
1176 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 > &)
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)