38#ifdef ENABLE_PYTHON_BINDING
39#include <boost/python.hpp>
40#include <boost/python/def.hpp>
41#include <boost/python/numpy.hpp>
42namespace bp = boost::python;
43namespace np = boost::python::numpy;
46#include <ContactOps.hpp>
52 double circumcenter[3],
double *xi,
double *
eta);
58 using ContactOps::SDFPython::SDFPython;
62 boost::shared_ptr<ContactSDFPython> sdf_python_ptr;
64#ifdef ENABLE_PYTHON_BINDING
66 auto file_exists = [](std::string myfile) {
67 std::ifstream file(myfile.c_str());
74 char sdf_file_name[255] =
"sdf.py";
76 sdf_file_name, 255, PETSC_NULLPTR);
78 if (file_exists(sdf_file_name)) {
79 MOFEM_LOG(
"EP", Sev::inform) << sdf_file_name <<
" file found";
80 sdf_python_ptr = boost::make_shared<ContactSDFPython>();
81 CHKERR sdf_python_ptr->sdfInit(sdf_file_name);
82 ContactOps::sdfPythonWeakPtr = sdf_python_ptr;
83 MOFEM_LOG(
"EP", Sev::inform) <<
"SdfPython initialized";
85 MOFEM_LOG(
"EP", Sev::warning) << sdf_file_name <<
" file NOT found";
90 return sdf_python_ptr;
100 boost::shared_ptr<moab::Core> core_mesh_ptr,
int max_order,
101 std::map<int, Range> &&body_map);
133 using MapFaceData = std::map<EntityHandle, std::vector<FaceData>>;
137 auto it = map_face_data.find(fe_ent);
138 if (it == map_face_data.end()) {
139 return (std::vector<FaceData> *)
nullptr;
141 return &(it->second);
145 std::vector<FaceData> *vec_ptr) {
147 if (it != vec_ptr->end()) {
148 if (it->gaussPtNb == gg) {
149 face_data_ptr = &(*it);
153 return face_data_ptr;
178 for (
auto &m_sdf : sdf_map_range) {
179 if (m_sdf.second.find(fe_ent) != m_sdf.second.end())
185template <
typename OP_PTR>
189 auto nb_gauss_pts = op_ptr->getGaussPts().size2();
191 auto ts_time = op_ptr->getTStime();
192 auto ts_time_step = op_ptr->getTStimeStep();
199 op_ptr->getFTensor1CoordsAtGaussPts(),
200 getFTensor1FromMat<3>(contact_disp), nb_gauss_pts);
202 op_ptr->getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
204 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
207 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
211 if (v_sdf.size() != nb_gauss_pts)
213 "Wrong number of integration pts");
214 if (m_grad_sdf.size2() != nb_gauss_pts)
216 "Wrong number of integration pts");
217 if (m_grad_sdf.size1() != 3)
223 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
225 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
229 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
245template <
typename T1>
248 std::array<double, 3> &unit_ray, std::array<double, 3> &point,
250 std::array<double, 9> &elem_point_nodes,
251 std::array<double, 9> &elem_traction_nodes,
256 auto t_unit_ray = getFTensor1FromPtr<3>(unit_ray.data());
257 auto t_point = getFTensor1FromPtr<3>(point.data());
259 auto get_normal = [](
auto &ele_coords) {
265 auto t_normal = get_normal(elem_point_nodes);
266 t_normal(
i) /= t_normal.l2();
268 auto sn = t_normal(
i) * t_point(
i);
269 auto nm = t_normal(
i) * t_spatial_coords(
i);
270 auto nr = t_normal(
i) * t_unit_ray(
i);
272 auto gamma = (sn - nm) / nr;
275 t_point_current(
i) = t_spatial_coords(
i) + gamma * t_unit_ray(
i);
277 auto get_local_point_shape_functions = [&](
auto &&t_elem_coords,
279 std::array<T1, 2> loc_coords;
282 &t_elem_coords(0, 0), &t_point(0), 1, loc_coords.data()),
285 N_MBTRI1(loc_coords[0], loc_coords[1]),
286 N_MBTRI2(loc_coords[0], loc_coords[1])};
289 auto eval_position = [&](
auto &&t_field,
auto &t_shape_fun) {
293 t_point_out(
i) = t_shape_fun(
j) * t_field(
j,
i);
297 auto t_shape_fun = get_local_point_shape_functions(
298 getFTensor2FromPtr<3, 3>(elem_point_nodes.data()), t_point_current);
299 auto t_slave_point_updated = eval_position(
300 getFTensor2FromPtr<3, 3>(elem_point_nodes.data()), t_shape_fun);
301 auto t_traction_updated = eval_position(
302 getFTensor2FromPtr<3, 3>(elem_traction_nodes.data()), t_shape_fun);
304 return std::make_tuple(t_slave_point_updated, t_traction_updated, t_normal);
307template <
typename T1>
320template <
typename T1>
336template <
typename T1>
340 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
342 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
348 t_normal(
i) = t_master_normal(
i) - t_slave_normal(
i);
351 auto gap = t_normal(
i) * (t_slave_point_current(
i) - t_spatial_coords(
i));
352 auto tn_master = t_master_traction_current(
i) * t_normal(
i);
353 auto tn_slave = t_slave_traction_current(
i) * t_normal(
i);
354 auto tn = std::max(-tn_master, tn_slave);
356 return std::make_tuple(gap, tn_master, tn_slave,
358 t_master_traction_current, t_slave_traction_current);
365template <
typename T1,
typename T2,
typename T3>
378 t_u(
i) = t_spatial_coords(
i) - t_coords(
i);
385 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
387 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
392 t_normal(
i) = t_master_normal(
i) - t_slave_normal(
i);
397 t_P(
i,
j) = t_normal(
i) * t_normal(
j);
399 t_Q(
i,
j) = kronecker_delta(
i,
j) - t_P(
i,
j);
401 constexpr double beta = 0.5;
408 auto f_min_gap = [
zeta](
auto d,
auto s) {
409 return 0.5 * (d + s - std::sqrt((d - s) * (d - s) +
zeta));
412 auto f_diff_min_gap = [
zeta](
auto d,
auto s) {
413 return 0.5 * (1 - (d - s) / std::sqrt((d - s) * (d - s) +
zeta));
417 auto f_barrier = [alpha1, alpha2, f_min_gap, f_diff_min_gap](
auto g,
421 0.5 * (tn + f_min_gap(d, tn) + f_diff_min_gap(d, tn) * (d - tn));
422 auto b2 = alpha2 * f_min_gap(
g, 0) *
g;
427 t_gap_vec(
i) = beta * t_spatial_coords(
i) +
428 (1 - beta) * t_slave_point_current(
i) - t_spatial_coords(
i);
431 -beta * t_master_traction(
i) + (beta - 1) * t_slave_traction_current(
i);
433 auto t_gap = t_normal(
i) * t_gap_vec(
i);
434 auto t_tn = t_normal(
i) * t_traction_vec(
i);
435 auto barrier = f_barrier(t_gap, t_tn);
439 t_traction_bar(
i) = t_normal(
i) * barrier;
443 t_Q(
i,
j) * t_master_traction(
j) +
445 t_P(
i,
j) * (t_master_traction(
j) - t_traction_bar(
j));
448 auto is_nan_or_inf = [](
double value) ->
bool {
449 return std::isnan(value) || std::isinf(value);
452 double v = std::complex<double>(t_rhs(
i) * t_rhs(
i)).real();
453 if (is_nan_or_inf(
v)) {
455 MOFEM_LOG(
"SELF", Sev::error) <<
"t_rhs " << t_rhs;
471 const std::string row_field_name,
472 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
473 boost::shared_ptr<ContactTree> contact_tree_ptr,
474 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr =
nullptr);
485 const std::string row_field_name,
486 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
487 boost::shared_ptr<ContactTree> contact_tree_ptr,
488 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
491 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
492 sdfMapRangePtr(sdf_map_range_ptr) {
501 "get alpha contact failed");
503 PETSC_NULLPTR,
"",
"-alpha_contact_quadratic",
505 "get alpha contact failed");
509 "get alpha contact failed");
529 const size_t nb_gauss_pts = getGaussPts().size2();
534 "Wrong number of integration pts %ld != %ld",
542 auto t_w = getFTensor0IntegrationWeight();
543 auto t_coords = getFTensor1CoordsAtGaussPts();
544 auto t_disp = getFTensor1FromMat<3>(
commonDataPtr->contactDisp);
545 auto t_traction = getFTensor1FromMat<3>(
commonDataPtr->contactTraction);
550 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
555 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
556 auto t_normalize_normal = getFTensor1FromMat<3>(m_normals_at_pts);
563 ++t_normalize_normal;
568 auto face_data_vec_ptr =
570 auto face_gauss_pts_it = face_data_vec_ptr->begin();
572 auto nb_base_functions = data.
getN().size2();
574 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
577 auto face_data_ptr =
contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
580 auto check_face_contact = [&]() {
590#ifdef ENABLE_PYTHON_BINDING
592 if (ContactOps::sdfPythonWeakPtr.lock()) {
593 auto tn = t_traction(
i) * t_grad_sdf_v(
i);
597 constexpr double c = 0;
600 if (!
c && check_face_contact()) {
602 t_spatial_coords(
i) = t_coords(
i) + t_disp(
i);
603 auto t_rhs_tmp =
multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
605 t_rhs(
i) = t_rhs_tmp(
i);
609#ifdef ENABLE_PYTHON_BINDING
612 if (ContactOps::sdfPythonWeakPtr.lock()) {
614 t_cP(
i,
j) = (
c * t_grad_sdf_v(
i)) * t_grad_sdf_v(
j);
615 t_cQ(
i,
j) = kronecker_delta(
i,
j) - t_cP(
i,
j);
616 t_rhs(
i) = t_cQ(
i,
j) * t_traction(
j) +
617 (
c * inv_cn * t_sdf_v) * t_grad_sdf_v(
i);
619 t_rhs(
i) = t_traction(
i);
622 t_rhs(
i) = t_traction(
i);
626 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
627 const double alpha = t_w * getMeasure();
630 for (; bb != nbRows / 3; ++bb) {
631 const double beta = alpha * t_base;
632 t_nf(
i) -= beta * t_rhs(
i);
636 for (; bb < nb_base_functions; ++bb)
647template <AssemblyType A>
655 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
656 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
657 boost::shared_ptr<ContactTree> contact_tree_ptr);
666template <AssemblyType A>
669 boost::shared_ptr<std::vector<BrokenBaseSideData>>
670 broken_base_side_data,
671 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
672 boost::shared_ptr<ContactTree> contact_tree_ptr)
673 :
OP(broken_base_side_data), commonDataPtr(common_data_ptr),
674 contactTreePtr(contact_tree_ptr) {}
676template <AssemblyType A>
686 const size_t nb_gauss_pts = OP::getGaussPts().size2();
689 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
691 "Wrong number of integration pts %ld != %ld",
692 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
699 auto t_w = OP::getFTensor0IntegrationWeight();
700 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
701 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
703 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
704 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
714 auto face_data_vec_ptr =
715 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
716 auto face_gauss_pts_it = face_data_vec_ptr->begin();
718 auto nb_base_functions = data.
getN().size2() / 3;
720 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
722 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
723 const double alpha = t_w / 2.;
726 for (; bb != OP::nbRows / 3; ++bb) {
727 const double beta = alpha * t_base(
i) * t_material_normal(
i);
728 t_nf(
i) += beta * t_disp(
i);
732 for (; bb < nb_base_functions; ++bb)
744 const std::string row_field_name,
const std::string col_field_name,
745 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
746 boost::shared_ptr<ContactTree> contact_tree_ptr,
747 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr =
nullptr);
759 const std::string row_field_name,
const std::string col_field_name,
760 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
761 boost::shared_ptr<ContactTree> contact_tree_ptr,
762 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
765 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
766 sdfMapRangePtr(sdf_map_range_ptr) {
785 auto &locMat = AssemblyBoundaryEleOp::locMat;
786 locMat.resize(nb_rows, nb_cols,
false);
789 if (nb_cols && nb_rows) {
791 auto nb_gauss_pts = getGaussPts().size2();
792 auto t_w = getFTensor0IntegrationWeight();
793 auto t_coords = getFTensor1CoordsAtGaussPts();
794 auto t_disp = getFTensor1FromMat<3>(
commonDataPtr->contactDisp);
795 auto t_traction = getFTensor1FromMat<3>(
commonDataPtr->contactTraction);
798 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
803 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
804 auto t_hess_sdf_v = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
805 auto t_normalized_normal = getFTensor1FromMat<3>(m_normals_at_pts);
815 ++t_normalized_normal;
818 auto face_data_vec_ptr =
820 auto face_gauss_pts_it = face_data_vec_ptr->begin();
823 auto nb_face_functions = row_data.
getN().size2() / 3;
826 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
828 auto face_data_ptr =
contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
831 auto check_face_contact = [&]() {
843#ifdef ENABLE_PYTHON_BINDING
845 if (ContactOps::sdfPythonWeakPtr.lock()) {
846 auto tn = t_traction(
i) * t_grad_sdf_v(
i);
850 constexpr double c = 0;
853 if (!
c && check_face_contact()) {
855 t_spatial_coords(
i) = t_coords(
i) + t_disp(
i);
856 constexpr double eps = std::numeric_limits<float>::epsilon();
857 for (
auto ii = 0; ii < 3; ++ii) {
859 t_spatial_coords(0), t_spatial_coords(1), t_spatial_coords(2)};
860 t_spatial_coords_cx(ii) +=
eps * 1
i;
864 for (
int jj = 0; jj != 3; ++jj) {
865 auto v = t_rhs_tmp(jj).imag();
866 t_res_dU(jj, ii) =
v /
eps;
872#ifdef ENABLE_PYTHON_BINDING
874 if (ContactOps::sdfPythonWeakPtr.lock()) {
878 (-
c) * (t_hess_sdf_v(
i,
j) * t_grad_sdf_v(
k) * t_traction(
k) +
879 t_grad_sdf_v(
i) * t_hess_sdf_v(
k,
j) * t_traction(
k))
881 + (
c * inv_cn) * (t_sdf_v * t_hess_sdf_v(
i,
j) +
883 t_grad_sdf_v(
j) * t_grad_sdf_v(
i));
892 auto alpha = t_w * getMeasure();
895 for (; rr != nb_rows / 3; ++rr) {
897 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
900 for (
size_t cc = 0; cc != nb_cols / 3; ++cc) {
901 auto beta = alpha * t_row_base * t_col_base;
902 t_mat(
i,
j) -= beta * t_res_dU(
i,
j);
909 for (; rr < nb_face_functions; ++rr)
921template <AssemblyType A>
929 std::string row_field_name,
930 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
931 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
932 boost::shared_ptr<ContactTree> contact_tree_ptr,
933 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr =
nullptr);
944template <AssemblyType A>
947 std::string row_field_name,
948 boost::shared_ptr<std::vector<BrokenBaseSideData>>
949 broken_base_side_data,
950 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
951 boost::shared_ptr<ContactTree> contact_tree_ptr,
952 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
953 :
OP(row_field_name, broken_base_side_data, false, false),
954 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
955 sdfMapRangePtr(sdf_map_range_ptr) {
959template <AssemblyType A>
975 auto &locMat = AssemblyBoundaryEleOp::locMat;
976 locMat.resize(nb_rows, nb_cols,
false);
979 if (nb_cols && nb_rows) {
981 const size_t nb_gauss_pts = OP::getGaussPts().size2();
983 auto t_w = OP::getFTensor0IntegrationWeight();
984 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
985 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
986 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
987 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
990 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
991 getSdf(
this, commonDataPtr->contactDisp,
992 checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr),
false);
995 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
1002 ++t_material_normal;
1007 auto face_data_vec_ptr =
1008 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1009 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1012 auto nb_face_functions = row_data.
getN().size2();
1016 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1018 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
1021 auto check_face_contact = [&]() {
1022 if (
checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr) != -1)
1025 if (face_data_ptr) {
1033#ifdef ENABLE_PYTHON_BINDING
1035 if (ContactOps::sdfPythonWeakPtr.lock()) {
1036 auto tn = t_traction(
i) * t_grad_sdf_v(
i);
1040 constexpr double c = 0;
1043 if (!
c && check_face_contact()) {
1045 t_spatial_coords(
i) = t_coords(
i) + t_disp(
i);
1046 constexpr double eps = std::numeric_limits<float>::epsilon();
1047 for (
auto ii = 0; ii != 3; ++ii) {
1049 t_traction(0), t_traction(1), t_traction(2)};
1050 t_traction_cx(ii) +=
eps * 1
i;
1054 for (
int jj = 0; jj != 3; ++jj) {
1055 auto v = t_rhs_tmp(jj).imag();
1056 t_res_dP(jj, ii) =
v /
eps;
1061#ifdef ENABLE_PYTHON_BINDING
1062 if (ContactOps::sdfPythonWeakPtr.lock()) {
1064 t_cP(
i,
j) = (
c * t_grad_sdf_v(
i)) * t_grad_sdf_v(
j);
1065 t_cQ(
i,
j) = kronecker_delta(
i,
j) - t_cP(
i,
j);
1066 t_res_dP(
i,
j) = t_cQ(
i,
j);
1075 const double alpha = t_w / 2.;
1077 for (; rr != nb_rows / 3; ++rr) {
1079 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1082 for (
size_t cc = 0; cc != nb_cols / 3; ++cc) {
1083 auto col_base = t_col_base(
i) * t_material_normal(
i);
1084 const double beta = alpha * t_row_base * col_base;
1085 t_mat(
i,
j) -= beta * t_res_dP(
i,
j);
1092 for (; rr < nb_face_functions; ++rr)
1102template <AssemblyType A, IntegrationType I>
1105template <AssemblyType A>
1113 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1114 std::string col_field_name,
1115 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1116 boost::shared_ptr<ContactTree> contact_tree_ptr);
1126template <AssemblyType A>
1129 boost::shared_ptr<std::vector<BrokenBaseSideData>>
1130 broken_base_side_data,
1131 std::string col_field_name,
1132 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1133 boost::shared_ptr<ContactTree> contact_tree_ptr)
1134 :
OP(col_field_name, broken_base_side_data, true, true),
1135 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr) {
1139template <AssemblyType A>
1158 auto &locMat = AssemblyBoundaryEleOp::locMat;
1159 locMat.resize(nb_rows, nb_cols,
false);
1162 if (nb_cols && nb_rows) {
1164 const size_t nb_gauss_pts = OP::getGaussPts().size2();
1166 auto t_w = OP::getFTensor0IntegrationWeight();
1167 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
1168 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
1169 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
1170 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
1177 ++t_material_normal;
1182 auto face_data_vec_ptr =
1183 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1184 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1187 auto nb_face_functions = row_data.
getN().size2() / 3;
1188 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1190 const auto alpha = t_w / 2.;
1193 for (; rr != nb_rows / 3; ++rr) {
1195 auto row_base = alpha * (t_row_base(
i) * t_material_normal(
i));
1197 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1200 for (
size_t cc = 0; cc != nb_cols / 3; ++cc) {
1201 const auto beta = row_base * t_col_base;
1209 for (; rr < nb_face_functions; ++rr)
1216 locMat = trans(locMat);
1222 boost::shared_ptr<moab::Core> core_mesh_ptr,
1223 int max_order, std::map<int, Range> &&body_map)
1224 :
Base(m_field, core_mesh_ptr,
"contact"), maxOrder(max_order),
1227 auto ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
1228 ref_ele_ptr->hoNodes =
1233 "Error when generating reference element");
1235 MOFEM_LOG(
"EP", Sev::inform) <<
"Contact hoNodes " << ref_ele_ptr->hoNodes
1239 <<
"Contact maxOrder " << ref_ele_ptr->defMaxLevel;
1243 int def_ele_id = -1;
1245 MB_TAG_CREAT | MB_TAG_DENSE,
1248 thBodyId, MB_TAG_CREAT | MB_TAG_DENSE,
1268 std::array<double, 3> def_small_x{0., 0., 0.};
1270 MB_TAG_CREAT | MB_TAG_DENSE,
1271 def_small_x.data());
1273 MB_TAG_CREAT | MB_TAG_DENSE,
1274 def_small_x.data());
1276 std::array<double, 3> def_tractions{0., 0., 0.};
1278 "TRACTION", 3, MB_TYPE_DOUBLE,
thTraction, MB_TAG_CREAT | MB_TAG_DENSE,
1279 &*def_tractions.begin());
1294 PetscBarrier(
nullptr);
1296 ParallelComm *pcomm_post_proc_mesh =
1298 if (pcomm_post_proc_mesh ==
nullptr)
1303 auto brodacts = [&](
auto &brodacts_ents) {
1307 pcomm_post_proc_mesh->broadcast_entities(
1314 Range brodacts_ents;
1316 CHKERR brodacts(brodacts_ents);
1335 treeSurfPtr = boost::shared_ptr<OrientedBoxTreeTool>(
1348 OpMoveNode(boost::shared_ptr<ContactTree> contact_tree_ptr,
1349 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1350 boost::shared_ptr<MatrixDouble> u_h1_ptr);
1361 boost::shared_ptr<ContactTree> contact_tree_ptr,
1362 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1363 boost::shared_ptr<MatrixDouble> u_h1_ptr)
1364 :
UOP(
NOSPACE,
UOP::OPSPACE), contactTreePtr(contact_tree_ptr),
1365 uH1Ptr(u_h1_ptr), commonDataPtr(common_data_ptr) {}
1371 auto get_body_id = [&](
auto fe_ent) {
1373 if (
m.second.find(fe_ent) !=
m.second.end()) {
1383 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1385 auto body_id = get_body_id(fe_ent);
1389 post_proc_ents, &fe_id);
1391 post_proc_ents, &body_id);
1393 auto nb_gauss_pts = getGaussPts().size2();
1394 auto t_u_h1 = getFTensor1FromMat<3>(*
uH1Ptr);
1395 auto t_u_l2 = getFTensor1FromMat<3>(
commonDataPtr->contactDisp);
1396 auto t_coords = getFTensor1CoordsAtGaussPts();
1399 auto t_x_h1 = getFTensor1FromPtr<3>(&x_h1(0, 0));
1401 auto t_x_l2 = getFTensor1FromPtr<3>(&x_l2(0, 0));
1408 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
1409 t_x_h1(
i) = t_coords(
i) + t_u_h1(
i);
1410 t_x_l2(
i) = t_coords(
i) + t_u_l2(
i);
1419 CHKERR moab_post_proc_mesh.set_coords(
1420 &*map_gauss_pts.begin(), map_gauss_pts.size(), &*x_h1.data().begin());
1421 CHKERR moab_post_proc_mesh.tag_set_data(
1422 contactTreePtr->thSmallX, &*map_gauss_pts.begin(), map_gauss_pts.size(),
1423 &*x_h1.data().begin());
1424 CHKERR moab_post_proc_mesh.tag_set_data(
1425 contactTreePtr->thLargeX, &*map_gauss_pts.begin(), map_gauss_pts.size(),
1426 &*coords.data().begin());
1427 CHKERR moab_post_proc_mesh.tag_set_data(
1428 contactTreePtr->thTraction, &*map_gauss_pts.begin(), map_gauss_pts.size(),
1429 &*tractions.data().begin());
1438 OpTreeSearch(boost::shared_ptr<ContactTree> contact_tree_ptr,
1439 boost::shared_ptr<MatrixDouble> u_h1_ptr,
1440 boost::shared_ptr<MatrixDouble> traction_ptr,
Range r,
1442 moab::Interface *post_proc_mesh_ptr,
1443 std::vector<EntityHandle> *map_gauss_pts_ptr
1460 boost::shared_ptr<MatrixDouble> u_h1_ptr,
1461 boost::shared_ptr<MatrixDouble> traction_ptr,
1464 moab::Interface *post_proc_mesh_ptr,
1465 std::vector<EntityHandle> *map_gauss_pts_ptr
1468 :
UOP(
NOSPACE,
UOP::OPSPACE), contactTreePtr(contact_tree_ptr),
1469 uH1Ptr(u_h1_ptr), tractionPtr(traction_ptr),
1470 postProcMeshPtr(post_proc_mesh_ptr), mapGaussPtsPtr(map_gauss_pts_ptr),
1477 auto &m_field = getPtrFE()->mField;
1478 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1487 const auto nb_gauss_pts = getGaussPts().size2();
1489 auto t_disp_h1 = getFTensor1FromMat<3>(*
uH1Ptr);
1490 auto t_coords = getFTensor1CoordsAtGaussPts();
1491 auto t_traction = getFTensor1FromMat<3>(*
tractionPtr);
1499 auto get_ele_centre = [
i](
auto t_ele_coords) {
1501 t_ele_center(
i) = 0;
1502 for (
int nn = 0; nn != 3; nn++) {
1503 t_ele_center(
i) += t_ele_coords(
i);
1506 t_ele_center(
i) /= 3;
1507 return t_ele_center;
1510 auto get_ele_radius = [
i](
auto t_ele_center,
auto t_ele_coords) {
1512 t_n0(
i) = t_ele_center(
i) - t_ele_coords(
i);
1516 auto get_face_conn = [
this](
auto face) {
1520 face, conn, num_nodes,
true),
1522 if (num_nodes != 3) {
1528 auto get_face_coords = [
this](
auto conn) {
1529 std::array<double, 9> coords;
1534 auto get_closet_face = [
this](
auto *point_ptr,
auto r) {
1536 std::vector<EntityHandle> faces_out;
1540 "get closest faces");
1544 auto get_faces_out = [
this](
auto *point_ptr,
auto *unit_ray_ptr,
auto radius,
1546 std::vector<double> distances_out;
1547 std::vector<EntityHandle> faces_out;
1552 point_ptr, unit_ray_ptr, &radius),
1554 "get closest faces");
1555 return std::make_pair(faces_out, distances_out);
1558 auto get_normal = [](
auto &ele_coords) {
1564 auto make_map = [&](
auto &face_out,
auto &face_dist,
auto &t_ray_point,
1565 auto &t_unit_ray,
auto &t_master_coord) {
1568 std::map<double, EntityHandle>
m;
1569 for (
auto ii = 0; ii != face_out.size(); ++ii) {
1570 auto face_conn = get_face_conn(face_out[ii]);
1573 t_face_normal.normalize();
1575 t_x(
i) = t_ray_point(
i) + t_unit_ray(
i) * face_dist[ii];
1578 t_x(
i) * t_face_normal(
j) - t_master_coord(
i) * t_unit_ray(
j);
1579 if (t_unit_ray(
i) * t_face_normal(
i) > std::cos(M_PI / 3)) {
1580 auto dot = std::sqrt(t_chi(
i,
j) * t_chi(
i,
j));
1581 m[dot] = face_out[ii];
1587 auto get_tag_data = [
this](
auto tag,
auto face,
auto &vec) {
1591 vec.resize(tag_length);
1597 auto create_tag = [
this](
const std::string tag_name,
const int size) {
1598 double def_VAL[] = {0, 0, 0, 0, 0, 0, 0, 0, 0};
1602 tag_name.c_str(), size, MB_TYPE_DOUBLE,
th,
1603 MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
1610 auto set_float_precision = [](
const double x) {
1611 if (std::abs(x) < std::numeric_limits<float>::epsilon())
1618 auto save_scal_tag = [&](
auto &
th,
auto v,
const int gg) {
1621 v = set_float_precision(
v);
1628 auto get_fe_adjacencies = [
this](
auto fe_ent) {
1631 &fe_ent, 1, 2,
false, adj_faces, moab::Interface::UNION),
1633 std::set<int> adj_ids;
1634 for (
auto f : adj_faces) {
1640 auto get_face_id = [
this](
auto face) {
1649 auto get_body_id = [
this](
auto face) {
1658 auto get_face_part = [
this](
auto face) {
1659 ParallelComm *pcomm_post_proc_mesh = ParallelComm::get_pcomm(
1663 pcomm_post_proc_mesh->part_tag(), &face, 1, &part) == MB_SUCCESS) {
1669 auto check_face = [&](
auto face,
auto fe_id,
auto part) {
1670 auto face_id = get_face_id(face);
1671 auto face_part = get_face_part(face);
1672 if (face_id == fe_id && face_part == part)
1680 auto save_vec_tag = [&](
auto &
th,
auto &t_d,
const int gg) {
1684 for (
auto &
a :
v.data())
1685 a = set_float_precision(
a);
1687 &*
v.data().begin());
1692 Tag th_mark = create_tag(
"contact_mark", 1);
1693 Tag th_mark_slave = create_tag(
"contact_mark_slave", 1);
1694 Tag th_body_id = create_tag(
"contact_body_id", 1);
1695 Tag th_gap = create_tag(
"contact_gap", 1);
1696 Tag th_tn_master = create_tag(
"contact_tn_master", 1);
1697 Tag th_tn_slave = create_tag(
"contact_tn_slave", 1);
1698 Tag th_contact_traction = create_tag(
"contact_traction", 3);
1699 Tag th_contact_traction_master = create_tag(
"contact_traction_master", 3);
1700 Tag th_contact_traction_slave = create_tag(
"contact_traction_slave", 3);
1701 Tag th_c = create_tag(
"contact_c", 1);
1702 Tag th_normal = create_tag(
"contact_normal", 3);
1703 Tag th_dist = create_tag(
"contact_dip", 3);
1705 auto t_ele_centre = get_ele_centre(getFTensor1Coords());
1706 auto ele_radius = get_ele_radius(t_ele_centre, getFTensor1Coords());
1712 auto adj_fe_ids = get_fe_adjacencies(fe_ent);
1714 for (
auto gg = 0; gg != nb_gauss_pts; ++gg) {
1717 t_spatial_coords(
i) = t_coords(
i) + t_disp_h1(
i);
1720 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1723 auto faces_close = get_closet_face(&t_spatial_coords(0), ele_radius);
1724 for (
auto face_close : faces_close) {
1725 if (check_face(face_close, fe_id, m_field.get_comm_rank())) {
1727 auto body_id = get_body_id(face_close);
1729 auto master_face_conn = get_face_conn(face_close);
1730 std::array<double, 9> master_coords;
1733 master_coords.data());
1734 std::array<double, 9> master_traction;
1737 master_traction.data());
1738 auto t_normal_face_close = get_normal(master_coords);
1739 t_normal_face_close.normalize();
1743 CHKERR save_scal_tag(th_mark,
m, gg);
1744 CHKERR save_scal_tag(th_body_id,
static_cast<double>(body_id), gg);
1745 CHKERR save_vec_tag(th_normal, t_normal_face_close, gg);
1746 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1750 t_unit_ray(
i) = -t_normal_face_close(
i);
1753 t_spatial_coords(
i) -
1756 constexpr double eps = 1e-3;
1757 auto [faces_out, faces_dist] =
1758 get_faces_out(&t_ray_point(0), &t_unit_ray(0),
1762 auto m = make_map(faces_out, faces_dist, t_ray_point, t_unit_ray,
1764 for (
auto m_it =
m.begin(); m_it !=
m.end(); ++m_it) {
1765 auto face = m_it->second;
1766 if (face != face_close) {
1770 (adj_fe_ids.find(get_face_id(face)) == adj_fe_ids.end() ||
1771 get_face_part(face) != m_field.get_comm_rank())
1776 shadow_vec.back().gaussPtNb = gg;
1778 auto slave_face_conn = get_face_conn(face);
1779 std::array<double, 9> slave_coords;
1782 slave_coords.data());
1783 auto t_normal_face = get_normal(slave_coords);
1784 std::array<double, 9> slave_tractions;
1787 slave_tractions.data());
1789 auto t_master_point =
1790 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data());
1791 auto t_slave_point =
1792 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data());
1793 auto t_ray_point_data =
1794 getFTensor1FromPtr<3>(shadow_vec.back().rayPoint.data());
1795 auto t_unit_ray_data =
1796 getFTensor1FromPtr<3>(shadow_vec.back().unitRay.data());
1798 t_slave_point(
i) = t_ray_point(
i) + m_it->first * t_unit_ray(
i);
1800 auto eval_position = [&](
auto &&t_elem_coords,
auto &&t_point) {
1801 std::array<double, 2> loc_coords;
1804 &t_elem_coords(0, 0), &t_point(0), 1,
1806 "get local coords");
1815 t_point_out(
i) = t_shape_fun(
j) * t_elem_coords(
j,
i);
1819 auto t_master_point_updated = eval_position(
1820 getFTensor2FromPtr<3, 3>(master_coords.data()),
1821 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data()));
1822 t_master_point(
i) = t_master_point_updated(
i);
1824 auto t_slave_point_updated = eval_position(
1825 getFTensor2FromPtr<3, 3>(slave_coords.data()),
1826 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data()));
1827 t_slave_point(
i) = t_slave_point_updated(
i);
1829 t_ray_point_data(
i) = t_ray_point(
i);
1830 t_unit_ray_data(
i) = t_unit_ray(
i);
1832 std::copy(master_coords.begin(), master_coords.end(),
1833 shadow_vec.back().masterPointNodes.begin());
1834 std::copy(master_traction.begin(), master_traction.end(),
1835 shadow_vec.back().masterTractionNodes.begin());
1836 std::copy(slave_coords.begin(), slave_coords.end(),
1837 shadow_vec.back().slavePointNodes.begin());
1838 std::copy(slave_tractions.begin(), slave_tractions.end(),
1839 shadow_vec.back().slaveTractionNodes.begin());
1841 shadow_vec.back().eleRadius = ele_radius;
1851 auto [gap, tn_master, tn_slave,
c, t_master_traction,
1853 multiGetGap(&(shadow_vec.back()), t_spatial_coords);
1855 t_gap_vec(
i) = t_slave_point(
i) - t_spatial_coords(
i);
1856 CHKERR save_scal_tag(th_gap, gap, gg);
1857 CHKERR save_scal_tag(th_tn_master, tn_master, gg);
1858 CHKERR save_scal_tag(th_tn_slave, tn_slave, gg);
1859 CHKERR save_scal_tag(th_c,
c, gg);
1861 CHKERR save_scal_tag(th_mark_slave,
m, gg);
1862 CHKERR save_vec_tag(th_dist, t_gap_vec, gg);
1863 CHKERR save_vec_tag(th_contact_traction_master,
1864 t_master_traction, gg);
1865 CHKERR save_vec_tag(th_contact_traction_slave, t_slave_traction,
1883 const std::string block_name,
int dim) {
1889 std::regex((boost::format(
"%s(.*)") % block_name).str())
1893 for (
auto bc : bcs) {
1897 "get meshset ents");
1904boost::shared_ptr<ForcesAndSourcesCore>
1907 auto &m_field = ep.
mField;
1909 boost::shared_ptr<ContactTree> fe_contact_tree;
1916 std::map<int, Range> map;
1922 (boost::format(
"%s(.*)") % name).str()
1929 m_field.get_moab(), dim, ents,
true),
1931 map[m_ptr->getMeshsetId()] = ents;
1932 MOFEM_LOG(
"EPSYNC", sev) <<
"Meshset: " << m_ptr->getMeshsetId() <<
" "
1933 << ents.size() <<
" entities";
1940 auto get_map_skin = [&](
auto &&map) {
1941 ParallelComm *pcomm =
1942 ParallelComm::get_pcomm(&m_field.get_moab(),
MYPCOMM_INDEX);
1944 Skinner skin(&m_field.get_moab());
1945 for (
auto &
m : map) {
1947 CHKERR skin.find_skin(0,
m.second,
false, skin_faces);
1949 skin_faces, PSTATUS_SHARED | PSTATUS_MULTISHARED,
1950 PSTATUS_NOT, -1,
nullptr),
1952 m.second.swap(skin_faces);
1962 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1964 auto calcs_side_traction = [&](
auto &pip) {
1968 using SideEleOp = EleOnSide::UserDataOperator;
1971 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1972 boost::make_shared<CGGUserPolynomialBase>();
1973 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
1974 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
1976 op_loop_domain_side->getOpPtrVector().push_back(
1978 ep.
piolaStress, contact_common_data_ptr->contactTractionPtr(),
1979 boost::make_shared<double>(1.0)));
1980 pip.push_back(op_loop_domain_side);
1984 auto add_contact_three = [&]() {
1986 auto tree_moab_ptr = boost::make_shared<moab::Core>();
1987 fe_contact_tree = boost::make_shared<ContactTree>(
1990 fe_contact_tree->getOpPtrVector().push_back(
1992 ep.
contactDisp, contact_common_data_ptr->contactDispPtr()));
1993 CHKERR calcs_side_traction(fe_contact_tree->getOpPtrVector());
1994 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
1995 fe_contact_tree->getOpPtrVector().push_back(
1997 fe_contact_tree->getOpPtrVector().push_back(
1998 new OpMoveNode(fe_contact_tree, contact_common_data_ptr, u_h1_ptr));
2002 CHKERR add_contact_three();
2009 struct exclude_sdf {
2010 exclude_sdf(
Range &&r) : map(r) {}
2011 bool operator()(
FEMethod *fe_method_ptr) {
2013 if (map.find(ent) != map.end()) {
2023 fe_contact_tree->exeTestHook =
2026 return fe_contact_tree;
2031 std::map<int, Range> map;
2036 (boost::format(
"%s(.*)") % name).str()
2045 map[m_ptr->getMeshsetId()] = ents;
2052 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2053 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2056 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2058 auto &m_field = ep.
mField;
2064 using SideEleOp = EleOnSide::UserDataOperator;
2065 using BdyEleOp = BoundaryEle::UserDataOperator;
2072 auto rule_contact = [](int, int,
int o) {
return -1; };
2075 auto set_rule_contact = [refine](
2078 int order_col,
int order_data
2082 auto rule = 2 * order_data;
2087 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2088 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2090 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2096 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2099 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2101 auto add_ops_domain_side = [&](
auto &pip) {
2106 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2107 boost::make_shared<CGGUserPolynomialBase>();
2109 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
2110 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
2112 op_loop_domain_side->getOpPtrVector().push_back(
2115 op_loop_domain_side->getOpPtrVector().push_back(
2117 ep.
piolaStress, contact_common_data_ptr->contactTractionPtr()));
2118 pip.push_back(op_loop_domain_side);
2122 auto add_ops_contact_rhs = [&](
auto &pip) {
2125 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2129 ep.
contactDisp, contact_common_data_ptr->contactDispPtr()));
2130 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2134 contact_tree_ptr, u_h1_ptr,
2135 contact_common_data_ptr->contactTractionPtr(),
2139 ep.
contactDisp, contact_common_data_ptr, contact_tree_ptr,
2140 contact_sfd_map_range_ptr));
2142 broken_data_ptr, contact_common_data_ptr, contact_tree_ptr));
2148 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2149 CHKERR add_ops_contact_rhs(op_loop_skeleton_side->getOpPtrVector());
2152 pip.push_back(op_loop_skeleton_side);
2158 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2159 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2162 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2163 auto &m_field = ep.
mField;
2169 using SideEleOp = EleOnSide::UserDataOperator;
2170 using BdyEleOp = BoundaryEle::UserDataOperator;
2177 auto rule_contact = [](int, int,
int o) {
return -1; };
2180 auto set_rule_contact = [refine](
2183 int order_col,
int order_data
2187 auto rule = 2 * order_data;
2192 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2193 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2195 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2201 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2204 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2206 auto add_ops_domain_side = [&](
auto &pip) {
2211 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2212 boost::make_shared<CGGUserPolynomialBase>();
2214 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
2215 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
2217 op_loop_domain_side->getOpPtrVector().push_back(
2220 op_loop_domain_side->getOpPtrVector().push_back(
2222 ep.
piolaStress, contact_common_data_ptr->contactTractionPtr()));
2223 pip.push_back(op_loop_domain_side);
2227 auto add_ops_contact_lhs = [&](
auto &pip) {
2230 ep.
contactDisp, contact_common_data_ptr->contactDispPtr()));
2231 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2235 contact_tree_ptr, u_h1_ptr,
2236 contact_common_data_ptr->contactTractionPtr(),
2241 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2246 contact_tree_ptr, contact_sfd_map_range_ptr));
2249 ep.
contactDisp, broken_data_ptr, contact_common_data_ptr,
2250 contact_tree_ptr, contact_sfd_map_range_ptr));
2253 broken_data_ptr, ep.
contactDisp, contact_common_data_ptr,
2260 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2261 CHKERR add_ops_contact_lhs(op_loop_skeleton_side->getOpPtrVector());
2264 pip.push_back(op_loop_skeleton_side);
2271 boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2272 boost::shared_ptr<MatrixDouble> u_h1_ptr,
2273 boost::shared_ptr<MatrixDouble> contact_traction_ptr,
2274 Range r, moab::Interface *post_proc_mesh_ptr,
2275 std::vector<EntityHandle> *map_gauss_pts_ptr) {
2277 auto &m_field = ep.
mField;
2278 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2280 contact_tree_ptr, u_h1_ptr, contact_traction_ptr,
2282 post_proc_mesh_ptr, map_gauss_pts_ptr);
Implementation of tonsorial bubble base div(v) = 0.
Eshelbian plasticity interface.
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define FTENSOR_INDEX(DIM, I)
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Tensor1< T, Tensor_Dim > normalize()
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MYPCOMM_INDEX
default communicator number PCOMM
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define CHK_MOAB_THROW(err, msg)
Check error code of MoAB function and throw MoFEM exception.
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static const double face_coords[4][9]
MoFEMErrorCode writeFile(const std::string file_name)
wrote results in (MOAB) format, use "file_name.h5m"
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
MoFEMErrorCode getCubitMeshsetPtr(const int ms_id, const CubitBCType cubit_bc_type, const CubitMeshSets **cubit_meshset_ptr) const
get cubit meshset
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
boost::shared_ptr< ContactSDFPython > setupContactSdf()
Read SDF file and setup contact SDF.
static auto get_range_from_block(MoFEM::Interface &m_field, const std::string block_name, int dim)
auto checkSdf(EntityHandle fe_ent, std::map< int, Range > &sdf_map_range)
ForcesAndSourcesCore::UserDataOperator * getOpContactDetection(EshelbianCore &ep, boost::shared_ptr< ForcesAndSourcesCore > contact_tree_ptr, boost::shared_ptr< MatrixDouble > u_h1_ptr, boost::shared_ptr< MatrixDouble > contact_traction_ptr, Range r, moab::Interface *post_proc_mesh_ptr, std::vector< EntityHandle > *map_gauss_pts_ptr)
Push operator for contact detection.
boost::shared_ptr< ForcesAndSourcesCore > createContactDetectionFiniteElement(EshelbianCore &ep)
Create a Contact Tree finite element.
auto multiMasterPoint(ContactTree::FaceData *face_data_ptr, FTensor::Tensor1< T1, 3 > &t_spatial_coords)
auto multiGetGap(ContactTree::FaceData *face_data_ptr, FTensor::Tensor1< T1, 3 > &t_spatial_coords)
auto multiPointRhs(ContactTree::FaceData *face_data_ptr, FTensor::Tensor1< T1, 3 > &t_coords, FTensor::Tensor1< T2, 3 > &t_spatial_coords, FTensor::Tensor1< T3, 3 > &t_master_traction, MultiPointRhsType type, bool debug=false)
MoFEMErrorCode pushContactOpsRhs(EshelbianCore &ep, boost::shared_ptr< ForcesAndSourcesCore > contact_tree_ptr, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
Push contact operations to the right-hand side.
auto multiSlavePoint(ContactTree::FaceData *face_data_ptr, FTensor::Tensor1< T1, 3 > &t_spatial_coords)
auto multiPoint(std::array< double, 3 > &unit_ray, std::array< double, 3 > &point, std::array< double, 9 > &elem_point_nodes, std::array< double, 9 > &elem_traction_nodes, FTensor::Tensor1< T1, 3 > &t_spatial_coords)
Calculate points data on contact surfaces.
MoFEMErrorCode pushContactOpsLhs(EshelbianCore &ep, boost::shared_ptr< ForcesAndSourcesCore > contact_tree_ptr, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
Push contact operations to the left-hand side.
static auto get_body_range(MoFEM::Interface &m_field, const std::string name, int dim)
auto getSdf(OP_PTR op_ptr, MatrixDouble &contact_disp, int block_id, bool eval_hessian)
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
VectorBoundedArray< double, 3 > VectorDouble3
UBlasMatrix< double > MatrixDouble
implementation of Data Operators for Forces and Sources
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
auto id_from_handle(const EntityHandle h)
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
void tricircumcenter3d_tp(double a[3], double b[3], double c[3], double circumcenter[3], double *xi, double *eta)
MoFEM::OpBrokenTopoBase GAUSS
Gaussian quadrature integration.
static auto getFTensor0FromVec(ublas::vector< T, A > &data)
Get tensor rank 0 (scalar) form data vector.
PetscErrorCode PetscOptionsGetString(PetscOptions *, const char pre[], const char name[], char str[], size_t size, PetscBool *set)
FTensor::Index< 'm', 3 > m
boost::shared_ptr< Range > frontAdjEdges
static double dynamicTime
MoFEM::Interface & mField
const std::string materialH1Positions
const std::string elementVolumeName
static PetscBool dynamicRelaxation
const std::string spatialH1Disp
const std::string piolaStress
int contactRefinementLevels
const std::string contactDisp
const std::string contactElement
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
boost::shared_ptr< ContactTree > contactTreePtr
boost::shared_ptr< ContactTree > contactTreePtr
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
boost::shared_ptr< std::map< int, Range > > sdfMapRangePtr
boost::shared_ptr< ContactTree > contactTreePtr
boost::shared_ptr< std::map< int, Range > > sdfMapRangePtr
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpConstrainBoundaryL2Lhs_dU(const std::string row_field_name, const std::string col_field_name, boost::shared_ptr< ContactOps::CommonData > common_data_ptr, boost::shared_ptr< ContactTree > contact_tree_ptr, boost::shared_ptr< std::map< int, Range > > sdf_map_range_ptr=nullptr)
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
boost::shared_ptr< ContactTree > contactTreePtr
boost::shared_ptr< std::map< int, Range > > sdfMapRangePtr
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
boost::shared_ptr< ContactTree > contactTreePtr
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
OpConstrainBoundaryL2Rhs(const std::string row_field_name, boost::shared_ptr< ContactOps::CommonData > common_data_ptr, boost::shared_ptr< ContactTree > contact_tree_ptr, boost::shared_ptr< std::map< int, Range > > sdf_map_range_ptr=nullptr)
boost::shared_ptr< ContactOps::CommonData > commonDataPtr
boost::shared_ptr< ContactTree > contactTreePtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
FaceElementForcesAndSourcesCore::UserDataOperator UOP
boost::shared_ptr< MatrixDouble > uH1Ptr
OpMoveNode(boost::shared_ptr< ContactTree > contact_tree_ptr, boost::shared_ptr< ContactOps::CommonData > common_data_ptr, boost::shared_ptr< MatrixDouble > u_h1_ptr)
moab::Interface * postProcMeshPtr
boost::shared_ptr< ContactTree > contactTreePtr
std::vector< EntityHandle > * mapGaussPtsPtr
boost::shared_ptr< MatrixDouble > uH1Ptr
OpTreeSearch(boost::shared_ptr< ContactTree > contact_tree_ptr, boost::shared_ptr< MatrixDouble > u_h1_ptr, boost::shared_ptr< MatrixDouble > traction_ptr, Range r, moab::Interface *post_proc_mesh_ptr, std::vector< EntityHandle > *map_gauss_pts_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< MatrixDouble > tractionPtr
FaceElementForcesAndSourcesCore::UserDataOperator UOP
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual int get_comm_rank() const =0
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
FTensor::Tensor0< FTensor::PackPtr< double *, 1 > > getFTensor0N(const FieldApproximationBase base)
Get base function as Tensor0.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
Structure for user loop methods on finite elements.
EntityHandle getFEEntityHandle() const
Get the entity handle of the current finite element.
Base face element used to integrate on skeleton.
default operator for TRI element
structure to get information from mofem into EntitiesFieldData
MatrixDouble gaussPts
Matrix of integration points.
Interface for managing meshsets containing materials and boundary conditions.
Operator for broken loop side.
Calculate trace of vector (Hdiv/Hcurl) space.
Specialization for double precision vector field values calculation.
Element used to execute operators on side of the element.
Template struct for dimension-specific finite element types.
MoFEMErrorCode preProcess()
Generate vertices and elements.
std::string optionsPrefix
Prefix for options.
auto & getPostProcMesh()
Get postprocessing mesh.
MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
std::map< EntityType, PostProcGenerateRefMeshPtr > refElementsMap
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
BoundaryEle::UserDataOperator BdyEleOp
double zeta
Viscous hardening.