v0.15.5
Loading...
Searching...
No Matches
EshelbianContact.cpp
Go to the documentation of this file.
1/**
2 * @file EshelbianContact.cpp
3 * @brief
4 * @date 2023-05-13
5 *
6 * @copyright Copyright (c) 2023
7 *
8 */
9
10#define SINGULARITY
11#include <MoFEM.hpp>
12using namespace MoFEM;
13
15
17
18namespace ContactOps {
19
22using BoundaryEleOp = BoundaryEle::UserDataOperator;
24
27
28double cn_contact = 1;
31double scale = 1;
32
33double airplane_ray_distance = 1; // thi is point from which plane send ray.
34 // This is multiple of elem radius.
35
36} // namespace ContactOps
37
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;
44#endif
45
46#include <ContactOps.hpp>
47
48#include <EshelbianContact.hpp>
49
50extern "C" {
51void tricircumcenter3d_tp(double a[3], double b[3], double c[3],
52 double circumcenter[3], double *xi, double *eta);
53}
54
55namespace EshelbianPlasticity {
56
57struct ContactSDFPython : public ContactOps::SDFPython {
58 using ContactOps::SDFPython::SDFPython;
59};
60
61boost::shared_ptr<ContactSDFPython> setupContactSdf() {
62 boost::shared_ptr<ContactSDFPython> sdf_python_ptr;
63
64#ifdef ENABLE_PYTHON_BINDING
65
66 auto file_exists = [](std::string myfile) {
67 std::ifstream file(myfile.c_str());
68 if (file) {
69 return true;
70 }
71 return false;
72 };
73
74 char sdf_file_name[255] = "sdf.py";
75 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR, "-sdf_file",
76 sdf_file_name, 255, PETSC_NULLPTR);
77
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";
84 } else {
85 MOFEM_LOG("EP", Sev::warning) << sdf_file_name << " file NOT found";
86 }
87
88#endif // ENABLE_PYTHON_BINDING
89
90 return sdf_python_ptr;
91};
92
94 : public PostProcBrokenMeshInMoabBase<FaceElementForcesAndSourcesCore> {
95
98
100 boost::shared_ptr<moab::Core> core_mesh_ptr, int max_order,
101 std::map<int, Range> &&body_map);
103 if (treeSurfPtr) {
104 treeSurfPtr->delete_tree(rootSetSurf);
105 }
106 }
109 inline boost::shared_ptr<OrientedBoxTreeTool> &getTreeSurfPtr() {
110 return treeSurfPtr;
111 }
112 inline auto getRootSetSurf() { return rootSetSurf; }
113 int getMaxLevel() const { return refElementsMap.at(MBTRI)->defMaxLevel; }
114
115 friend struct OpMoveNode;
116 friend struct OpTreeSearch;
117
118 struct FaceData {
119 int gaussPtNb; //< integration points number
120 std::array<double, 3> slavePoint;
121 std::array<double, 3> masterPoint;
122 std::array<double, 3> rayPoint;
123 std::array<double, 3> unitRay;
124 double eleRadius;
125
126 // std::vector<int> dofsSlaveIds;
127 // std::vector<double> dofsSlaveCoeff;
128 // std::vector<double> baseSlaveFuncs;
129
130 std::array<double, 9> masterPointNodes;
131 std::array<double, 9> masterTractionNodes;
132 std::array<double, 9> slavePointNodes;
133 std::array<double, 9> slaveTractionNodes;
134
135 FaceData() = default;
136 };
137
138 using MapFaceData = std::map<EntityHandle, std::vector<FaceData>>;
139
140 inline auto findFaceDataVecPtr(EntityHandle fe_ent) {
141 auto &map_face_data = shadowDataMap;
142 auto it = map_face_data.find(fe_ent);
143 if (it == map_face_data.end()) {
144 return (std::vector<FaceData> *)nullptr;
145 }
146 return &(it->second);
147 }
148
149 inline auto getFaceDataPtr(std::vector<FaceData>::iterator &it, int gg,
150 std::vector<FaceData> *vec_ptr) {
151 FaceData *face_data_ptr = nullptr;
152 if (it != vec_ptr->end()) {
153 if (it->gaussPtNb == gg) {
154 face_data_ptr = &(*it);
155 ++it;
156 }
157 }
158 return face_data_ptr;
159 }
160
161protected:
163 boost::shared_ptr<OrientedBoxTreeTool> treeSurfPtr;
165
167 // Tag thCoeff;
168 // Tag thIds;
169 // Tag thBases;
174
176
177 std::map<int, Range> bodyMap;
178
179 const int maxOrder;
180};
181
182auto checkSdf(EntityHandle fe_ent, std::map<int, Range> &sdf_map_range) {
183 for (auto &m_sdf : sdf_map_range) {
184 if (m_sdf.second.find(fe_ent) != m_sdf.second.end())
185 return m_sdf.first;
186 }
187 return -1;
188}
189
190template <typename OP_PTR>
191auto getSdf(OP_PTR op_ptr, MatrixDouble &contact_disp, int block_id,
192 bool eval_hessian) {
193
194 auto nb_gauss_pts = op_ptr->getGaussPts().size2();
195
196 auto ts_time = op_ptr->getTStime();
197 auto ts_time_step = op_ptr->getTStimeStep();
200 ts_time_step = EshelbianCore::dynamicStep;
201 }
202
203 auto m_spatial_coords = ContactOps::get_spatial_coords(
204 op_ptr->getFTensor1CoordsAtGaussPts(),
205 getFTensor1FromMat<3>(contact_disp), nb_gauss_pts);
206 auto m_normals_at_pts = ContactOps::get_normalize_normals(
207 op_ptr->getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
209 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
210 block_id);
212 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
213 block_id);
214
215#ifndef NDEBUG
216 if (v_sdf.size() != nb_gauss_pts)
218 "Wrong number of integration pts");
219 if (m_grad_sdf.size2() != nb_gauss_pts)
221 "Wrong number of integration pts");
222 if (m_grad_sdf.size1() != 3)
223 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Should be size of 3");
224#endif // NDEBUG
225
226 if (eval_hessian) {
228 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
229 block_id);
230 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
231 m_hess_sdf);
232 } else {
233
234 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
235 MatrixDouble(6, nb_gauss_pts, 0.));
236 }
237};
238
239/**
240 * @brief Calculate points data on contact surfaces
241 *
242 * @tparam T1
243 * @param unit_ray
244 * @param point
245 * @param elem_point_nodes
246 * @param elem_traction_nodes
247 * @param t_spatial_coords
248 * @return auto
249 */
250template <typename T1>
252
253 std::array<double, 3> &unit_ray, std::array<double, 3> &point,
254
255 std::array<double, 9> &elem_point_nodes,
256 std::array<double, 9> &elem_traction_nodes,
257
258 FTensor::Tensor1<T1, 3> &t_spatial_coords) {
259 FTensor::Index<'i', 3> i;
260
261 auto t_unit_ray = getFTensor1FromPtr<3>(unit_ray.data());
262 auto t_point = getFTensor1FromPtr<3>(point.data());
263
264 auto get_normal = [](auto &ele_coords) {
266 Tools::getTriNormal(ele_coords.data(), &t_normal(0));
267 return t_normal;
268 };
269
270 auto t_normal = get_normal(elem_point_nodes);
271 t_normal(i) /= t_normal.l2();
272
273 auto sn = t_normal(i) * t_point(i);
274 auto nm = t_normal(i) * t_spatial_coords(i);
275 auto nr = t_normal(i) * t_unit_ray(i);
276
277 auto gamma = (sn - nm) / nr;
278
279 FTensor::Tensor1<T1, 3> t_point_current;
280 t_point_current(i) = t_spatial_coords(i) + gamma * t_unit_ray(i);
281
282 auto get_local_point_shape_functions = [&](auto &&t_elem_coords,
283 auto &t_point) {
284 std::array<T1, 2> loc_coords;
287 &t_elem_coords(0, 0), &t_point(0), 1, loc_coords.data()),
288 "get local coords");
289 return FTensor::Tensor1<T1, 3>{N_MBTRI0(loc_coords[0], loc_coords[1]),
290 N_MBTRI1(loc_coords[0], loc_coords[1]),
291 N_MBTRI2(loc_coords[0], loc_coords[1])};
292 };
293
294 auto eval_position = [&](auto &&t_field, auto &t_shape_fun) {
295 FTensor::Index<'i', 3> i;
296 FTensor::Index<'j', 3> j;
297 FTensor::Tensor1<T1, 3> t_point_out;
298 t_point_out(i) = t_shape_fun(j) * t_field(j, i);
299 return t_point_out;
300 };
301
302 auto t_shape_fun = get_local_point_shape_functions(
303 getFTensor2FromPtr<3, 3>(elem_point_nodes.data()), t_point_current);
304 auto t_slave_point_updated = eval_position(
305 getFTensor2FromPtr<3, 3>(elem_point_nodes.data()), t_shape_fun);
306 auto t_traction_updated = eval_position(
307 getFTensor2FromPtr<3, 3>(elem_traction_nodes.data()), t_shape_fun);
308
309 return std::make_tuple(t_slave_point_updated, t_traction_updated, t_normal);
310};
311
312template <typename T1>
314
315 ContactTree::FaceData *face_data_ptr,
316 FTensor::Tensor1<T1, 3> &t_spatial_coords
317
318) {
319
320 return multiPoint(face_data_ptr->unitRay, face_data_ptr->masterPoint,
321 face_data_ptr->masterPointNodes,
322 face_data_ptr->masterTractionNodes, t_spatial_coords);
323};
324
325template <typename T1>
327
328 ContactTree::FaceData *face_data_ptr,
329 FTensor::Tensor1<T1, 3> &t_spatial_coords
330
331) {
332
333 return multiPoint(face_data_ptr->unitRay, face_data_ptr->slavePoint,
334 face_data_ptr->slavePointNodes,
335 face_data_ptr->slaveTractionNodes, t_spatial_coords);
336};
337
338/**
339 * Evaluate gap and tractions between master and slave points
340 */
341template <typename T1>
343 FTensor::Tensor1<T1, 3> &t_spatial_coords) {
344
345 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
346 multiSlavePoint(face_data_ptr, t_spatial_coords);
347 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
348 multiMasterPoint(face_data_ptr, t_spatial_coords);
349
350 FTensor::Index<'i', 3> i;
351
353 t_normal(i) = t_master_normal(i) - t_slave_normal(i);
354 t_normal.normalize();
355
356 auto gap = t_normal(i) * (t_slave_point_current(i) - t_spatial_coords(i));
357 auto tn_master = t_master_traction_current(i) * t_normal(i);
358 auto tn_slave = t_slave_traction_current(i) * t_normal(i);
359 auto tn = std::max(-tn_master, tn_slave);
360
361 return std::make_tuple(gap, tn_master, tn_slave,
362 ContactOps::constrain(gap, tn),
363 t_master_traction_current, t_slave_traction_current);
364};
365
367
368/** Caluclate rhs term for contact
369 */
370template <typename T1, typename T2, typename T3>
372
373 ContactTree::FaceData *face_data_ptr, FTensor::Tensor1<T1, 3> &t_coords,
374 FTensor::Tensor1<T2, 3> &t_spatial_coords,
375 FTensor::Tensor1<T3, 3> &t_master_traction, MultiPointRhsType type,
376 bool debug = false
377
378) {
379 FTensor::Index<'i', 3> i;
380 FTensor::Index<'j', 3> j;
381
383 t_u(i) = t_spatial_coords(i) - t_coords(i);
384
386
387 switch (type) {
389
390 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
391 multiSlavePoint(face_data_ptr, t_spatial_coords);
392 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
393 multiMasterPoint(face_data_ptr, t_spatial_coords);
394
395 // average normal surface
397 t_normal(i) = t_master_normal(i) - t_slave_normal(i);
398 t_normal.normalize();
399
400 // get projection operators
402 t_P(i, j) = t_normal(i) * t_normal(j);
404 t_Q(i, j) = kronecker_delta(i, j) - t_P(i, j);
405
406 constexpr double beta = 0.5;
407
408 auto zeta = 1e-8 * face_data_ptr->eleRadius;
411
412 // this is regularised std::min(d,s)
413 auto f_min_gap = [zeta](auto d, auto s) {
414 return 0.5 * (d + s - std::sqrt((d - s) * (d - s) + zeta));
415 };
416 // this derivative of regularised std::min(d,s)
417 auto f_diff_min_gap = [zeta](auto d, auto s) {
418 return 0.5 * (1 - (d - s) / std::sqrt((d - s) * (d - s) + zeta));
419 };
420
421 // add penalty on penetration side
422 auto f_barrier = [alpha1, alpha2, f_min_gap, f_diff_min_gap](auto g,
423 auto tn) {
424 auto d = alpha1 * g;
425 auto b1 =
426 0.5 * (tn + f_min_gap(d, tn) + f_diff_min_gap(d, tn) * (d - tn));
427 auto b2 = alpha2 * f_min_gap(g, 0) * g;
428 return b1 - b2;
429 };
430
431 FTensor::Tensor1<T2, 3> t_gap_vec;
432 t_gap_vec(i) = beta * t_spatial_coords(i) +
433 (1 - beta) * t_slave_point_current(i) - t_spatial_coords(i);
435 t_traction_vec(i) =
436 -beta * t_master_traction(i) + (beta - 1) * t_slave_traction_current(i);
437
438 auto t_gap = t_normal(i) * t_gap_vec(i);
439 auto t_tn = t_normal(i) * t_traction_vec(i);
440 auto barrier = f_barrier(t_gap, t_tn);
441
443 // add penalty on penetration side
444 t_traction_bar(i) = t_normal(i) * barrier;
445
446 t_rhs(i) =
447
448 t_Q(i, j) * t_master_traction(j) +
449
450 t_P(i, j) * (t_master_traction(j) - t_traction_bar(j));
451
452 if (debug) {
453 auto is_nan_or_inf = [](double value) -> bool {
454 return std::isnan(value) || std::isinf(value);
455 };
456
457 double v = std::complex<double>(t_rhs(i) * t_rhs(i)).real();
458 if (is_nan_or_inf(v)) {
459 MOFEM_LOG_CHANNEL("SELF");
460 MOFEM_LOG("SELF", Sev::error) << "t_rhs " << t_rhs;
461 CHK_MOAB_THROW(MOFEM_DATA_INCONSISTENCY, "rhs is nan or inf");
462 }
463 }
464
465 } break;
468 break;
469 }
470
471 return t_rhs;
472};
473
476 const std::string row_field_name,
477 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
478 boost::shared_ptr<ContactTree> contact_tree_ptr,
479 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr = nullptr);
480
482
483private:
484 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
485 boost::shared_ptr<ContactTree> contactTreePtr;
486 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
487};
488
490 const std::string row_field_name,
491 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
492 boost::shared_ptr<ContactTree> contact_tree_ptr,
493 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
494 : ContactOps::AssemblyBoundaryEleOp(row_field_name, row_field_name,
495 ContactOps::BoundaryEleOp::OPROW),
496 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
497 sdfMapRangePtr(sdf_map_range_ptr) {
498 CHK_THROW_MESSAGE(PetscOptionsGetScalar(PETSC_NULLPTR, "", "-cn",
500 PETSC_NULLPTR),
501 "get cn failed");
504 PetscOptionsGetScalar(PETSC_NULLPTR, "", "-alpha_contact_const",
505 &ContactOps::alpha_contact_const, PETSC_NULLPTR),
506 "get alpha contact failed");
508 PETSC_NULLPTR, "", "-alpha_contact_quadratic",
510 "get alpha contact failed");
512 PetscOptionsGetScalar(PETSC_NULLPTR, "", "-airplane_ray_distance",
513 &ContactOps::airplane_ray_distance, PETSC_NULLPTR),
514 "get alpha contact failed");
515
516 MOFEM_LOG("EP", Sev::inform) << "cn " << ContactOps::cn_contact;
517 MOFEM_LOG("EP", Sev::inform)
518 << "alpha_contact_const " << ContactOps::alpha_contact_const;
519 MOFEM_LOG("EP", Sev::inform)
520 << "alpha_contact_quadratic " << ContactOps::alpha_contact_quadratic;
521 MOFEM_LOG("EP", Sev::inform)
522 << "airplane_ray_distance " << ContactOps::airplane_ray_distance;
523}
524
528
529 FTensor::Index<'i', 3> i;
530 FTensor::Index<'j', 3> j;
531 FTensor::Index<'k', 3> k;
532 FTensor::Index<'l', 3> l;
533
534 const size_t nb_gauss_pts = getGaussPts().size2();
535
536#ifndef NDEBUG
537 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
538 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
539 "Wrong number of integration pts %ld != %ld",
540 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
541 }
542#endif // !NDEBUG
543
544 auto &nf = locF;
545 locF.clear();
546
547 auto t_w = getFTensor0IntegrationWeight();
548 auto t_coords = getFTensor1CoordsAtGaussPts();
549 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
550 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
551
552 // placeholder to pass boundary block id to python. default SDF is set on
553 // block = -1, one can choose different block by making block "CONTACT_SDF",
554 // then specific SDF can be set to that block.
555 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
556 getSdf(this, commonDataPtr->contactDisp,
557 checkSdf(getFEEntityHandle(), *sdfMapRangePtr), false);
558
559 auto t_sdf_v = getFTensor0FromVec(v_sdf);
560 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
561 auto t_normalize_normal = getFTensor1FromMat<3>(m_normals_at_pts);
562
563 auto next = [&]() {
564 ++t_w;
565 ++t_coords;
566 ++t_disp;
567 ++t_traction;
568 ++t_normalize_normal;
569 ++t_sdf_v;
570 ++t_grad_sdf_v;
571 };
572
573 auto face_data_vec_ptr =
574 contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
575 auto face_gauss_pts_it = face_data_vec_ptr->begin();
576
577 auto nb_base_functions = data.getN().size2();
578 auto t_base = data.getFTensor0N();
579 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
580
582 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
583 face_data_vec_ptr);
584
585 auto check_face_contact = [&]() {
586 if (checkSdf(getFEEntityHandle(), *sdfMapRangePtr) != -1)
587 return false;
588
589 if (face_data_ptr) {
590 return true;
591 }
592 return false;
593 };
594
595#ifdef ENABLE_PYTHON_BINDING
596 double c = 0.;
597 if (ContactOps::sdfPythonWeakPtr.lock()) {
598 auto tn = t_traction(i) * t_grad_sdf_v(i);
599 c = ContactOps::constrain(t_sdf_v, tn);
600 }
601#else
602 constexpr double c = 0;
603#endif
604
605 if (!c && check_face_contact()) {
606 FTensor::Tensor1<double, 3> t_spatial_coords;
607 t_spatial_coords(i) = t_coords(i) + t_disp(i);
608 auto t_rhs_tmp = multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
609 t_traction, MultiPointRhsType::U, true);
610 t_rhs(i) = t_rhs_tmp(i);
611
612 } else {
613
614#ifdef ENABLE_PYTHON_BINDING
615 auto inv_cn = 1. / ContactOps::cn_contact;
616
617 if (ContactOps::sdfPythonWeakPtr.lock()) {
619 t_cP(i, j) = (c * t_grad_sdf_v(i)) * t_grad_sdf_v(j);
620 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
621 t_rhs(i) = t_cQ(i, j) * t_traction(j) +
622 (c * inv_cn * t_sdf_v) * t_grad_sdf_v(i);
623 } else {
624 t_rhs(i) = t_traction(i);
625 }
626#else
627 t_rhs(i) = t_traction(i);
628#endif
629 }
630
631 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
632 const double alpha = t_w * getMeasure();
633
634 size_t bb = 0;
635 for (; bb != nbRows / 3; ++bb) {
636 const double beta = alpha * t_base;
637 t_nf(i) -= beta * t_rhs(i);
638 ++t_nf;
639 ++t_base;
640 }
641 for (; bb < nb_base_functions; ++bb)
642 ++t_base;
643
644 next();
645 }
646
648}
649
650template <AssemblyType A, IntegrationType I> struct OpConstrainBoundaryHDivRhs;
651
652template <AssemblyType A>
654 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
655
658
660 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
661 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
662 boost::shared_ptr<ContactTree> contact_tree_ptr);
663
664 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
665
666private:
667 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
668 boost::shared_ptr<ContactTree> contactTreePtr;
669};
670
671template <AssemblyType A>
674 boost::shared_ptr<std::vector<BrokenBaseSideData>>
675 broken_base_side_data,
676 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
677 boost::shared_ptr<ContactTree> contact_tree_ptr)
678 : OP(broken_base_side_data), commonDataPtr(common_data_ptr),
679 contactTreePtr(contact_tree_ptr) {}
680
681template <AssemblyType A>
685
690
691 const size_t nb_gauss_pts = OP::getGaussPts().size2();
692
693#ifndef NDEBUG
694 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
695 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
696 "Wrong number of integration pts %ld != %ld",
697 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
698 }
699#endif // !NDEBUG
700
701 auto &nf = OP::locF;
702 OP::locF.clear();
703
704 auto t_w = OP::getFTensor0IntegrationWeight();
705 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
706 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
707
708 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
709 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
710
711 auto next = [&]() {
712 ++t_w;
713 ++t_disp;
714 ++t_traction;
715 ++t_coords;
716 ++t_material_normal;
717 };
718
719 auto face_data_vec_ptr =
720 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
721 auto face_gauss_pts_it = face_data_vec_ptr->begin();
722
723 auto nb_base_functions = data.getN().size2() / 3;
724 auto t_base = data.getFTensor1N<3>();
725 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
726
727 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
728 const double alpha = t_w / 2.;
729
730 size_t bb = 0;
731 for (; bb != OP::nbRows / 3; ++bb) {
732 const double beta = alpha * t_base(i) * t_material_normal(i);
733 t_nf(i) += beta * t_disp(i);
734 ++t_nf;
735 ++t_base;
736 }
737 for (; bb < nb_base_functions; ++bb)
738 ++t_base;
739
740 next();
741 }
742
744}
745
747
749 const std::string row_field_name, const std::string col_field_name,
750 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
751 boost::shared_ptr<ContactTree> contact_tree_ptr,
752 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr = nullptr);
753
756
757private:
758 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
759 boost::shared_ptr<ContactTree> contactTreePtr;
760 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
761};
762
764 const std::string row_field_name, const std::string col_field_name,
765 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
766 boost::shared_ptr<ContactTree> contact_tree_ptr,
767 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
768 : ContactOps::AssemblyBoundaryEleOp(row_field_name, col_field_name,
769 ContactOps::BoundaryEleOp::OPROWCOL),
770 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
771 sdfMapRangePtr(sdf_map_range_ptr) {
772
773 sYmm = false;
774}
775
778 EntitiesFieldData::EntData &col_data) {
780
781 using namespace ContactOps;
782
783 FTensor::Index<'i', 3> i;
784 FTensor::Index<'j', 3> j;
785 FTensor::Index<'k', 3> k;
786
787 auto nb_rows = row_data.getIndices().size();
788 auto nb_cols = col_data.getIndices().size();
789
790 auto &locMat = AssemblyBoundaryEleOp::locMat;
791 locMat.resize(nb_rows, nb_cols, false);
792 locMat.clear();
793
794 if (nb_cols && nb_rows) {
795
796 auto nb_gauss_pts = getGaussPts().size2();
797 auto t_w = getFTensor0IntegrationWeight();
798 auto t_coords = getFTensor1CoordsAtGaussPts();
799 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
800 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
801
802 // placeholder to pass boundary block id to python
803 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
804 getSdf(this, commonDataPtr->contactDisp,
805 checkSdf(getFEEntityHandle(), *sdfMapRangePtr), true);
806
807 auto t_sdf_v = getFTensor0FromVec(v_sdf);
808 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
809 auto t_hess_sdf_v = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
810 auto t_normalized_normal = getFTensor1FromMat<3>(m_normals_at_pts);
811
812 auto next = [&]() {
813 ++t_w;
814 ++t_coords;
815 ++t_disp;
816 ++t_traction;
817 ++t_sdf_v;
818 ++t_grad_sdf_v;
819 ++t_hess_sdf_v;
820 ++t_normalized_normal;
821 };
822
823 auto face_data_vec_ptr =
824 contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
825 auto face_gauss_pts_it = face_data_vec_ptr->begin();
826
827 auto t_row_base = row_data.getFTensor0N();
828 auto nb_face_functions = row_data.getN().size2() / 3;
829 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
830
831 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
832
833 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
834 face_data_vec_ptr);
835
836 auto check_face_contact = [&]() {
837 if (checkSdf(getFEEntityHandle(), *sdfMapRangePtr) != -1)
838 return false;
839
840 if (face_data_ptr) {
841 return true;
842 }
843 return false;
844 };
845
847
848#ifdef ENABLE_PYTHON_BINDING
849 double c = 0.;
850 if (ContactOps::sdfPythonWeakPtr.lock()) {
851 auto tn = t_traction(i) * t_grad_sdf_v(i);
852 c = ContactOps::constrain(t_sdf_v, tn);
853 }
854#else
855 constexpr double c = 0;
856#endif
857
858 if (!c && check_face_contact()) {
859 FTensor::Tensor1<double, 3> t_spatial_coords;
860 t_spatial_coords(i) = t_coords(i) + t_disp(i);
861 constexpr double eps = std::numeric_limits<float>::epsilon();
862 for (auto ii = 0; ii < 3; ++ii) {
863 FTensor::Tensor1<std::complex<double>, 3> t_spatial_coords_cx{
864 t_spatial_coords(0), t_spatial_coords(1), t_spatial_coords(2)};
865 t_spatial_coords_cx(ii) += eps * 1i;
866 auto t_rhs_tmp =
867 multiPointRhs(face_data_ptr, t_coords, t_spatial_coords_cx,
868 t_traction, MultiPointRhsType::U);
869 for (int jj = 0; jj != 3; ++jj) {
870 auto v = t_rhs_tmp(jj).imag();
871 t_res_dU(jj, ii) = v / eps;
872 }
873 }
874
875 } else {
876
877#ifdef ENABLE_PYTHON_BINDING
878
879 if (ContactOps::sdfPythonWeakPtr.lock()) {
880 auto inv_cn = 1. / ContactOps::cn_contact;
881 t_res_dU(i, j) =
882
883 (-c) * (t_hess_sdf_v(i, j) * t_grad_sdf_v(k) * t_traction(k) +
884 t_grad_sdf_v(i) * t_hess_sdf_v(k, j) * t_traction(k))
885
886 + (c * inv_cn) * (t_sdf_v * t_hess_sdf_v(i, j) +
887
888 t_grad_sdf_v(j) * t_grad_sdf_v(i));
889 } else {
890 t_res_dU(i, j) = 0;
891 }
892#else
893 t_res_dU(i, j) = 0;
894#endif
895 }
896
897 auto alpha = t_w * getMeasure();
898
899 size_t rr = 0;
900 for (; rr != nb_rows / 3; ++rr) {
901
902 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
903 auto t_col_base = col_data.getFTensor0N(gg, 0);
904
905 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
906 auto beta = alpha * t_row_base * t_col_base;
907 t_mat(i, j) -= beta * t_res_dU(i, j);
908 ++t_col_base;
909 ++t_mat;
910 }
911
912 ++t_row_base;
913 }
914 for (; rr < nb_face_functions; ++rr)
915 ++t_row_base;
916
917 next();
918 }
919 }
920
922}
923
924template <AssemblyType A, IntegrationType I> struct OpConstrainBoundaryL2Lhs_dP;
925
926template <AssemblyType A>
928 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
929
932
934 std::string row_field_name,
935 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
936 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
937 boost::shared_ptr<ContactTree> contact_tree_ptr,
938 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr = nullptr);
939
940 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
942
943private:
944 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
945 boost::shared_ptr<ContactTree> contactTreePtr;
946 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
947};
948
949template <AssemblyType A>
952 std::string row_field_name,
953 boost::shared_ptr<std::vector<BrokenBaseSideData>>
954 broken_base_side_data,
955 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
956 boost::shared_ptr<ContactTree> contact_tree_ptr,
957 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
958 : OP(row_field_name, broken_base_side_data, false, false),
959 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
960 sdfMapRangePtr(sdf_map_range_ptr) {
961 OP::sYmm = false;
962}
963
964template <AssemblyType A>
968 EntitiesFieldData::EntData &col_data) {
970
971 using namespace ContactOps;
972
973 FTensor::Index<'i', 3> i;
974 FTensor::Index<'j', 3> j;
975 FTensor::Index<'k', 3> k;
976
977 auto nb_rows = row_data.getIndices().size();
978 auto nb_cols = col_data.getIndices().size();
979
980 auto &locMat = AssemblyBoundaryEleOp::locMat;
981 locMat.resize(nb_rows, nb_cols, false);
982 locMat.clear();
983
984 if (nb_cols && nb_rows) {
985
986 const size_t nb_gauss_pts = OP::getGaussPts().size2();
987
988 auto t_w = OP::getFTensor0IntegrationWeight();
989 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
990 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
991 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
992 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
993
994 // placeholder to pass boundary block id to python
995 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
996 getSdf(this, commonDataPtr->contactDisp,
997 checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr), false);
998
999 auto t_sdf_v = getFTensor0FromVec(v_sdf);
1000 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
1001
1002 auto next = [&]() {
1003 ++t_w;
1004 ++t_disp;
1005 ++t_traction;
1006 ++t_coords;
1007 ++t_material_normal;
1008 ++t_sdf_v;
1009 ++t_grad_sdf_v;
1010 };
1011
1012 auto face_data_vec_ptr =
1013 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1014 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1015
1016 auto t_row_base = row_data.getFTensor0N();
1017 auto nb_face_functions = row_data.getN().size2();
1018
1019 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
1020
1021 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1022
1023 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
1024 face_data_vec_ptr);
1025
1026 auto check_face_contact = [&]() {
1027 if (checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr) != -1)
1028 return false;
1029
1030 if (face_data_ptr) {
1031 return true;
1032 }
1033 return false;
1034 };
1035
1037
1038#ifdef ENABLE_PYTHON_BINDING
1039 double c = 0.;
1040 if (ContactOps::sdfPythonWeakPtr.lock()) {
1041 auto tn = t_traction(i) * t_grad_sdf_v(i);
1042 c = ContactOps::constrain(t_sdf_v, tn);
1043 }
1044#else
1045 constexpr double c = 0;
1046#endif
1047
1048 if (!c && check_face_contact()) {
1049 FTensor::Tensor1<double, 3> t_spatial_coords;
1050 t_spatial_coords(i) = t_coords(i) + t_disp(i);
1051 constexpr double eps = std::numeric_limits<float>::epsilon();
1052 for (auto ii = 0; ii != 3; ++ii) {
1053 FTensor::Tensor1<std::complex<double>, 3> t_traction_cx{
1054 t_traction(0), t_traction(1), t_traction(2)};
1055 t_traction_cx(ii) += eps * 1i;
1056 auto t_rhs_tmp =
1057 multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
1058 t_traction_cx, MultiPointRhsType::U);
1059 for (int jj = 0; jj != 3; ++jj) {
1060 auto v = t_rhs_tmp(jj).imag();
1061 t_res_dP(jj, ii) = v / eps;
1062 }
1063 }
1064 } else {
1065
1066#ifdef ENABLE_PYTHON_BINDING
1067 if (ContactOps::sdfPythonWeakPtr.lock()) {
1069 t_cP(i, j) = (c * t_grad_sdf_v(i)) * t_grad_sdf_v(j);
1070 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
1071 t_res_dP(i, j) = t_cQ(i, j);
1072 } else {
1073 t_res_dP(i, j) = t_kd(i, j);
1074 }
1075#else
1076 t_res_dP(i, j) = t_kd(i, j);
1077#endif
1078 }
1079
1080 const double alpha = t_w / 2.;
1081 size_t rr = 0;
1082 for (; rr != nb_rows / 3; ++rr) {
1083
1084 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1085 auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
1086
1087 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
1088 auto col_base = t_col_base(i) * t_material_normal(i);
1089 const double beta = alpha * t_row_base * col_base;
1090 t_mat(i, j) -= beta * t_res_dP(i, j);
1091 ++t_col_base;
1092 ++t_mat;
1093 }
1094
1095 ++t_row_base;
1096 }
1097 for (; rr < nb_face_functions; ++rr)
1098 ++t_row_base;
1099
1100 next();
1101 }
1102 }
1103
1105}
1106
1107template <AssemblyType A, IntegrationType I>
1109
1110template <AssemblyType A>
1112 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
1113
1116
1118 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1119 std::string col_field_name,
1120 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1121 boost::shared_ptr<ContactTree> contact_tree_ptr);
1122
1123 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
1124 EntitiesFieldData::EntData &col_data);
1125
1126private:
1127 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
1128 boost::shared_ptr<ContactTree> contactTreePtr;
1129};
1130
1131template <AssemblyType A>
1134 boost::shared_ptr<std::vector<BrokenBaseSideData>>
1135 broken_base_side_data,
1136 std::string col_field_name,
1137 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1138 boost::shared_ptr<ContactTree> contact_tree_ptr)
1139 : OP(col_field_name, broken_base_side_data, true, true),
1140 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr) {
1141 OP::sYmm = false;
1142}
1143
1144template <AssemblyType A>
1148 EntitiesFieldData::EntData &row_data) {
1150
1151 // Note: col_data and row_data are swapped in this function, we going to
1152 // transpose locMat at the end
1153
1154 using namespace ContactOps;
1155
1159
1160 auto nb_rows = row_data.getIndices().size();
1161 auto nb_cols = col_data.getIndices().size();
1162
1163 auto &locMat = AssemblyBoundaryEleOp::locMat;
1164 locMat.resize(nb_rows, nb_cols, false);
1165 locMat.clear();
1166
1167 if (nb_cols && nb_rows) {
1168
1169 const size_t nb_gauss_pts = OP::getGaussPts().size2();
1170
1171 auto t_w = OP::getFTensor0IntegrationWeight();
1172 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
1173 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
1174 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
1175 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
1176
1177 auto next = [&]() {
1178 ++t_w;
1179 ++t_disp;
1180 ++t_traction;
1181 ++t_coords;
1182 ++t_material_normal;
1183 };
1184
1185 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
1186
1187 auto face_data_vec_ptr =
1188 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1189 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1190
1191 auto t_row_base = row_data.getFTensor1N<3>();
1192 auto nb_face_functions = row_data.getN().size2() / 3;
1193 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1194
1195 const auto alpha = t_w / 2.;
1196
1197 size_t rr = 0;
1198 for (; rr != nb_rows / 3; ++rr) {
1199
1200 auto row_base = alpha * (t_row_base(i) * t_material_normal(i));
1201
1202 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1203 auto t_col_base = col_data.getFTensor0N(gg, 0);
1204
1205 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
1206 const auto beta = row_base * t_col_base;
1207 t_mat(i, j) += beta * t_kd(i, j);
1208 ++t_col_base;
1209 ++t_mat;
1210 }
1211
1212 ++t_row_base;
1213 }
1214 for (; rr < nb_face_functions; ++rr)
1215 ++t_row_base;
1216
1217 next();
1218 }
1219 }
1220
1221 locMat = trans(locMat);
1222
1224}
1225
1227 boost::shared_ptr<moab::Core> core_mesh_ptr,
1228 int max_order, std::map<int, Range> &&body_map)
1229 : Base(m_field, core_mesh_ptr, "contact"), maxOrder(max_order),
1230 bodyMap(body_map) {
1231
1232 auto ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
1233 ref_ele_ptr->hoNodes =
1234 PETSC_FALSE; ///< So far only linear geometry is implemented for contact
1235
1236 CHK_THROW_MESSAGE(ref_ele_ptr->getOptions(optionsPrefix), "getOptions");
1237 CHK_THROW_MESSAGE(ref_ele_ptr->generateReferenceElementMesh(),
1238 "Error when generating reference element");
1239
1240 MOFEM_LOG("EP", Sev::inform) << "Contact hoNodes " << ref_ele_ptr->hoNodes
1241 ? "true"
1242 : "false";
1243 MOFEM_LOG("EP", Sev::inform)
1244 << "Contact maxOrder " << ref_ele_ptr->defMaxLevel;
1245
1246 refElementsMap[MBTRI] = ref_ele_ptr;
1247
1248 int def_ele_id = -1;
1249 CHKERR getPostProcMesh().tag_get_handle("ELE_ID", 1, MB_TYPE_INTEGER, thEleId,
1250 MB_TAG_CREAT | MB_TAG_DENSE,
1251 &def_ele_id);
1252 CHKERR getPostProcMesh().tag_get_handle("BODY_ID", 1, MB_TYPE_INTEGER,
1253 thBodyId, MB_TAG_CREAT | MB_TAG_DENSE,
1254 &def_ele_id);
1255
1256 // std::vector<int> def_ids(3 * nb_bases, 0);
1257 // CHKERR getPostProcMesh().tag_get_handle("IDS", 3 * nb_bases,
1258 // MB_TYPE_INTEGER,
1259 // thIds, MB_TAG_CREAT | MB_TAG_DENSE,
1260 // &*def_ids.begin());
1261 // std::vector<double> def_coeffs(3 * nb_bases, 0);
1262 // CHKERR getPostProcMesh().tag_get_handle("COEFF", 3 * nb_bases,
1263 // MB_TYPE_DOUBLE,
1264 // thCoeff, MB_TAG_CREAT |
1265 // MB_TAG_DENSE,
1266 // &*def_coeffs.begin());
1267 // std::vector<double> def_basses(nb_bases, 0);
1268 // CHKERR getPostProcMesh().tag_get_handle("BASES", nb_bases, MB_TYPE_DOUBLE,
1269 // thBases, MB_TAG_CREAT |
1270 // MB_TAG_DENSE,
1271 // &*def_basses.begin());
1272
1273 std::array<double, 3> def_small_x{0., 0., 0.};
1274 CHKERR getPostProcMesh().tag_get_handle("x", 3, MB_TYPE_DOUBLE, thSmallX,
1275 MB_TAG_CREAT | MB_TAG_DENSE,
1276 def_small_x.data());
1277 CHKERR getPostProcMesh().tag_get_handle("X", 3, MB_TYPE_DOUBLE, thLargeX,
1278 MB_TAG_CREAT | MB_TAG_DENSE,
1279 def_small_x.data());
1280
1281 std::array<double, 3> def_tractions{0., 0., 0.};
1282 CHKERR getPostProcMesh().tag_get_handle(
1283 "TRACTION", 3, MB_TYPE_DOUBLE, thTraction, MB_TAG_CREAT | MB_TAG_DENSE,
1284 &*def_tractions.begin());
1285}
1286
1292
1295
1297 shadowDataMap.clear();
1298
1299 PetscBarrier(nullptr);
1300
1301 auto pcomm_post_proc_mesh = this->getPostProcMeshPcommPtr();
1302 if (pcomm_post_proc_mesh == nullptr)
1303 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
1304
1305 auto brodacts = [&](auto &brodacts_ents) {
1307 Range recived_ents;
1308 for (auto rr = 0; rr != mField.get_comm_size(); ++rr) {
1309 pcomm_post_proc_mesh->broadcast_entities(
1310 rr, rr == mField.get_comm_rank() ? brodacts_ents : recived_ents,
1311 false, true);
1312 }
1314 };
1315
1316 Range brodacts_ents;
1317 CHKERR getPostProcMesh().get_entities_by_dimension(0, 2, brodacts_ents, true);
1318 CHKERR brodacts(brodacts_ents);
1319
1320 Range ents;
1321 CHKERR getPostProcMesh().get_entities_by_dimension(0, 2, ents, true);
1322 CHKERR buildTree(ents);
1323
1324#ifndef NDEBUG
1325 CHKERR writeFile("debug_tree.h5m");
1326#endif // NDEBUG
1327
1329}
1330
1333 // MOFEM_LOG("SELF", Sev::inform) << ents << endl;
1334 if (treeSurfPtr) {
1335 treeSurfPtr->delete_tree(rootSetSurf);
1336 }
1337 treeSurfPtr = boost::shared_ptr<OrientedBoxTreeTool>(
1338 new OrientedBoxTreeTool(&getPostProcMesh(), "ROOTSETSURF", true));
1339 CHKERR treeSurfPtr->build(ents, rootSetSurf);
1340 // CHKERR treeSurfPtr->stats(rootSetSurf, std::cerr);
1342}
1343
1345
1347
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);
1351 MoFEMErrorCode doWork(int side, EntityType type,
1353
1354protected:
1355 boost::weak_ptr<ContactTree> contactTreePtr;
1356 boost::shared_ptr<MatrixDouble> uH1Ptr;
1357 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
1358};
1359
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) {}
1366
1367MoFEMErrorCode OpMoveNode::doWork(int side, EntityType type,
1370
1371 if (auto contact_tree_ptr = contactTreePtr.lock()) {
1372
1373 auto get_body_id = [&](auto fe_ent) {
1374 for (auto &m : contact_tree_ptr->bodyMap) {
1375 if (m.second.find(fe_ent) != m.second.end()) {
1376 return m.first;
1377 }
1378 }
1379 return -1;
1380 };
1381
1382 auto &moab_post_proc_mesh = contact_tree_ptr->getPostProcMesh();
1383 auto &post_proc_ents = contact_tree_ptr->getPostProcElements();
1384
1385 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1386 auto fe_id = id_from_handle(fe_ent);
1387 auto body_id = get_body_id(fe_ent);
1388 auto &map_gauss_pts = contact_tree_ptr->getMapGaussPts();
1389
1390 CHKERR moab_post_proc_mesh.tag_clear_data(contact_tree_ptr->thEleId,
1391 post_proc_ents, &fe_id);
1392 CHKERR moab_post_proc_mesh.tag_clear_data(contact_tree_ptr->thBodyId,
1393 post_proc_ents, &body_id);
1394
1395 auto nb_gauss_pts = getGaussPts().size2();
1396 auto t_u_h1 = getFTensor1FromMat<3>(*uH1Ptr);
1397 auto t_u_l2 = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
1398 auto t_coords = getFTensor1CoordsAtGaussPts();
1399
1400 MatrixDouble x_h1(nb_gauss_pts, 3);
1401 auto t_x_h1 = getFTensor1FromPtr<3>(&x_h1(0, 0));
1402 MatrixDouble x_l2(nb_gauss_pts, 3);
1403 auto t_x_l2 = getFTensor1FromPtr<3>(&x_l2(0, 0));
1404 MatrixDouble tractions = trans(commonDataPtr->contactTraction);
1405 MatrixDouble coords = getCoordsAtGaussPts();
1406
1407 FTensor::Index<'i', 3> i;
1408
1409 // VectorDouble bases(nb_bases, 0);
1410 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1411 t_x_h1(i) = t_coords(i) + t_u_h1(i);
1412 t_x_l2(i) = t_coords(i) + t_u_l2(i);
1413
1414 ++t_coords;
1415 ++t_u_h1;
1416 ++t_u_l2;
1417 ++t_x_h1;
1418 ++t_x_l2;
1419 }
1420
1421 CHKERR moab_post_proc_mesh.set_coords(
1422 &*map_gauss_pts.begin(), map_gauss_pts.size(), &*x_h1.data().begin());
1423 CHKERR moab_post_proc_mesh.tag_set_data(
1424 contact_tree_ptr->thSmallX, &*map_gauss_pts.begin(),
1425 map_gauss_pts.size(), &*x_h1.data().begin());
1426 CHKERR moab_post_proc_mesh.tag_set_data(
1427 contact_tree_ptr->thLargeX, &*map_gauss_pts.begin(),
1428 map_gauss_pts.size(), &*coords.data().begin());
1429 CHKERR moab_post_proc_mesh.tag_set_data(
1430 contact_tree_ptr->thTraction, &*map_gauss_pts.begin(),
1431 map_gauss_pts.size(), &*tractions.data().begin());
1432
1433 } else {
1434 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY,
1435 "ContactTree pointer expired in OpMoveNode");
1436 }
1437
1439}
1440
1442
1444
1445 OpTreeSearch(boost::shared_ptr<ContactTree> contact_tree_ptr,
1446 boost::shared_ptr<MatrixDouble> u_h1_ptr,
1447 boost::shared_ptr<MatrixDouble> traction_ptr, Range r,
1448
1449 moab::Interface *post_proc_mesh_ptr,
1450 std::vector<EntityHandle> *map_gauss_pts_ptr
1451
1452 );
1453 MoFEMErrorCode doWork(int side, EntityType type,
1455
1456protected:
1457 boost::shared_ptr<ContactTree> contactTreePtr;
1458 boost::shared_ptr<MatrixDouble> uH1Ptr;
1459 boost::shared_ptr<MatrixDouble> tractionPtr;
1460 moab::Interface *postProcMeshPtr = nullptr;
1461 std::vector<EntityHandle> *mapGaussPtsPtr = nullptr;
1462
1464};
1465
1466OpTreeSearch::OpTreeSearch(boost::shared_ptr<ContactTree> contact_tree_ptr,
1467 boost::shared_ptr<MatrixDouble> u_h1_ptr,
1468 boost::shared_ptr<MatrixDouble> traction_ptr,
1469 Range r,
1470
1471 moab::Interface *post_proc_mesh_ptr,
1472 std::vector<EntityHandle> *map_gauss_pts_ptr
1473
1474 )
1475 : UOP(NOSPACE, UOP::OPSPACE), contactTreePtr(contact_tree_ptr),
1476 uH1Ptr(u_h1_ptr), tractionPtr(traction_ptr),
1477 postProcMeshPtr(post_proc_mesh_ptr), mapGaussPtsPtr(map_gauss_pts_ptr),
1478 contactRange(r) {}
1479
1480MoFEMErrorCode OpTreeSearch::doWork(int side, EntityType type,
1483
1484 auto &m_field = getPtrFE()->mField;
1485 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1486 auto fe_id = id_from_handle(fe_ent);
1487
1488 if (contactRange.find(fe_ent) == contactRange.end())
1490
1491 FTensor::Index<'i', 3> i;
1492 FTensor::Index<'j', 3> j;
1493
1494 const auto nb_gauss_pts = getGaussPts().size2();
1495
1496 auto t_disp_h1 = getFTensor1FromMat<3>(*uH1Ptr);
1497 auto t_coords = getFTensor1CoordsAtGaussPts();
1498 auto t_traction = getFTensor1FromMat<3>(*tractionPtr);
1499
1500 auto next = [&]() {
1501 ++t_disp_h1;
1502 ++t_traction;
1503 ++t_coords;
1504 };
1505
1506 auto get_ele_centre = [i](auto t_ele_coords) {
1507 FTensor::Tensor1<double, 3> t_ele_center;
1508 t_ele_center(i) = 0;
1509 for (int nn = 0; nn != 3; nn++) {
1510 t_ele_center(i) += t_ele_coords(i);
1511 ++t_ele_coords;
1512 }
1513 t_ele_center(i) /= 3;
1514 return t_ele_center;
1515 };
1516
1517 auto get_ele_radius = [i](auto t_ele_center, auto t_ele_coords) {
1519 t_n0(i) = t_ele_center(i) - t_ele_coords(i);
1520 return t_n0.l2();
1521 };
1522
1523 auto get_face_conn = [this](auto face) {
1524 const EntityHandle *conn;
1525 int num_nodes;
1526 CHK_MOAB_THROW(contactTreePtr->getPostProcMesh().get_connectivity(
1527 face, conn, num_nodes, true),
1528 "get conn");
1529 if (num_nodes != 3) {
1530 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "face is not a triangle");
1531 }
1532 return conn;
1533 };
1534
1535 auto get_face_coords = [this](auto conn) {
1536 std::array<double, 9> coords;
1537 CHKERR contactTreePtr->getPostProcMesh().get_coords(conn, 3, coords.data());
1538 return coords;
1539 };
1540
1541 auto get_closet_face = [this](auto *point_ptr, auto r) {
1542 FTensor::Tensor1<double, 3> t_point_out;
1543 std::vector<EntityHandle> faces_out;
1545 contactTreePtr->getTreeSurfPtr()->sphere_intersect_triangles(
1546 point_ptr, r / 8, contactTreePtr->getRootSetSurf(), faces_out),
1547 "get closest faces");
1548 return faces_out;
1549 };
1550
1551 auto get_faces_out = [this](auto *point_ptr, auto *unit_ray_ptr, auto radius,
1552 auto eps) {
1553 std::vector<double> distances_out;
1554 std::vector<EntityHandle> faces_out;
1556
1557 contactTreePtr->getTreeSurfPtr()->ray_intersect_triangles(
1558 distances_out, faces_out, contactTreePtr->getRootSetSurf(), eps,
1559 point_ptr, unit_ray_ptr, &radius),
1560
1561 "get closest faces");
1562 return std::make_pair(faces_out, distances_out);
1563 };
1564
1565 auto get_normal = [](auto &ele_coords) {
1567 Tools::getTriNormal(ele_coords.data(), &t_normal(0));
1568 return t_normal;
1569 };
1570
1571 auto make_map = [&](auto &face_out, auto &face_dist, auto &t_ray_point,
1572 auto &t_unit_ray, auto &t_master_coord) {
1573 FTensor::Index<'i', 3> i;
1574 FTensor::Index<'j', 3> j;
1575 std::map<double, EntityHandle> m;
1576 for (auto ii = 0; ii != face_out.size(); ++ii) {
1577 auto face_conn = get_face_conn(face_out[ii]);
1578 auto face_coords = get_face_coords(face_conn);
1579 auto t_face_normal = get_normal(face_coords);
1580 t_face_normal.normalize();
1582 t_x(i) = t_ray_point(i) + t_unit_ray(i) * face_dist[ii];
1584 t_chi(i, j) =
1585 t_x(i) * t_face_normal(j) - t_master_coord(i) * t_unit_ray(j);
1586 if (t_unit_ray(i) * t_face_normal(i) > std::cos(M_PI / 3)) {
1587 auto dot = std::sqrt(t_chi(i, j) * t_chi(i, j));
1588 m[dot] = face_out[ii];
1589 }
1590 }
1591 return m;
1592 };
1593
1594 auto get_tag_data = [this](auto tag, auto face, auto &vec) {
1596 int tag_length;
1597 CHKERR contactTreePtr->getPostProcMesh().tag_get_length(tag, tag_length);
1598 vec.resize(tag_length);
1599 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(tag, &face, 1,
1600 &*vec.begin());
1602 };
1603
1604 auto create_tag = [this](const std::string tag_name, const int size) {
1605 double def_VAL[] = {0, 0, 0, 0, 0, 0, 0, 0, 0};
1606 Tag th;
1607 if (postProcMeshPtr) {
1608 CHKERR postProcMeshPtr->tag_get_handle(
1609 tag_name.c_str(), size, MB_TYPE_DOUBLE, th,
1610 MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
1611 CHKERR postProcMeshPtr->tag_clear_data(th, &*mapGaussPtsPtr->begin(),
1612 mapGaussPtsPtr->size(), def_VAL);
1613 }
1614 return th;
1615 };
1616
1617 auto set_float_precision = [](const double x) {
1618 if (std::abs(x) < std::numeric_limits<float>::epsilon())
1619 return 0.;
1620 else
1621 return x;
1622 };
1623
1624 // scalars
1625 auto save_scal_tag = [&](auto &th, auto v, const int gg) {
1627 if (postProcMeshPtr) {
1628 v = set_float_precision(v);
1629 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1, &v);
1630 }
1632 };
1633
1634 // adjacencies
1635 auto get_fe_adjacencies = [this](auto fe_ent) {
1636 Range adj_faces;
1637 CHK_MOAB_THROW(getFaceFE()->mField.get_moab().get_adjacencies(
1638 &fe_ent, 1, 2, false, adj_faces, moab::Interface::UNION),
1639 "get adj");
1640 std::set<int> adj_ids;
1641 for (auto f : adj_faces) {
1642 adj_ids.insert(id_from_handle(f));
1643 }
1644 return adj_ids;
1645 };
1646
1647 auto get_face_id = [this](auto face) {
1648 int id;
1649 if (contactTreePtr->getPostProcMesh().tag_get_data(
1650 contactTreePtr->thEleId, &face, 1, &id) == MB_SUCCESS) {
1651 return id;
1652 }
1653 return -1;
1654 };
1655
1656 auto get_body_id = [this](auto face) {
1657 int id;
1658 if (contactTreePtr->getPostProcMesh().tag_get_data(
1659 contactTreePtr->thBodyId, &face, 1, &id) == MB_SUCCESS) {
1660 return id;
1661 }
1662 return -1;
1663 };
1664
1665 auto get_face_part = [this](auto face) {
1666 auto pcomm_post_proc_mesh = contactTreePtr->getPostProcMeshPcommPtr();
1667 if (pcomm_post_proc_mesh == nullptr)
1668 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
1669 int part;
1670 if (contactTreePtr->getPostProcMesh().tag_get_data(
1671 pcomm_post_proc_mesh->part_tag(), &face, 1, &part) == MB_SUCCESS) {
1672 return part;
1673 }
1674 return -1;
1675 };
1676
1677 auto check_face = [&](auto face, auto fe_id, auto part) {
1678 auto face_id = get_face_id(face);
1679 auto face_part = get_face_part(face);
1680 if (face_id == fe_id && face_part == part)
1681 return true;
1682 return false;
1683 };
1684
1685 // vectors
1686 VectorDouble3 v(3);
1687 FTensor::Tensor1<FTensor::PackPtr<double *, 0>, 3> t_v(&v[0], &v[1], &v[2]);
1688 auto save_vec_tag = [&](auto &th, auto &t_d, const int gg) {
1690 if (postProcMeshPtr) {
1691 t_v(i) = t_d(i);
1692 for (auto &a : v.data())
1693 a = set_float_precision(a);
1694 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1,
1695 &*v.data().begin());
1696 }
1698 };
1699
1700 Tag th_mark = create_tag("contact_mark", 1);
1701 Tag th_mark_slave = create_tag("contact_mark_slave", 1);
1702 Tag th_body_id = create_tag("contact_body_id", 1);
1703 Tag th_gap = create_tag("contact_gap", 1);
1704 Tag th_tn_master = create_tag("contact_tn_master", 1);
1705 Tag th_tn_slave = create_tag("contact_tn_slave", 1);
1706 Tag th_contact_traction = create_tag("contact_traction", 3);
1707 Tag th_contact_traction_master = create_tag("contact_traction_master", 3);
1708 Tag th_contact_traction_slave = create_tag("contact_traction_slave", 3);
1709 Tag th_c = create_tag("contact_c", 1);
1710 Tag th_normal = create_tag("contact_normal", 3);
1711 Tag th_dist = create_tag("contact_dip", 3);
1712
1713 auto t_ele_centre = get_ele_centre(getFTensor1Coords());
1714 auto ele_radius = get_ele_radius(t_ele_centre, getFTensor1Coords());
1715
1716 contactTreePtr->shadowDataMap.clear();
1717 auto &shadow_vec = contactTreePtr->shadowDataMap[fe_ent];
1718 shadow_vec.clear();
1719
1720 auto adj_fe_ids = get_fe_adjacencies(fe_ent);
1721
1722 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1723
1724 FTensor::Tensor1<double, 3> t_spatial_coords;
1725 t_spatial_coords(i) = t_coords(i) + t_disp_h1(i);
1726
1727 if (postProcMeshPtr) {
1728 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1729 }
1730
1731 auto faces_close = get_closet_face(&t_spatial_coords(0), ele_radius);
1732 for (auto face_close : faces_close) {
1733 if (check_face(face_close, fe_id, m_field.get_comm_rank())) {
1734
1735 auto body_id = get_body_id(face_close);
1736
1737 auto master_face_conn = get_face_conn(face_close);
1738 std::array<double, 9> master_coords;
1739 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1740 contactTreePtr->thSmallX, master_face_conn, 3,
1741 master_coords.data());
1742 std::array<double, 9> master_traction;
1743 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1744 contactTreePtr->thTraction, master_face_conn, 3,
1745 master_traction.data());
1746 auto t_normal_face_close = get_normal(master_coords);
1747 t_normal_face_close.normalize();
1748
1749 if (postProcMeshPtr) {
1750 double m = 1;
1751 CHKERR save_scal_tag(th_mark, m, gg);
1752 CHKERR save_scal_tag(th_body_id, static_cast<double>(body_id), gg);
1753 CHKERR save_vec_tag(th_normal, t_normal_face_close, gg);
1754 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1755 }
1756
1757 FTensor::Tensor1<double, 3> t_unit_ray;
1758 t_unit_ray(i) = -t_normal_face_close(i);
1759 FTensor::Tensor1<double, 3> t_ray_point;
1760 t_ray_point(i) =
1761 t_spatial_coords(i) -
1762 t_unit_ray(i) * ContactOps::airplane_ray_distance * ele_radius;
1763
1764 constexpr double eps = 1e-3;
1765 auto [faces_out, faces_dist] =
1766 get_faces_out(&t_ray_point(0), &t_unit_ray(0),
1767 2 * ContactOps::airplane_ray_distance * ele_radius,
1768 eps * ele_radius);
1769
1770 auto m = make_map(faces_out, faces_dist, t_ray_point, t_unit_ray,
1771 t_spatial_coords);
1772 for (auto m_it = m.begin(); m_it != m.end(); ++m_it) {
1773 auto face = m_it->second;
1774 if (face != face_close) {
1775
1776 if (
1777
1778 (adj_fe_ids.find(get_face_id(face)) == adj_fe_ids.end() ||
1779 get_face_part(face) != m_field.get_comm_rank())
1780
1781 ) {
1782
1783 shadow_vec.push_back(ContactTree::FaceData());
1784 shadow_vec.back().gaussPtNb = gg;
1785
1786 auto slave_face_conn = get_face_conn(face);
1787 std::array<double, 9> slave_coords;
1788 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1789 contactTreePtr->thSmallX, slave_face_conn, 3,
1790 slave_coords.data());
1791 auto t_normal_face = get_normal(slave_coords);
1792 std::array<double, 9> slave_tractions;
1793 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1794 contactTreePtr->thTraction, slave_face_conn, 3,
1795 slave_tractions.data());
1796
1797 auto t_master_point =
1798 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data());
1799 auto t_slave_point =
1800 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data());
1801 auto t_ray_point_data =
1802 getFTensor1FromPtr<3>(shadow_vec.back().rayPoint.data());
1803 auto t_unit_ray_data =
1804 getFTensor1FromPtr<3>(shadow_vec.back().unitRay.data());
1805
1806 t_slave_point(i) = t_ray_point(i) + m_it->first * t_unit_ray(i);
1807
1808 auto eval_position = [&](auto &&t_elem_coords, auto &&t_point) {
1809 std::array<double, 2> loc_coords;
1812 &t_elem_coords(0, 0), &t_point(0), 1,
1813 loc_coords.data()),
1814 "get local coords");
1815 FTensor::Tensor1<double, 3> t_shape_fun;
1816 CHK_THROW_MESSAGE(Tools::shapeFunMBTRI<0>(&t_shape_fun(0),
1817 &loc_coords[0],
1818 &loc_coords[1], 1),
1819 "calc shape fun");
1820 FTensor::Index<'i', 3> i;
1821 FTensor::Index<'j', 3> j;
1822 FTensor::Tensor1<double, 3> t_point_out;
1823 t_point_out(i) = t_shape_fun(j) * t_elem_coords(j, i);
1824 return t_point_out;
1825 };
1826
1827 auto t_master_point_updated = eval_position(
1828 getFTensor2FromPtr<3, 3>(master_coords.data()),
1829 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data()));
1830 t_master_point(i) = t_master_point_updated(i);
1831
1832 auto t_slave_point_updated = eval_position(
1833 getFTensor2FromPtr<3, 3>(slave_coords.data()),
1834 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data()));
1835 t_slave_point(i) = t_slave_point_updated(i);
1836
1837 t_ray_point_data(i) = t_ray_point(i);
1838 t_unit_ray_data(i) = t_unit_ray(i);
1839
1840 std::copy(master_coords.begin(), master_coords.end(),
1841 shadow_vec.back().masterPointNodes.begin());
1842 std::copy(master_traction.begin(), master_traction.end(),
1843 shadow_vec.back().masterTractionNodes.begin());
1844 std::copy(slave_coords.begin(), slave_coords.end(),
1845 shadow_vec.back().slavePointNodes.begin());
1846 std::copy(slave_tractions.begin(), slave_tractions.end(),
1847 shadow_vec.back().slaveTractionNodes.begin());
1848
1849 shadow_vec.back().eleRadius = ele_radius;
1850
1851 // CHKERR get_tag_data(contactTreePtr->thIds, face,
1852 // shadow_vec.back().dofsSlaveIds);
1853 // CHKERR get_tag_data(contactTreePtr->thCoeff, face,
1854 // shadow_vec.back().dofsSlaveCoeff);
1855 // CHKERR get_tag_data(contactTreePtr->thBases, face,
1856 // shadow_vec.back().baseSlaveFuncs);
1857
1858 if (postProcMeshPtr) {
1859 auto [gap, tn_master, tn_slave, c, t_master_traction,
1860 t_slave_traction] =
1861 multiGetGap(&(shadow_vec.back()), t_spatial_coords);
1863 t_gap_vec(i) = t_slave_point(i) - t_spatial_coords(i);
1864 CHKERR save_scal_tag(th_gap, gap, gg);
1865 CHKERR save_scal_tag(th_tn_master, tn_master, gg);
1866 CHKERR save_scal_tag(th_tn_slave, tn_slave, gg);
1867 CHKERR save_scal_tag(th_c, c, gg);
1868 double m = 1;
1869 CHKERR save_scal_tag(th_mark_slave, m, gg);
1870 CHKERR save_vec_tag(th_dist, t_gap_vec, gg);
1871 CHKERR save_vec_tag(th_contact_traction_master,
1872 t_master_traction, gg);
1873 CHKERR save_vec_tag(th_contact_traction_slave, t_slave_traction,
1874 gg);
1875 }
1876
1877 break;
1878 }
1879 }
1880 }
1881 break;
1882 }
1883 }
1884 next();
1885 }
1886
1888}
1889
1891 const std::string block_name, int dim) {
1892 Range r;
1893
1894 auto mesh_mng = m_field.getInterface<MeshsetsManager>();
1895 auto bcs = mesh_mng->getCubitMeshsetPtr(
1896
1897 std::regex((boost::format("%s(.*)") % block_name).str())
1898
1899 );
1900
1901 for (auto bc : bcs) {
1902 Range faces;
1903 CHK_MOAB_THROW(bc->getMeshsetIdEntitiesByDimension(m_field.get_moab(), dim,
1904 faces, true),
1905 "get meshset ents");
1906 r.merge(faces);
1907 }
1908
1909 return r;
1910};
1911
1912boost::shared_ptr<ForcesAndSourcesCore>
1914
1915 auto &m_field = ep.mField;
1916
1917 boost::shared_ptr<ContactTree> fe_contact_tree;
1918
1919 auto impl = [&]() {
1921
1922 /** Contact requires that body is marked */
1923 auto get_body_range = [&](auto name, int dim, auto sev) {
1924 std::map<int, Range> map;
1925
1926 for (auto m_ptr :
1927 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
1928 std::regex(
1929
1930 (boost::format("%s(.*)") % name).str()
1931
1932 ))
1933
1934 ) {
1935 Range ents;
1936 CHK_MOAB_THROW(m_ptr->getMeshsetIdEntitiesByDimension(
1937 m_field.get_moab(), dim, ents, true),
1938 "by dim");
1939 map[m_ptr->getMeshsetId()] = ents;
1940 MOFEM_LOG("EPSYNC", sev) << "Meshset: " << m_ptr->getMeshsetId() << " "
1941 << ents.size() << " entities";
1942 }
1943
1944 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), sev);
1945 return map;
1946 };
1947
1948 auto get_map_skin = [&](auto &&map) {
1949 ParallelComm *pcomm =
1950 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
1951
1952 Skinner skin(&m_field.get_moab());
1953 for (auto &m : map) {
1954 Range skin_faces;
1955 CHKERR skin.find_skin(0, m.second, false, skin_faces);
1956 CHK_MOAB_THROW(pcomm->filter_pstatus(
1957 skin_faces, PSTATUS_SHARED | PSTATUS_MULTISHARED,
1958 PSTATUS_NOT, -1, nullptr),
1959 "filter");
1960 m.second.swap(skin_faces);
1961 }
1962 return map;
1963 };
1964
1965 /* The above code is written in C++ and it appears to be defining and using
1966 various operations on boundary elements and side elements. */
1968 using BoundaryEleOp = BoundaryEle::UserDataOperator;
1969
1970 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1971
1972 auto calcs_side_traction = [&](auto &pip) {
1974 using EleOnSide =
1976 using SideEleOp = EleOnSide::UserDataOperator;
1977 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1978 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
1979 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1980 boost::make_shared<CGGUserPolynomialBase>(nullptr, true);
1981 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
1982 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
1984 op_loop_domain_side->getOpPtrVector().push_back(
1986 ep.piolaStress, contact_common_data_ptr->contactTractionPtr(),
1987 boost::make_shared<double>(1.0)));
1988 pip.push_back(op_loop_domain_side);
1990 };
1991
1992 auto add_contact_three = [&]() {
1994 auto tree_moab_ptr = boost::make_shared<moab::Core>();
1995 fe_contact_tree = boost::make_shared<ContactTree>(
1996 m_field, tree_moab_ptr, ep.spaceOrder,
1997 get_body_range("CONTACT", SPACE_DIM - 1, Sev::inform));
1998 fe_contact_tree->getOpPtrVector().push_back(
2000 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2001 CHKERR calcs_side_traction(fe_contact_tree->getOpPtrVector());
2002 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2003 fe_contact_tree->getOpPtrVector().push_back(
2005 fe_contact_tree->getOpPtrVector().push_back(
2006 new OpMoveNode(fe_contact_tree, contact_common_data_ptr, u_h1_ptr));
2008 };
2009
2010 CHKERR add_contact_three();
2011
2013 };
2014
2015 CHK_THROW_MESSAGE(impl(), "createContactDetectionFiniteElement");
2016
2017 struct exclude_sdf {
2018 exclude_sdf(Range &&r) : map(r) {}
2019 bool operator()(FEMethod *fe_method_ptr) {
2020 auto ent = fe_method_ptr->getFEEntityHandle();
2021 if (map.find(ent) != map.end()) {
2022 return false;
2023 }
2024 return true;
2025 }
2026
2027 private:
2028 Range map;
2029 };
2030
2031 fe_contact_tree->exeTestHook =
2032 exclude_sdf(get_range_from_block(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2033
2034 return fe_contact_tree;
2035};
2036
2037static auto get_body_range(MoFEM::Interface &m_field, const std::string name,
2038 int dim) {
2039 std::map<int, Range> map;
2040
2041 for (auto m_ptr :
2042 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
2043
2044 (boost::format("%s(.*)") % name).str()
2045
2046 ))
2047
2048 ) {
2049 Range ents;
2050 CHK_MOAB_THROW(m_ptr->getMeshsetIdEntitiesByDimension(m_field.get_moab(),
2051 dim, ents, true),
2052 "by dim");
2053 map[m_ptr->getMeshsetId()] = ents;
2054 }
2055
2056 return map;
2057};
2058
2060 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2061 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2063
2064 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2065
2066 auto &m_field = ep.mField;
2067
2068 using BoundaryEle =
2070 using EleOnSide =
2072 using SideEleOp = EleOnSide::UserDataOperator;
2073 using BdyEleOp = BoundaryEle::UserDataOperator;
2074
2075 // First: Iterate over skeleton FEs adjacent to Domain FEs
2076 // Note: BoundaryEle, i.e. uses skeleton interation rule
2077 auto op_loop_skeleton_side = new OpLoopSide<BoundaryEle>(
2078 m_field, ep.contactElement, SPACE_DIM - 1, Sev::noisy);
2079
2080 auto rule_contact = [](int, int, int o) { return -1; };
2082
2083 auto set_rule_contact = [refine](
2084
2085 ForcesAndSourcesCore *fe_raw_ptr, int order_row,
2086 int order_col, int order_data
2087
2088 ) {
2090 auto rule = 2 * order_data;
2091 fe_raw_ptr->gaussPts = Tools::refineTriangleIntegrationPts(rule, refine);
2093 };
2094
2095 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2096 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2097 CHKERR
2098 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2099 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
2100 ep.frontAdjEdges);
2101
2102 // Second: Iterate over domain FEs adjacent to skelton, particularly
2103 // one domain element.
2104 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2105
2106 // Data storing contact fields
2107 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2108
2109 auto add_ops_domain_side = [&](auto &pip) {
2111 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
2112 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
2113 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
2114 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2115 boost::make_shared<CGGUserPolynomialBase>(nullptr, true);
2116 CHKERR
2117 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
2118 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
2120 op_loop_domain_side->getOpPtrVector().push_back(
2122 broken_data_ptr));
2123 op_loop_domain_side->getOpPtrVector().push_back(
2125 ep.piolaStress, contact_common_data_ptr->contactTractionPtr()));
2126 pip.push_back(op_loop_domain_side);
2128 };
2129
2130 auto add_ops_contact_rhs = [&](auto &pip) {
2132 // get body id and SDF range
2133 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2134 get_body_range(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2135
2136 pip.push_back(new OpCalculateVectorFieldValues<3>(
2137 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2138 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2139 pip.push_back(
2141 pip.push_back(new OpTreeSearch(
2142 contact_tree_ptr, u_h1_ptr,
2143 contact_common_data_ptr->contactTractionPtr(),
2144 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1), nullptr,
2145 nullptr));
2147 ep.contactDisp, contact_common_data_ptr, contact_tree_ptr,
2148 contact_sfd_map_range_ptr));
2150 broken_data_ptr, contact_common_data_ptr, contact_tree_ptr));
2151
2153 };
2154
2155 // push ops to face/side pipeline
2156 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2157 CHKERR add_ops_contact_rhs(op_loop_skeleton_side->getOpPtrVector());
2158
2159 // Add skeleton to domain pipeline
2160 pip.push_back(op_loop_skeleton_side);
2161
2163};
2164
2166 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2167 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2169
2170 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2171 auto &m_field = ep.mField;
2172
2173 using BoundaryEle =
2175 using EleOnSide =
2177 using SideEleOp = EleOnSide::UserDataOperator;
2178 using BdyEleOp = BoundaryEle::UserDataOperator;
2179
2180 // First: Iterate over skeleton FEs adjacent to Domain FEs
2181 // Note: BoundaryEle, i.e. uses skeleton interation rule
2182 auto op_loop_skeleton_side = new OpLoopSide<BoundaryEle>(
2183 m_field, ep.contactElement, SPACE_DIM - 1, Sev::noisy);
2184
2185 auto rule_contact = [](int, int, int o) { return -1; };
2187
2188 auto set_rule_contact = [refine](
2189
2190 ForcesAndSourcesCore *fe_raw_ptr, int order_row,
2191 int order_col, int order_data
2192
2193 ) {
2195 auto rule = 2 * order_data;
2196 fe_raw_ptr->gaussPts = Tools::refineTriangleIntegrationPts(rule, refine);
2198 };
2199
2200 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2201 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2202 CHKERR
2203 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2204 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
2205 ep.frontAdjEdges);
2206
2207 // Second: Iterate over domain FEs adjacent to skelton, particularly
2208 // one domain element.
2209 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2210
2211 // Data storing contact fields
2212 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2213
2214 auto add_ops_domain_side = [&](auto &pip) {
2216 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
2217 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
2218 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
2219 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2220 boost::make_shared<CGGUserPolynomialBase>(nullptr, true);
2221 CHKERR
2222 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
2223 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
2225 op_loop_domain_side->getOpPtrVector().push_back(
2227 broken_data_ptr));
2228 op_loop_domain_side->getOpPtrVector().push_back(
2230 ep.piolaStress, contact_common_data_ptr->contactTractionPtr()));
2231 pip.push_back(op_loop_domain_side);
2233 };
2234
2235 auto add_ops_contact_lhs = [&](auto &pip) {
2237 pip.push_back(new OpCalculateVectorFieldValues<3>(
2238 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2239 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2240 pip.push_back(
2242 pip.push_back(new OpTreeSearch(
2243 contact_tree_ptr, u_h1_ptr,
2244 contact_common_data_ptr->contactTractionPtr(),
2245 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1), nullptr,
2246 nullptr));
2247
2248 // get body id and SDF range
2249 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2250 get_body_range(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2251
2253 ep.contactDisp, ep.contactDisp, contact_common_data_ptr,
2254 contact_tree_ptr, contact_sfd_map_range_ptr));
2255 pip.push_back(
2257 ep.contactDisp, broken_data_ptr, contact_common_data_ptr,
2258 contact_tree_ptr, contact_sfd_map_range_ptr));
2259 pip.push_back(
2261 broken_data_ptr, ep.contactDisp, contact_common_data_ptr,
2262 contact_tree_ptr));
2263
2265 };
2266
2267 // push ops to face/side pipeline
2268 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2269 CHKERR add_ops_contact_lhs(op_loop_skeleton_side->getOpPtrVector());
2270
2271 // Add skeleton to domain pipeline
2272 pip.push_back(op_loop_skeleton_side);
2273
2275};
2276
2279 boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2280 boost::shared_ptr<MatrixDouble> u_h1_ptr,
2281 boost::shared_ptr<MatrixDouble> contact_traction_ptr,
2282 Range r, moab::Interface *post_proc_mesh_ptr,
2283 std::vector<EntityHandle> *map_gauss_pts_ptr) {
2284
2285 auto &m_field = ep.mField;
2286 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2287 return new OpTreeSearch(
2288 contact_tree_ptr, u_h1_ptr, contact_traction_ptr,
2289 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1),
2290 post_proc_mesh_ptr, map_gauss_pts_ptr);
2291}
2292
2293} // namespace EshelbianPlasticity
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)
constexpr double a
static const double eps
constexpr int SPACE_DIM
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Kronecker Delta class.
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()
@ NOSPACE
Definition definitions.h:83
#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
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
static const bool debug
#define N_MBTRI1(x, y)
triangle shape function
Definition fem_tools.h:47
#define N_MBTRI0(x, y)
triangle shape function
Definition fem_tools.h:46
#define N_MBTRI2(x, y)
triangle shape function
Definition fem_tools.h:48
static const double face_coords[4][9]
constexpr auto t_kd
double eta
IntegrationType
Form integrator integration types.
@ GAUSS
Gaussian quadrature integration.
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
auto get_normalize_normals(FTensor::Tensor1< T1, DIM1 > &&t_normal_at_pts, size_t nb_gauss_pts)
double cn_contact
Definition contact.cpp:97
double alpha_contact_quadratic
FormsIntegrators< BoundaryEleOp >::Assembly< A >::OpBase AssemblyBoundaryEleOp
MatrixDouble grad_surface_distance_function(double delta_t, double t, int nb_gauss_pts, MatrixDouble &m_spatial_coords, MatrixDouble &m_normals_at_pts, int block_id)
MatrixDouble hess_surface_distance_function(double delta_t, double t, int nb_gauss_pts, MatrixDouble &m_spatial_coords, MatrixDouble &m_normals_at_pts, int block_id)
PostProcEleByDim< SPACE_DIM >::SideEle SideEle
double alpha_contact_const
double airplane_ray_distance
auto get_spatial_coords(FTensor::Tensor1< T1, DIM1 > &&t_coords, FTensor::Tensor1< T2, DIM2 > &&t_disp, size_t nb_gauss_pts)
double constrain(double sdf, double tn)
constrain function
VectorDouble surface_distance_function(double delta_t, double t, int nb_gauss_pts, MatrixDouble &m_spatial_coords, MatrixDouble &m_normals_at_pts, int block_id)
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
Definition Types.hpp:92
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
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)
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)
constexpr AssemblyType A
constexpr double g
FTensor::Index< 'm', 3 > m
boost::shared_ptr< Range > frontAdjEdges
static double dynamicTime
MoFEM::Interface & mField
const std::string materialH1Positions
static int dynamicStep
const std::string elementVolumeName
static PetscBool dynamicRelaxation
const std::string spatialH1Disp
const std::string piolaStress
const std::string contactDisp
const std::string contactElement
boost::shared_ptr< OrientedBoxTreeTool > & getTreeSurfPtr()
std::map< EntityHandle, std::vector< FaceData > > MapFaceData
boost::shared_ptr< OrientedBoxTreeTool > treeSurfPtr
auto findFaceDataVecPtr(EntityHandle fe_ent)
auto getFaceDataPtr(std::vector< FaceData >::iterator &it, int gg, std::vector< FaceData > *vec_ptr)
MoFEMErrorCode buildTree(Range &ents)
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
MoFEMErrorCode postProcess()
Post-processing function executed at loop completion.
ContactTree(MoFEM::Interface &m_field, boost::shared_ptr< moab::Core > core_mesh_ptr, int max_order, std::map< int, Range > &&body_map)
int getMaxLevel() const
Determine refinement level based on fields approx ordre.
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< 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
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::weak_ptr< ContactTree > contactTreePtr
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)
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.
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 postProcess()
Post-processing function executed at loop completion.
static RefineTrianglesReturn refineTriangle(int nb_levels)
create uniform triangle mesh of refined elements
Definition Tools.cpp:724
static MatrixDouble refineTriangleIntegrationPts(MatrixDouble pts, RefineTrianglesReturn refined)
generate integration points for refined triangle mesh for last level
Definition Tools.cpp:791
static MoFEMErrorCode getLocalCoordinatesOnReferenceThreeNodeTri(const double *elem_coords, const double *glob_coords, const int nb_nodes, double *local_coords)
Get the local coordinates on reference three node tri object.
Definition Tools.cpp:188
static MoFEMErrorCode getTriNormal(const double *coords, double *normal, double *d_normal=nullptr)
Get the Tri Normal objectGet triangle normal.
Definition Tools.cpp:353
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
BoundaryEle::UserDataOperator BdyEleOp
double zeta
Viscous hardening.
Definition plastic.cpp:130