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);
104 inline boost::shared_ptr<OrientedBoxTreeTool> &getTreeSurfPtr() {
105 return treeSurfPtr;
106 }
107 inline auto getRootSetSurf() { return rootSetSurf; }
108 int getMaxLevel() const { return refElementsMap.at(MBTRI)->defMaxLevel; }
109
110 friend struct OpMoveNode;
111 friend struct OpTreeSearch;
112
113 struct FaceData {
114 int gaussPtNb; //< integration points number
115 std::array<double, 3> slavePoint;
116 std::array<double, 3> masterPoint;
117 std::array<double, 3> rayPoint;
118 std::array<double, 3> unitRay;
119 double eleRadius;
120
121 // std::vector<int> dofsSlaveIds;
122 // std::vector<double> dofsSlaveCoeff;
123 // std::vector<double> baseSlaveFuncs;
124
125 std::array<double, 9> masterPointNodes;
126 std::array<double, 9> masterTractionNodes;
127 std::array<double, 9> slavePointNodes;
128 std::array<double, 9> slaveTractionNodes;
129
130 FaceData() = default;
131 };
132
133 using MapFaceData = std::map<EntityHandle, std::vector<FaceData>>;
134
135 inline auto findFaceDataVecPtr(EntityHandle fe_ent) {
136 auto &map_face_data = shadowDataMap;
137 auto it = map_face_data.find(fe_ent);
138 if (it == map_face_data.end()) {
139 return (std::vector<FaceData> *)nullptr;
140 }
141 return &(it->second);
142 }
143
144 inline auto getFaceDataPtr(std::vector<FaceData>::iterator &it, int gg,
145 std::vector<FaceData> *vec_ptr) {
146 FaceData *face_data_ptr = nullptr;
147 if (it != vec_ptr->end()) {
148 if (it->gaussPtNb == gg) {
149 face_data_ptr = &(*it);
150 ++it;
151 }
152 }
153 return face_data_ptr;
154 }
155
156protected:
158 boost::shared_ptr<OrientedBoxTreeTool> treeSurfPtr;
160
162 // Tag thCoeff;
163 // Tag thIds;
164 // Tag thBases;
169
171
172 std::map<int, Range> bodyMap;
173
174 const int maxOrder;
175};
176
177auto checkSdf(EntityHandle fe_ent, std::map<int, Range> &sdf_map_range) {
178 for (auto &m_sdf : sdf_map_range) {
179 if (m_sdf.second.find(fe_ent) != m_sdf.second.end())
180 return m_sdf.first;
181 }
182 return -1;
183}
184
185template <typename OP_PTR>
186auto getSdf(OP_PTR op_ptr, MatrixDouble &contact_disp, int block_id,
187 bool eval_hessian) {
188
189 auto nb_gauss_pts = op_ptr->getGaussPts().size2();
190
191 auto ts_time = op_ptr->getTStime();
192 auto ts_time_step = op_ptr->getTStimeStep();
195 ts_time_step = EshelbianCore::dynamicStep;
196 }
197
198 auto m_spatial_coords = ContactOps::get_spatial_coords(
199 op_ptr->getFTensor1CoordsAtGaussPts(),
200 getFTensor1FromMat<3>(contact_disp), nb_gauss_pts);
201 auto m_normals_at_pts = ContactOps::get_normalize_normals(
202 op_ptr->getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
204 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
205 block_id);
207 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
208 block_id);
209
210#ifndef NDEBUG
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)
218 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Should be size of 3");
219#endif // NDEBUG
220
221 if (eval_hessian) {
223 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
224 block_id);
225 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
226 m_hess_sdf);
227 } else {
228
229 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
230 MatrixDouble(6, nb_gauss_pts, 0.));
231 }
232};
233
234/**
235 * @brief Calculate points data on contact surfaces
236 *
237 * @tparam T1
238 * @param unit_ray
239 * @param point
240 * @param elem_point_nodes
241 * @param elem_traction_nodes
242 * @param t_spatial_coords
243 * @return auto
244 */
245template <typename T1>
247
248 std::array<double, 3> &unit_ray, std::array<double, 3> &point,
249
250 std::array<double, 9> &elem_point_nodes,
251 std::array<double, 9> &elem_traction_nodes,
252
253 FTensor::Tensor1<T1, 3> &t_spatial_coords) {
254 FTensor::Index<'i', 3> i;
255
256 auto t_unit_ray = getFTensor1FromPtr<3>(unit_ray.data());
257 auto t_point = getFTensor1FromPtr<3>(point.data());
258
259 auto get_normal = [](auto &ele_coords) {
261 Tools::getTriNormal(ele_coords.data(), &t_normal(0));
262 return t_normal;
263 };
264
265 auto t_normal = get_normal(elem_point_nodes);
266 t_normal(i) /= t_normal.l2();
267
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);
271
272 auto gamma = (sn - nm) / nr;
273
274 FTensor::Tensor1<T1, 3> t_point_current;
275 t_point_current(i) = t_spatial_coords(i) + gamma * t_unit_ray(i);
276
277 auto get_local_point_shape_functions = [&](auto &&t_elem_coords,
278 auto &t_point) {
279 std::array<T1, 2> loc_coords;
282 &t_elem_coords(0, 0), &t_point(0), 1, loc_coords.data()),
283 "get local coords");
284 return FTensor::Tensor1<T1, 3>{N_MBTRI0(loc_coords[0], loc_coords[1]),
285 N_MBTRI1(loc_coords[0], loc_coords[1]),
286 N_MBTRI2(loc_coords[0], loc_coords[1])};
287 };
288
289 auto eval_position = [&](auto &&t_field, auto &t_shape_fun) {
290 FTensor::Index<'i', 3> i;
291 FTensor::Index<'j', 3> j;
292 FTensor::Tensor1<T1, 3> t_point_out;
293 t_point_out(i) = t_shape_fun(j) * t_field(j, i);
294 return t_point_out;
295 };
296
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);
303
304 return std::make_tuple(t_slave_point_updated, t_traction_updated, t_normal);
305};
306
307template <typename T1>
309
310 ContactTree::FaceData *face_data_ptr,
311 FTensor::Tensor1<T1, 3> &t_spatial_coords
312
313) {
314
315 return multiPoint(face_data_ptr->unitRay, face_data_ptr->masterPoint,
316 face_data_ptr->masterPointNodes,
317 face_data_ptr->masterTractionNodes, t_spatial_coords);
318};
319
320template <typename T1>
322
323 ContactTree::FaceData *face_data_ptr,
324 FTensor::Tensor1<T1, 3> &t_spatial_coords
325
326) {
327
328 return multiPoint(face_data_ptr->unitRay, face_data_ptr->slavePoint,
329 face_data_ptr->slavePointNodes,
330 face_data_ptr->slaveTractionNodes, t_spatial_coords);
331};
332
333/**
334 * Evaluate gap and tractions between master and slave points
335 */
336template <typename T1>
338 FTensor::Tensor1<T1, 3> &t_spatial_coords) {
339
340 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
341 multiSlavePoint(face_data_ptr, t_spatial_coords);
342 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
343 multiMasterPoint(face_data_ptr, t_spatial_coords);
344
345 FTensor::Index<'i', 3> i;
346
348 t_normal(i) = t_master_normal(i) - t_slave_normal(i);
349 t_normal.normalize();
350
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);
355
356 return std::make_tuple(gap, tn_master, tn_slave,
357 ContactOps::constrain(gap, tn),
358 t_master_traction_current, t_slave_traction_current);
359};
360
362
363/** Caluclate rhs term for contact
364 */
365template <typename T1, typename T2, typename T3>
367
368 ContactTree::FaceData *face_data_ptr, FTensor::Tensor1<T1, 3> &t_coords,
369 FTensor::Tensor1<T2, 3> &t_spatial_coords,
370 FTensor::Tensor1<T3, 3> &t_master_traction, MultiPointRhsType type,
371 bool debug = false
372
373) {
374 FTensor::Index<'i', 3> i;
375 FTensor::Index<'j', 3> j;
376
378 t_u(i) = t_spatial_coords(i) - t_coords(i);
379
381
382 switch (type) {
384
385 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
386 multiSlavePoint(face_data_ptr, t_spatial_coords);
387 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
388 multiMasterPoint(face_data_ptr, t_spatial_coords);
389
390 // average normal surface
392 t_normal(i) = t_master_normal(i) - t_slave_normal(i);
393 t_normal.normalize();
394
395 // get projection operators
397 t_P(i, j) = t_normal(i) * t_normal(j);
399 t_Q(i, j) = kronecker_delta(i, j) - t_P(i, j);
400
401 constexpr double beta = 0.5;
402
403 auto zeta = 1e-8 * face_data_ptr->eleRadius;
406
407 // this is regularised std::min(d,s)
408 auto f_min_gap = [zeta](auto d, auto s) {
409 return 0.5 * (d + s - std::sqrt((d - s) * (d - s) + zeta));
410 };
411 // this derivative of regularised std::min(d,s)
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));
414 };
415
416 // add penalty on penetration side
417 auto f_barrier = [alpha1, alpha2, f_min_gap, f_diff_min_gap](auto g,
418 auto tn) {
419 auto d = alpha1 * g;
420 auto b1 =
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;
423 return b1 - b2;
424 };
425
426 FTensor::Tensor1<T2, 3> t_gap_vec;
427 t_gap_vec(i) = beta * t_spatial_coords(i) +
428 (1 - beta) * t_slave_point_current(i) - t_spatial_coords(i);
430 t_traction_vec(i) =
431 -beta * t_master_traction(i) + (beta - 1) * t_slave_traction_current(i);
432
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);
436
438 // add penalty on penetration side
439 t_traction_bar(i) = t_normal(i) * barrier;
440
441 t_rhs(i) =
442
443 t_Q(i, j) * t_master_traction(j) +
444
445 t_P(i, j) * (t_master_traction(j) - t_traction_bar(j));
446
447 if (debug) {
448 auto is_nan_or_inf = [](double value) -> bool {
449 return std::isnan(value) || std::isinf(value);
450 };
451
452 double v = std::complex<double>(t_rhs(i) * t_rhs(i)).real();
453 if (is_nan_or_inf(v)) {
454 MOFEM_LOG_CHANNEL("SELF");
455 MOFEM_LOG("SELF", Sev::error) << "t_rhs " << t_rhs;
456 CHK_MOAB_THROW(MOFEM_DATA_INCONSISTENCY, "rhs is nan or inf");
457 }
458 }
459
460 } break;
463 break;
464 }
465
466 return t_rhs;
467};
468
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);
475
477
478private:
479 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
480 boost::shared_ptr<ContactTree> contactTreePtr;
481 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
482};
483
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)
489 : ContactOps::AssemblyBoundaryEleOp(row_field_name, row_field_name,
490 ContactOps::BoundaryEleOp::OPROW),
491 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
492 sdfMapRangePtr(sdf_map_range_ptr) {
493 CHK_THROW_MESSAGE(PetscOptionsGetScalar(PETSC_NULLPTR, "", "-cn",
495 PETSC_NULLPTR),
496 "get cn failed");
499 PetscOptionsGetScalar(PETSC_NULLPTR, "", "-alpha_contact_const",
500 &ContactOps::alpha_contact_const, PETSC_NULLPTR),
501 "get alpha contact failed");
503 PETSC_NULLPTR, "", "-alpha_contact_quadratic",
505 "get alpha contact failed");
507 PetscOptionsGetScalar(PETSC_NULLPTR, "", "-airplane_ray_distance",
508 &ContactOps::airplane_ray_distance, PETSC_NULLPTR),
509 "get alpha contact failed");
510
511 MOFEM_LOG("EP", Sev::inform) << "cn " << ContactOps::cn_contact;
512 MOFEM_LOG("EP", Sev::inform)
513 << "alpha_contact_const " << ContactOps::alpha_contact_const;
514 MOFEM_LOG("EP", Sev::inform)
515 << "alpha_contact_quadratic " << ContactOps::alpha_contact_quadratic;
516 MOFEM_LOG("EP", Sev::inform)
517 << "airplane_ray_distance " << ContactOps::airplane_ray_distance;
518}
519
523
524 FTensor::Index<'i', 3> i;
525 FTensor::Index<'j', 3> j;
526 FTensor::Index<'k', 3> k;
527 FTensor::Index<'l', 3> l;
528
529 const size_t nb_gauss_pts = getGaussPts().size2();
530
531#ifndef NDEBUG
532 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
533 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
534 "Wrong number of integration pts %ld != %ld",
535 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
536 }
537#endif // !NDEBUG
538
539 auto &nf = locF;
540 locF.clear();
541
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);
546
547 // placeholder to pass boundary block id to python. default SDF is set on
548 // block = -1, one can choose different block by making block "CONTACT_SDF",
549 // then specific SDF can be set to that block.
550 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
551 getSdf(this, commonDataPtr->contactDisp,
552 checkSdf(getFEEntityHandle(), *sdfMapRangePtr), false);
553
554 auto t_sdf_v = getFTensor0FromVec(v_sdf);
555 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
556 auto t_normalize_normal = getFTensor1FromMat<3>(m_normals_at_pts);
557
558 auto next = [&]() {
559 ++t_w;
560 ++t_coords;
561 ++t_disp;
562 ++t_traction;
563 ++t_normalize_normal;
564 ++t_sdf_v;
565 ++t_grad_sdf_v;
566 };
567
568 auto face_data_vec_ptr =
569 contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
570 auto face_gauss_pts_it = face_data_vec_ptr->begin();
571
572 auto nb_base_functions = data.getN().size2();
573 auto t_base = data.getFTensor0N();
574 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
575
577 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
578 face_data_vec_ptr);
579
580 auto check_face_contact = [&]() {
581 if (checkSdf(getFEEntityHandle(), *sdfMapRangePtr) != -1)
582 return false;
583
584 if (face_data_ptr) {
585 return true;
586 }
587 return false;
588 };
589
590#ifdef ENABLE_PYTHON_BINDING
591 double c = 0.;
592 if (ContactOps::sdfPythonWeakPtr.lock()) {
593 auto tn = t_traction(i) * t_grad_sdf_v(i);
594 c = ContactOps::constrain(t_sdf_v, tn);
595 }
596#else
597 constexpr double c = 0;
598#endif
599
600 if (!c && check_face_contact()) {
601 FTensor::Tensor1<double, 3> t_spatial_coords;
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,
604 t_traction, MultiPointRhsType::U, true);
605 t_rhs(i) = t_rhs_tmp(i);
606
607 } else {
608
609#ifdef ENABLE_PYTHON_BINDING
610 auto inv_cn = 1. / ContactOps::cn_contact;
611
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);
618 } else {
619 t_rhs(i) = t_traction(i);
620 }
621#else
622 t_rhs(i) = t_traction(i);
623#endif
624 }
625
626 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
627 const double alpha = t_w * getMeasure();
628
629 size_t bb = 0;
630 for (; bb != nbRows / 3; ++bb) {
631 const double beta = alpha * t_base;
632 t_nf(i) -= beta * t_rhs(i);
633 ++t_nf;
634 ++t_base;
635 }
636 for (; bb < nb_base_functions; ++bb)
637 ++t_base;
638
639 next();
640 }
641
643}
644
645template <AssemblyType A, IntegrationType I> struct OpConstrainBoundaryHDivRhs;
646
647template <AssemblyType A>
649 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
650
653
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);
658
660
661private:
662 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
663 boost::shared_ptr<ContactTree> contactTreePtr;
664};
665
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) {}
675
676template <AssemblyType A>
680
685
686 const size_t nb_gauss_pts = OP::getGaussPts().size2();
687
688#ifndef NDEBUG
689 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
690 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
691 "Wrong number of integration pts %ld != %ld",
692 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
693 }
694#endif // !NDEBUG
695
696 auto &nf = OP::locF;
697 OP::locF.clear();
698
699 auto t_w = OP::getFTensor0IntegrationWeight();
700 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
701 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
702
703 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
704 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
705
706 auto next = [&]() {
707 ++t_w;
708 ++t_disp;
709 ++t_traction;
710 ++t_coords;
711 ++t_material_normal;
712 };
713
714 auto face_data_vec_ptr =
715 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
716 auto face_gauss_pts_it = face_data_vec_ptr->begin();
717
718 auto nb_base_functions = data.getN().size2() / 3;
719 auto t_base = data.getFTensor1N<3>();
720 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
721
722 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
723 const double alpha = t_w / 2.;
724
725 size_t bb = 0;
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);
729 ++t_nf;
730 ++t_base;
731 }
732 for (; bb < nb_base_functions; ++bb)
733 ++t_base;
734
735 next();
736 }
737
739}
740
742
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);
748
751
752private:
753 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
754 boost::shared_ptr<ContactTree> contactTreePtr;
755 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
756};
757
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)
763 : ContactOps::AssemblyBoundaryEleOp(row_field_name, col_field_name,
764 ContactOps::BoundaryEleOp::OPROWCOL),
765 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
766 sdfMapRangePtr(sdf_map_range_ptr) {
767
768 sYmm = false;
769}
770
773 EntitiesFieldData::EntData &col_data) {
775
776 using namespace ContactOps;
777
778 FTensor::Index<'i', 3> i;
779 FTensor::Index<'j', 3> j;
780 FTensor::Index<'k', 3> k;
781
782 auto nb_rows = row_data.getIndices().size();
783 auto nb_cols = col_data.getIndices().size();
784
785 auto &locMat = AssemblyBoundaryEleOp::locMat;
786 locMat.resize(nb_rows, nb_cols, false);
787 locMat.clear();
788
789 if (nb_cols && nb_rows) {
790
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);
796
797 // placeholder to pass boundary block id to python
798 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
799 getSdf(this, commonDataPtr->contactDisp,
800 checkSdf(getFEEntityHandle(), *sdfMapRangePtr), true);
801
802 auto t_sdf_v = getFTensor0FromVec(v_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);
806
807 auto next = [&]() {
808 ++t_w;
809 ++t_coords;
810 ++t_disp;
811 ++t_traction;
812 ++t_sdf_v;
813 ++t_grad_sdf_v;
814 ++t_hess_sdf_v;
815 ++t_normalized_normal;
816 };
817
818 auto face_data_vec_ptr =
819 contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
820 auto face_gauss_pts_it = face_data_vec_ptr->begin();
821
822 auto t_row_base = row_data.getFTensor0N();
823 auto nb_face_functions = row_data.getN().size2() / 3;
824 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
825
826 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
827
828 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
829 face_data_vec_ptr);
830
831 auto check_face_contact = [&]() {
832 if (checkSdf(getFEEntityHandle(), *sdfMapRangePtr) != -1)
833 return false;
834
835 if (face_data_ptr) {
836 return true;
837 }
838 return false;
839 };
840
842
843#ifdef ENABLE_PYTHON_BINDING
844 double c = 0.;
845 if (ContactOps::sdfPythonWeakPtr.lock()) {
846 auto tn = t_traction(i) * t_grad_sdf_v(i);
847 c = ContactOps::constrain(t_sdf_v, tn);
848 }
849#else
850 constexpr double c = 0;
851#endif
852
853 if (!c && check_face_contact()) {
854 FTensor::Tensor1<double, 3> t_spatial_coords;
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) {
858 FTensor::Tensor1<std::complex<double>, 3> t_spatial_coords_cx{
859 t_spatial_coords(0), t_spatial_coords(1), t_spatial_coords(2)};
860 t_spatial_coords_cx(ii) += eps * 1i;
861 auto t_rhs_tmp =
862 multiPointRhs(face_data_ptr, t_coords, t_spatial_coords_cx,
863 t_traction, MultiPointRhsType::U);
864 for (int jj = 0; jj != 3; ++jj) {
865 auto v = t_rhs_tmp(jj).imag();
866 t_res_dU(jj, ii) = v / eps;
867 }
868 }
869
870 } else {
871
872#ifdef ENABLE_PYTHON_BINDING
873
874 if (ContactOps::sdfPythonWeakPtr.lock()) {
875 auto inv_cn = 1. / ContactOps::cn_contact;
876 t_res_dU(i, j) =
877
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))
880
881 + (c * inv_cn) * (t_sdf_v * t_hess_sdf_v(i, j) +
882
883 t_grad_sdf_v(j) * t_grad_sdf_v(i));
884 } else {
885 t_res_dU(i, j) = 0;
886 }
887#else
888 t_res_dU(i, j) = 0;
889#endif
890 }
891
892 auto alpha = t_w * getMeasure();
893
894 size_t rr = 0;
895 for (; rr != nb_rows / 3; ++rr) {
896
897 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
898 auto t_col_base = col_data.getFTensor0N(gg, 0);
899
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);
903 ++t_col_base;
904 ++t_mat;
905 }
906
907 ++t_row_base;
908 }
909 for (; rr < nb_face_functions; ++rr)
910 ++t_row_base;
911
912 next();
913 }
914 }
915
917}
918
919template <AssemblyType A, IntegrationType I> struct OpConstrainBoundaryL2Lhs_dP;
920
921template <AssemblyType A>
923 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
924
927
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);
934
937
938private:
939 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
940 boost::shared_ptr<ContactTree> contactTreePtr;
941 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
942};
943
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) {
956 OP::sYmm = false;
957}
958
959template <AssemblyType A>
963 EntitiesFieldData::EntData &col_data) {
965
966 using namespace ContactOps;
967
968 FTensor::Index<'i', 3> i;
969 FTensor::Index<'j', 3> j;
970 FTensor::Index<'k', 3> k;
971
972 auto nb_rows = row_data.getIndices().size();
973 auto nb_cols = col_data.getIndices().size();
974
975 auto &locMat = AssemblyBoundaryEleOp::locMat;
976 locMat.resize(nb_rows, nb_cols, false);
977 locMat.clear();
978
979 if (nb_cols && nb_rows) {
980
981 const size_t nb_gauss_pts = OP::getGaussPts().size2();
982
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();
988
989 // placeholder to pass boundary block id to python
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);
993
994 auto t_sdf_v = getFTensor0FromVec(v_sdf);
995 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
996
997 auto next = [&]() {
998 ++t_w;
999 ++t_disp;
1000 ++t_traction;
1001 ++t_coords;
1002 ++t_material_normal;
1003 ++t_sdf_v;
1004 ++t_grad_sdf_v;
1005 };
1006
1007 auto face_data_vec_ptr =
1008 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1009 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1010
1011 auto t_row_base = row_data.getFTensor0N();
1012 auto nb_face_functions = row_data.getN().size2();
1013
1014 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
1015
1016 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1017
1018 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
1019 face_data_vec_ptr);
1020
1021 auto check_face_contact = [&]() {
1022 if (checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr) != -1)
1023 return false;
1024
1025 if (face_data_ptr) {
1026 return true;
1027 }
1028 return false;
1029 };
1030
1032
1033#ifdef ENABLE_PYTHON_BINDING
1034 double c = 0.;
1035 if (ContactOps::sdfPythonWeakPtr.lock()) {
1036 auto tn = t_traction(i) * t_grad_sdf_v(i);
1037 c = ContactOps::constrain(t_sdf_v, tn);
1038 }
1039#else
1040 constexpr double c = 0;
1041#endif
1042
1043 if (!c && check_face_contact()) {
1044 FTensor::Tensor1<double, 3> t_spatial_coords;
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) {
1048 FTensor::Tensor1<std::complex<double>, 3> t_traction_cx{
1049 t_traction(0), t_traction(1), t_traction(2)};
1050 t_traction_cx(ii) += eps * 1i;
1051 auto t_rhs_tmp =
1052 multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
1053 t_traction_cx, MultiPointRhsType::U);
1054 for (int jj = 0; jj != 3; ++jj) {
1055 auto v = t_rhs_tmp(jj).imag();
1056 t_res_dP(jj, ii) = v / eps;
1057 }
1058 }
1059 } else {
1060
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);
1067 } else {
1068 t_res_dP(i, j) = t_kd(i, j);
1069 }
1070#else
1071 t_res_dP(i, j) = t_kd(i, j);
1072#endif
1073 }
1074
1075 const double alpha = t_w / 2.;
1076 size_t rr = 0;
1077 for (; rr != nb_rows / 3; ++rr) {
1078
1079 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1080 auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
1081
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);
1086 ++t_col_base;
1087 ++t_mat;
1088 }
1089
1090 ++t_row_base;
1091 }
1092 for (; rr < nb_face_functions; ++rr)
1093 ++t_row_base;
1094
1095 next();
1096 }
1097 }
1098
1100}
1101
1102template <AssemblyType A, IntegrationType I>
1104
1105template <AssemblyType A>
1107 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
1108
1111
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);
1117
1119 EntitiesFieldData::EntData &col_data);
1120
1121private:
1122 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
1123 boost::shared_ptr<ContactTree> contactTreePtr;
1124};
1125
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) {
1136 OP::sYmm = false;
1137}
1138
1139template <AssemblyType A>
1143 EntitiesFieldData::EntData &row_data) {
1145
1146 // Note: col_data and row_data are swapped in this function, we going to
1147 // transpose locMat at the end
1148
1149 using namespace ContactOps;
1150
1154
1155 auto nb_rows = row_data.getIndices().size();
1156 auto nb_cols = col_data.getIndices().size();
1157
1158 auto &locMat = AssemblyBoundaryEleOp::locMat;
1159 locMat.resize(nb_rows, nb_cols, false);
1160 locMat.clear();
1161
1162 if (nb_cols && nb_rows) {
1163
1164 const size_t nb_gauss_pts = OP::getGaussPts().size2();
1165
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();
1171
1172 auto next = [&]() {
1173 ++t_w;
1174 ++t_disp;
1175 ++t_traction;
1176 ++t_coords;
1177 ++t_material_normal;
1178 };
1179
1180 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
1181
1182 auto face_data_vec_ptr =
1183 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1184 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1185
1186 auto t_row_base = row_data.getFTensor1N<3>();
1187 auto nb_face_functions = row_data.getN().size2() / 3;
1188 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1189
1190 const auto alpha = t_w / 2.;
1191
1192 size_t rr = 0;
1193 for (; rr != nb_rows / 3; ++rr) {
1194
1195 auto row_base = alpha * (t_row_base(i) * t_material_normal(i));
1196
1197 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1198 auto t_col_base = col_data.getFTensor0N(gg, 0);
1199
1200 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
1201 const auto beta = row_base * t_col_base;
1202 t_mat(i, j) += beta * t_kd(i, j);
1203 ++t_col_base;
1204 ++t_mat;
1205 }
1206
1207 ++t_row_base;
1208 }
1209 for (; rr < nb_face_functions; ++rr)
1210 ++t_row_base;
1211
1212 next();
1213 }
1214 }
1215
1216 locMat = trans(locMat);
1217
1219}
1220
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),
1225 bodyMap(body_map) {
1226
1227 auto ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
1228 ref_ele_ptr->hoNodes =
1229 PETSC_FALSE; ///< So far only linear geometry is implemented for contact
1230
1231 CHK_THROW_MESSAGE(ref_ele_ptr->getOptions(optionsPrefix), "getOptions");
1232 CHK_THROW_MESSAGE(ref_ele_ptr->generateReferenceElementMesh(),
1233 "Error when generating reference element");
1234
1235 MOFEM_LOG("EP", Sev::inform) << "Contact hoNodes " << ref_ele_ptr->hoNodes
1236 ? "true"
1237 : "false";
1238 MOFEM_LOG("EP", Sev::inform)
1239 << "Contact maxOrder " << ref_ele_ptr->defMaxLevel;
1240
1241 refElementsMap[MBTRI] = ref_ele_ptr;
1242
1243 int def_ele_id = -1;
1244 CHKERR getPostProcMesh().tag_get_handle("ELE_ID", 1, MB_TYPE_INTEGER, thEleId,
1245 MB_TAG_CREAT | MB_TAG_DENSE,
1246 &def_ele_id);
1247 CHKERR getPostProcMesh().tag_get_handle("BODY_ID", 1, MB_TYPE_INTEGER,
1248 thBodyId, MB_TAG_CREAT | MB_TAG_DENSE,
1249 &def_ele_id);
1250
1251 // std::vector<int> def_ids(3 * nb_bases, 0);
1252 // CHKERR getPostProcMesh().tag_get_handle("IDS", 3 * nb_bases,
1253 // MB_TYPE_INTEGER,
1254 // thIds, MB_TAG_CREAT | MB_TAG_DENSE,
1255 // &*def_ids.begin());
1256 // std::vector<double> def_coeffs(3 * nb_bases, 0);
1257 // CHKERR getPostProcMesh().tag_get_handle("COEFF", 3 * nb_bases,
1258 // MB_TYPE_DOUBLE,
1259 // thCoeff, MB_TAG_CREAT |
1260 // MB_TAG_DENSE,
1261 // &*def_coeffs.begin());
1262 // std::vector<double> def_basses(nb_bases, 0);
1263 // CHKERR getPostProcMesh().tag_get_handle("BASES", nb_bases, MB_TYPE_DOUBLE,
1264 // thBases, MB_TAG_CREAT |
1265 // MB_TAG_DENSE,
1266 // &*def_basses.begin());
1267
1268 std::array<double, 3> def_small_x{0., 0., 0.};
1269 CHKERR getPostProcMesh().tag_get_handle("x", 3, MB_TYPE_DOUBLE, thSmallX,
1270 MB_TAG_CREAT | MB_TAG_DENSE,
1271 def_small_x.data());
1272 CHKERR getPostProcMesh().tag_get_handle("X", 3, MB_TYPE_DOUBLE, thLargeX,
1273 MB_TAG_CREAT | MB_TAG_DENSE,
1274 def_small_x.data());
1275
1276 std::array<double, 3> def_tractions{0., 0., 0.};
1277 CHKERR getPostProcMesh().tag_get_handle(
1278 "TRACTION", 3, MB_TYPE_DOUBLE, thTraction, MB_TAG_CREAT | MB_TAG_DENSE,
1279 &*def_tractions.begin());
1280}
1281
1287
1290
1292 shadowDataMap.clear();
1293
1294 PetscBarrier(nullptr);
1295
1296 ParallelComm *pcomm_post_proc_mesh =
1297 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
1298 if (pcomm_post_proc_mesh == nullptr)
1299 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
1300
1301 // MOFEM_LOG("EP", Sev::verbose) << "ContactTree broadcast entities";
1302
1303 auto brodacts = [&](auto &brodacts_ents) {
1305 Range recived_ents;
1306 for (auto rr = 0; rr != mField.get_comm_size(); ++rr) {
1307 pcomm_post_proc_mesh->broadcast_entities(
1308 rr, rr == mField.get_comm_rank() ? brodacts_ents : recived_ents,
1309 false, true);
1310 }
1312 };
1313
1314 Range brodacts_ents;
1315 CHKERR getPostProcMesh().get_entities_by_dimension(0, 2, brodacts_ents, true);
1316 CHKERR brodacts(brodacts_ents);
1317
1318 // MOFEM_LOG("EP", Sev::verbose) << "ContactTree build tree";
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
1335 treeSurfPtr = boost::shared_ptr<OrientedBoxTreeTool>(
1336 new OrientedBoxTreeTool(&getPostProcMesh(), "ROOTSETSURF", true));
1337 CHKERR treeSurfPtr->build(ents, rootSetSurf);
1338
1339 // CHKERR treeSurfPtr->stats(rootSetSurf, std::cerr);
1340
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::shared_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 auto get_body_id = [&](auto fe_ent) {
1372 for (auto &m : contactTreePtr->bodyMap) {
1373 if (m.second.find(fe_ent) != m.second.end()) {
1374 return m.first;
1375 }
1376 }
1377 return -1;
1378 };
1379
1380 auto &moab_post_proc_mesh = contactTreePtr->getPostProcMesh();
1381 auto &post_proc_ents = contactTreePtr->getPostProcElements();
1382
1383 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1384 auto fe_id = id_from_handle(fe_ent);
1385 auto body_id = get_body_id(fe_ent);
1386 auto &map_gauss_pts = contactTreePtr->getMapGaussPts();
1387
1388 CHKERR moab_post_proc_mesh.tag_clear_data(contactTreePtr->thEleId,
1389 post_proc_ents, &fe_id);
1390 CHKERR moab_post_proc_mesh.tag_clear_data(contactTreePtr->thBodyId,
1391 post_proc_ents, &body_id);
1392
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();
1397
1398 MatrixDouble x_h1(nb_gauss_pts, 3);
1399 auto t_x_h1 = getFTensor1FromPtr<3>(&x_h1(0, 0));
1400 MatrixDouble x_l2(nb_gauss_pts, 3);
1401 auto t_x_l2 = getFTensor1FromPtr<3>(&x_l2(0, 0));
1402 MatrixDouble tractions = trans(commonDataPtr->contactTraction);
1403 MatrixDouble coords = getCoordsAtGaussPts();
1404
1405 FTensor::Index<'i', 3> i;
1406
1407 // VectorDouble bases(nb_bases, 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);
1411
1412 ++t_coords;
1413 ++t_u_h1;
1414 ++t_u_l2;
1415 ++t_x_h1;
1416 ++t_x_l2;
1417 }
1418
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());
1430
1432}
1433
1435
1437
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,
1441
1442 moab::Interface *post_proc_mesh_ptr,
1443 std::vector<EntityHandle> *map_gauss_pts_ptr
1444
1445 );
1446 MoFEMErrorCode doWork(int side, EntityType type,
1448
1449protected:
1450 boost::shared_ptr<ContactTree> contactTreePtr;
1451 boost::shared_ptr<MatrixDouble> uH1Ptr;
1452 boost::shared_ptr<MatrixDouble> tractionPtr;
1453 moab::Interface *postProcMeshPtr = nullptr;
1454 std::vector<EntityHandle> *mapGaussPtsPtr = nullptr;
1455
1457};
1458
1459OpTreeSearch::OpTreeSearch(boost::shared_ptr<ContactTree> contact_tree_ptr,
1460 boost::shared_ptr<MatrixDouble> u_h1_ptr,
1461 boost::shared_ptr<MatrixDouble> traction_ptr,
1462 Range r,
1463
1464 moab::Interface *post_proc_mesh_ptr,
1465 std::vector<EntityHandle> *map_gauss_pts_ptr
1466
1467 )
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),
1471 contactRange(r) {}
1472
1473MoFEMErrorCode OpTreeSearch::doWork(int side, EntityType type,
1476
1477 auto &m_field = getPtrFE()->mField;
1478 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1479 auto fe_id = id_from_handle(fe_ent);
1480
1481 if (contactRange.find(fe_ent) == contactRange.end())
1483
1484 FTensor::Index<'i', 3> i;
1485 FTensor::Index<'j', 3> j;
1486
1487 const auto nb_gauss_pts = getGaussPts().size2();
1488
1489 auto t_disp_h1 = getFTensor1FromMat<3>(*uH1Ptr);
1490 auto t_coords = getFTensor1CoordsAtGaussPts();
1491 auto t_traction = getFTensor1FromMat<3>(*tractionPtr);
1492
1493 auto next = [&]() {
1494 ++t_disp_h1;
1495 ++t_traction;
1496 ++t_coords;
1497 };
1498
1499 auto get_ele_centre = [i](auto t_ele_coords) {
1500 FTensor::Tensor1<double, 3> t_ele_center;
1501 t_ele_center(i) = 0;
1502 for (int nn = 0; nn != 3; nn++) {
1503 t_ele_center(i) += t_ele_coords(i);
1504 ++t_ele_coords;
1505 }
1506 t_ele_center(i) /= 3;
1507 return t_ele_center;
1508 };
1509
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);
1513 return t_n0.l2();
1514 };
1515
1516 auto get_face_conn = [this](auto face) {
1517 const EntityHandle *conn;
1518 int num_nodes;
1519 CHK_MOAB_THROW(contactTreePtr->getPostProcMesh().get_connectivity(
1520 face, conn, num_nodes, true),
1521 "get conn");
1522 if (num_nodes != 3) {
1523 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "face is not a triangle");
1524 }
1525 return conn;
1526 };
1527
1528 auto get_face_coords = [this](auto conn) {
1529 std::array<double, 9> coords;
1530 CHKERR contactTreePtr->getPostProcMesh().get_coords(conn, 3, coords.data());
1531 return coords;
1532 };
1533
1534 auto get_closet_face = [this](auto *point_ptr, auto r) {
1535 FTensor::Tensor1<double, 3> t_point_out;
1536 std::vector<EntityHandle> faces_out;
1538 contactTreePtr->getTreeSurfPtr()->sphere_intersect_triangles(
1539 point_ptr, r / 8, contactTreePtr->getRootSetSurf(), faces_out),
1540 "get closest faces");
1541 return faces_out;
1542 };
1543
1544 auto get_faces_out = [this](auto *point_ptr, auto *unit_ray_ptr, auto radius,
1545 auto eps) {
1546 std::vector<double> distances_out;
1547 std::vector<EntityHandle> faces_out;
1549
1550 contactTreePtr->getTreeSurfPtr()->ray_intersect_triangles(
1551 distances_out, faces_out, contactTreePtr->getRootSetSurf(), eps,
1552 point_ptr, unit_ray_ptr, &radius),
1553
1554 "get closest faces");
1555 return std::make_pair(faces_out, distances_out);
1556 };
1557
1558 auto get_normal = [](auto &ele_coords) {
1560 Tools::getTriNormal(ele_coords.data(), &t_normal(0));
1561 return t_normal;
1562 };
1563
1564 auto make_map = [&](auto &face_out, auto &face_dist, auto &t_ray_point,
1565 auto &t_unit_ray, auto &t_master_coord) {
1566 FTensor::Index<'i', 3> i;
1567 FTensor::Index<'j', 3> j;
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]);
1571 auto face_coords = get_face_coords(face_conn);
1572 auto t_face_normal = get_normal(face_coords);
1573 t_face_normal.normalize();
1575 t_x(i) = t_ray_point(i) + t_unit_ray(i) * face_dist[ii];
1577 t_chi(i, j) =
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];
1582 }
1583 }
1584 return m;
1585 };
1586
1587 auto get_tag_data = [this](auto tag, auto face, auto &vec) {
1589 int tag_length;
1590 CHKERR contactTreePtr->getPostProcMesh().tag_get_length(tag, tag_length);
1591 vec.resize(tag_length);
1592 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(tag, &face, 1,
1593 &*vec.begin());
1595 };
1596
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};
1599 Tag th;
1600 if (postProcMeshPtr) {
1601 CHKERR postProcMeshPtr->tag_get_handle(
1602 tag_name.c_str(), size, MB_TYPE_DOUBLE, th,
1603 MB_TAG_CREAT | MB_TAG_SPARSE, def_VAL);
1604 CHKERR postProcMeshPtr->tag_clear_data(th, &*mapGaussPtsPtr->begin(),
1605 mapGaussPtsPtr->size(), def_VAL);
1606 }
1607 return th;
1608 };
1609
1610 auto set_float_precision = [](const double x) {
1611 if (std::abs(x) < std::numeric_limits<float>::epsilon())
1612 return 0.;
1613 else
1614 return x;
1615 };
1616
1617 // scalars
1618 auto save_scal_tag = [&](auto &th, auto v, const int gg) {
1620 if (postProcMeshPtr) {
1621 v = set_float_precision(v);
1622 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1, &v);
1623 }
1625 };
1626
1627 // adjacencies
1628 auto get_fe_adjacencies = [this](auto fe_ent) {
1629 Range adj_faces;
1630 CHK_MOAB_THROW(getFaceFE()->mField.get_moab().get_adjacencies(
1631 &fe_ent, 1, 2, false, adj_faces, moab::Interface::UNION),
1632 "get adj");
1633 std::set<int> adj_ids;
1634 for (auto f : adj_faces) {
1635 adj_ids.insert(id_from_handle(f));
1636 }
1637 return adj_ids;
1638 };
1639
1640 auto get_face_id = [this](auto face) {
1641 int id;
1642 if (contactTreePtr->getPostProcMesh().tag_get_data(
1643 contactTreePtr->thEleId, &face, 1, &id) == MB_SUCCESS) {
1644 return id;
1645 }
1646 return -1;
1647 };
1648
1649 auto get_body_id = [this](auto face) {
1650 int id;
1651 if (contactTreePtr->getPostProcMesh().tag_get_data(
1652 contactTreePtr->thBodyId, &face, 1, &id) == MB_SUCCESS) {
1653 return id;
1654 }
1655 return -1;
1656 };
1657
1658 auto get_face_part = [this](auto face) {
1659 ParallelComm *pcomm_post_proc_mesh = ParallelComm::get_pcomm(
1660 &(contactTreePtr->getPostProcMesh()), MYPCOMM_INDEX);
1661 int part;
1662 if (contactTreePtr->getPostProcMesh().tag_get_data(
1663 pcomm_post_proc_mesh->part_tag(), &face, 1, &part) == MB_SUCCESS) {
1664 return part;
1665 }
1666 return -1;
1667 };
1668
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)
1673 return true;
1674 return false;
1675 };
1676
1677 // vectors
1678 VectorDouble3 v(3);
1679 FTensor::Tensor1<FTensor::PackPtr<double *, 0>, 3> t_v(&v[0], &v[1], &v[2]);
1680 auto save_vec_tag = [&](auto &th, auto &t_d, const int gg) {
1682 if (postProcMeshPtr) {
1683 t_v(i) = t_d(i);
1684 for (auto &a : v.data())
1685 a = set_float_precision(a);
1686 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1,
1687 &*v.data().begin());
1688 }
1690 };
1691
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);
1704
1705 auto t_ele_centre = get_ele_centre(getFTensor1Coords());
1706 auto ele_radius = get_ele_radius(t_ele_centre, getFTensor1Coords());
1707
1708 contactTreePtr->shadowDataMap.clear();
1709 auto &shadow_vec = contactTreePtr->shadowDataMap[fe_ent];
1710 shadow_vec.clear();
1711
1712 auto adj_fe_ids = get_fe_adjacencies(fe_ent);
1713
1714 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1715
1716 FTensor::Tensor1<double, 3> t_spatial_coords;
1717 t_spatial_coords(i) = t_coords(i) + t_disp_h1(i);
1718
1719 if (postProcMeshPtr) {
1720 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1721 }
1722
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())) {
1726
1727 auto body_id = get_body_id(face_close);
1728
1729 auto master_face_conn = get_face_conn(face_close);
1730 std::array<double, 9> master_coords;
1731 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1732 contactTreePtr->thSmallX, master_face_conn, 3,
1733 master_coords.data());
1734 std::array<double, 9> master_traction;
1735 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1736 contactTreePtr->thTraction, master_face_conn, 3,
1737 master_traction.data());
1738 auto t_normal_face_close = get_normal(master_coords);
1739 t_normal_face_close.normalize();
1740
1741 if (postProcMeshPtr) {
1742 double m = 1;
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);
1747 }
1748
1749 FTensor::Tensor1<double, 3> t_unit_ray;
1750 t_unit_ray(i) = -t_normal_face_close(i);
1751 FTensor::Tensor1<double, 3> t_ray_point;
1752 t_ray_point(i) =
1753 t_spatial_coords(i) -
1754 t_unit_ray(i) * ContactOps::airplane_ray_distance * ele_radius;
1755
1756 constexpr double eps = 1e-3;
1757 auto [faces_out, faces_dist] =
1758 get_faces_out(&t_ray_point(0), &t_unit_ray(0),
1759 2 * ContactOps::airplane_ray_distance * ele_radius,
1760 eps * ele_radius);
1761
1762 auto m = make_map(faces_out, faces_dist, t_ray_point, t_unit_ray,
1763 t_spatial_coords);
1764 for (auto m_it = m.begin(); m_it != m.end(); ++m_it) {
1765 auto face = m_it->second;
1766 if (face != face_close) {
1767
1768 if (
1769
1770 (adj_fe_ids.find(get_face_id(face)) == adj_fe_ids.end() ||
1771 get_face_part(face) != m_field.get_comm_rank())
1772
1773 ) {
1774
1775 shadow_vec.push_back(ContactTree::FaceData());
1776 shadow_vec.back().gaussPtNb = gg;
1777
1778 auto slave_face_conn = get_face_conn(face);
1779 std::array<double, 9> slave_coords;
1780 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1781 contactTreePtr->thSmallX, slave_face_conn, 3,
1782 slave_coords.data());
1783 auto t_normal_face = get_normal(slave_coords);
1784 std::array<double, 9> slave_tractions;
1785 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1786 contactTreePtr->thTraction, slave_face_conn, 3,
1787 slave_tractions.data());
1788
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());
1797
1798 t_slave_point(i) = t_ray_point(i) + m_it->first * t_unit_ray(i);
1799
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,
1805 loc_coords.data()),
1806 "get local coords");
1807 FTensor::Tensor1<double, 3> t_shape_fun;
1808 CHK_THROW_MESSAGE(Tools::shapeFunMBTRI<0>(&t_shape_fun(0),
1809 &loc_coords[0],
1810 &loc_coords[1], 1),
1811 "calc shape fun");
1812 FTensor::Index<'i', 3> i;
1813 FTensor::Index<'j', 3> j;
1814 FTensor::Tensor1<double, 3> t_point_out;
1815 t_point_out(i) = t_shape_fun(j) * t_elem_coords(j, i);
1816 return t_point_out;
1817 };
1818
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);
1823
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);
1828
1829 t_ray_point_data(i) = t_ray_point(i);
1830 t_unit_ray_data(i) = t_unit_ray(i);
1831
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());
1840
1841 shadow_vec.back().eleRadius = ele_radius;
1842
1843 // CHKERR get_tag_data(contactTreePtr->thIds, face,
1844 // shadow_vec.back().dofsSlaveIds);
1845 // CHKERR get_tag_data(contactTreePtr->thCoeff, face,
1846 // shadow_vec.back().dofsSlaveCoeff);
1847 // CHKERR get_tag_data(contactTreePtr->thBases, face,
1848 // shadow_vec.back().baseSlaveFuncs);
1849
1850 if (postProcMeshPtr) {
1851 auto [gap, tn_master, tn_slave, c, t_master_traction,
1852 t_slave_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);
1860 double m = 1;
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,
1866 gg);
1867 }
1868
1869 break;
1870 }
1871 }
1872 }
1873 break;
1874 }
1875 }
1876 next();
1877 }
1878
1880}
1881
1883 const std::string block_name, int dim) {
1884 Range r;
1885
1886 auto mesh_mng = m_field.getInterface<MeshsetsManager>();
1887 auto bcs = mesh_mng->getCubitMeshsetPtr(
1888
1889 std::regex((boost::format("%s(.*)") % block_name).str())
1890
1891 );
1892
1893 for (auto bc : bcs) {
1894 Range faces;
1895 CHK_MOAB_THROW(bc->getMeshsetIdEntitiesByDimension(m_field.get_moab(), dim,
1896 faces, true),
1897 "get meshset ents");
1898 r.merge(faces);
1899 }
1900
1901 return r;
1902};
1903
1904boost::shared_ptr<ForcesAndSourcesCore>
1906
1907 auto &m_field = ep.mField;
1908
1909 boost::shared_ptr<ContactTree> fe_contact_tree;
1910
1911 auto impl = [&]() {
1913
1914 /** Contact requires that body is marked */
1915 auto get_body_range = [&](auto name, int dim, auto sev) {
1916 std::map<int, Range> map;
1917
1918 for (auto m_ptr :
1919 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
1920 std::regex(
1921
1922 (boost::format("%s(.*)") % name).str()
1923
1924 ))
1925
1926 ) {
1927 Range ents;
1928 CHK_MOAB_THROW(m_ptr->getMeshsetIdEntitiesByDimension(
1929 m_field.get_moab(), dim, ents, true),
1930 "by dim");
1931 map[m_ptr->getMeshsetId()] = ents;
1932 MOFEM_LOG("EPSYNC", sev) << "Meshset: " << m_ptr->getMeshsetId() << " "
1933 << ents.size() << " entities";
1934 }
1935
1936 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), sev);
1937 return map;
1938 };
1939
1940 auto get_map_skin = [&](auto &&map) {
1941 ParallelComm *pcomm =
1942 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
1943
1944 Skinner skin(&m_field.get_moab());
1945 for (auto &m : map) {
1946 Range skin_faces;
1947 CHKERR skin.find_skin(0, m.second, false, skin_faces);
1948 CHK_MOAB_THROW(pcomm->filter_pstatus(
1949 skin_faces, PSTATUS_SHARED | PSTATUS_MULTISHARED,
1950 PSTATUS_NOT, -1, nullptr),
1951 "filter");
1952 m.second.swap(skin_faces);
1953 }
1954 return map;
1955 };
1956
1957 /* The above code is written in C++ and it appears to be defining and using
1958 various operations on boundary elements and side elements. */
1960 using BoundaryEleOp = BoundaryEle::UserDataOperator;
1961
1962 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1963
1964 auto calcs_side_traction = [&](auto &pip) {
1966 using EleOnSide =
1968 using SideEleOp = EleOnSide::UserDataOperator;
1969 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1970 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
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);
1982 };
1983
1984 auto add_contact_three = [&]() {
1986 auto tree_moab_ptr = boost::make_shared<moab::Core>();
1987 fe_contact_tree = boost::make_shared<ContactTree>(
1988 m_field, tree_moab_ptr, ep.spaceOrder,
1989 get_body_range("CONTACT", SPACE_DIM - 1, Sev::inform));
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));
2000 };
2001
2002 CHKERR add_contact_three();
2003
2005 };
2006
2007 CHK_THROW_MESSAGE(impl(), "createContactDetectionFiniteElement");
2008
2009 struct exclude_sdf {
2010 exclude_sdf(Range &&r) : map(r) {}
2011 bool operator()(FEMethod *fe_method_ptr) {
2012 auto ent = fe_method_ptr->getFEEntityHandle();
2013 if (map.find(ent) != map.end()) {
2014 return false;
2015 }
2016 return true;
2017 }
2018
2019 private:
2020 Range map;
2021 };
2022
2023 fe_contact_tree->exeTestHook =
2024 exclude_sdf(get_range_from_block(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2025
2026 return fe_contact_tree;
2027};
2028
2029static auto get_body_range(MoFEM::Interface &m_field, const std::string name,
2030 int dim) {
2031 std::map<int, Range> map;
2032
2033 for (auto m_ptr :
2034 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
2035
2036 (boost::format("%s(.*)") % name).str()
2037
2038 ))
2039
2040 ) {
2041 Range ents;
2042 CHK_MOAB_THROW(m_ptr->getMeshsetIdEntitiesByDimension(m_field.get_moab(),
2043 dim, ents, true),
2044 "by dim");
2045 map[m_ptr->getMeshsetId()] = ents;
2046 }
2047
2048 return map;
2049};
2050
2052 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2053 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2055
2056 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2057
2058 auto &m_field = ep.mField;
2059
2060 using BoundaryEle =
2062 using EleOnSide =
2064 using SideEleOp = EleOnSide::UserDataOperator;
2065 using BdyEleOp = BoundaryEle::UserDataOperator;
2066
2067 // First: Iterate over skeleton FEs adjacent to Domain FEs
2068 // Note: BoundaryEle, i.e. uses skeleton interation rule
2069 auto op_loop_skeleton_side = new OpLoopSide<BoundaryEle>(
2070 m_field, ep.contactElement, SPACE_DIM - 1, Sev::noisy);
2071
2072 auto rule_contact = [](int, int, int o) { return -1; };
2074
2075 auto set_rule_contact = [refine](
2076
2077 ForcesAndSourcesCore *fe_raw_ptr, int order_row,
2078 int order_col, int order_data
2079
2080 ) {
2082 auto rule = 2 * order_data;
2083 fe_raw_ptr->gaussPts = Tools::refineTriangleIntegrationPts(rule, refine);
2085 };
2086
2087 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2088 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2089 CHKERR
2090 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2091 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
2092 ep.frontAdjEdges);
2093
2094 // Second: Iterate over domain FEs adjacent to skelton, particularly
2095 // one domain element.
2096 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2097
2098 // Data storing contact fields
2099 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2100
2101 auto add_ops_domain_side = [&](auto &pip) {
2103 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
2104 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
2105 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
2106 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2107 boost::make_shared<CGGUserPolynomialBase>();
2108 CHKERR
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(
2114 broken_data_ptr));
2115 op_loop_domain_side->getOpPtrVector().push_back(
2117 ep.piolaStress, contact_common_data_ptr->contactTractionPtr()));
2118 pip.push_back(op_loop_domain_side);
2120 };
2121
2122 auto add_ops_contact_rhs = [&](auto &pip) {
2124 // get body id and SDF range
2125 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2126 get_body_range(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2127
2128 pip.push_back(new OpCalculateVectorFieldValues<3>(
2129 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2130 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2131 pip.push_back(
2133 pip.push_back(new OpTreeSearch(
2134 contact_tree_ptr, u_h1_ptr,
2135 contact_common_data_ptr->contactTractionPtr(),
2136 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1), nullptr,
2137 nullptr));
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));
2143
2145 };
2146
2147 // push ops to face/side pipeline
2148 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2149 CHKERR add_ops_contact_rhs(op_loop_skeleton_side->getOpPtrVector());
2150
2151 // Add skeleton to domain pipeline
2152 pip.push_back(op_loop_skeleton_side);
2153
2155};
2156
2158 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2159 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2161
2162 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2163 auto &m_field = ep.mField;
2164
2165 using BoundaryEle =
2167 using EleOnSide =
2169 using SideEleOp = EleOnSide::UserDataOperator;
2170 using BdyEleOp = BoundaryEle::UserDataOperator;
2171
2172 // First: Iterate over skeleton FEs adjacent to Domain FEs
2173 // Note: BoundaryEle, i.e. uses skeleton interation rule
2174 auto op_loop_skeleton_side = new OpLoopSide<BoundaryEle>(
2175 m_field, ep.contactElement, SPACE_DIM - 1, Sev::noisy);
2176
2177 auto rule_contact = [](int, int, int o) { return -1; };
2179
2180 auto set_rule_contact = [refine](
2181
2182 ForcesAndSourcesCore *fe_raw_ptr, int order_row,
2183 int order_col, int order_data
2184
2185 ) {
2187 auto rule = 2 * order_data;
2188 fe_raw_ptr->gaussPts = Tools::refineTriangleIntegrationPts(rule, refine);
2190 };
2191
2192 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2193 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2194 CHKERR
2195 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2196 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
2197 ep.frontAdjEdges);
2198
2199 // Second: Iterate over domain FEs adjacent to skelton, particularly
2200 // one domain element.
2201 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2202
2203 // Data storing contact fields
2204 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2205
2206 auto add_ops_domain_side = [&](auto &pip) {
2208 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
2209 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
2210 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
2211 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2212 boost::make_shared<CGGUserPolynomialBase>();
2213 CHKERR
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(
2219 broken_data_ptr));
2220 op_loop_domain_side->getOpPtrVector().push_back(
2222 ep.piolaStress, contact_common_data_ptr->contactTractionPtr()));
2223 pip.push_back(op_loop_domain_side);
2225 };
2226
2227 auto add_ops_contact_lhs = [&](auto &pip) {
2229 pip.push_back(new OpCalculateVectorFieldValues<3>(
2230 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2231 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2232 pip.push_back(
2234 pip.push_back(new OpTreeSearch(
2235 contact_tree_ptr, u_h1_ptr,
2236 contact_common_data_ptr->contactTractionPtr(),
2237 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1), nullptr,
2238 nullptr));
2239
2240 // get body id and SDF range
2241 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2242 get_body_range(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2243
2245 ep.contactDisp, ep.contactDisp, contact_common_data_ptr,
2246 contact_tree_ptr, contact_sfd_map_range_ptr));
2247 pip.push_back(
2249 ep.contactDisp, broken_data_ptr, contact_common_data_ptr,
2250 contact_tree_ptr, contact_sfd_map_range_ptr));
2251 pip.push_back(
2253 broken_data_ptr, ep.contactDisp, contact_common_data_ptr,
2254 contact_tree_ptr));
2255
2257 };
2258
2259 // push ops to face/side pipeline
2260 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2261 CHKERR add_ops_contact_lhs(op_loop_skeleton_side->getOpPtrVector());
2262
2263 // Add skeleton to domain pipeline
2264 pip.push_back(op_loop_skeleton_side);
2265
2267};
2268
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) {
2276
2277 auto &m_field = ep.mField;
2278 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2279 return new OpTreeSearch(
2280 contact_tree_ptr, u_h1_ptr, contact_traction_ptr,
2281 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1),
2282 post_proc_mesh_ptr, map_gauss_pts_ptr);
2283}
2284
2285} // 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.
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
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)
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
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)
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