v0.15.4
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
45 #pragma message "With ENABLE_PYTHON_BINDING"
46
47#else
48
49 #pragma message "Without ENABLE_PYTHON_BINDING"
50
51#endif
52
53#include <ContactOps.hpp>
54
55#include <EshelbianContact.hpp>
56
57extern "C" {
58void tricircumcenter3d_tp(double a[3], double b[3], double c[3],
59 double circumcenter[3], double *xi, double *eta);
60}
61
62namespace EshelbianPlasticity {
63
64struct ContactSDFPython : public ContactOps::SDFPython {
65 using ContactOps::SDFPython::SDFPython;
66};
67
68boost::shared_ptr<ContactSDFPython> setupContactSdf() {
69 boost::shared_ptr<ContactSDFPython> sdf_python_ptr;
70
71#ifdef ENABLE_PYTHON_BINDING
72
73 auto file_exists = [](std::string myfile) {
74 std::ifstream file(myfile.c_str());
75 if (file) {
76 return true;
77 }
78 return false;
79 };
80
81 char sdf_file_name[255] = "sdf.py";
82 CHKERR PetscOptionsGetString(PETSC_NULLPTR, PETSC_NULLPTR, "-sdf_file",
83 sdf_file_name, 255, PETSC_NULLPTR);
84
85 if (file_exists(sdf_file_name)) {
86 MOFEM_LOG("EP", Sev::inform) << sdf_file_name << " file found";
87 sdf_python_ptr = boost::make_shared<ContactSDFPython>();
88 CHKERR sdf_python_ptr->sdfInit(sdf_file_name);
89 ContactOps::sdfPythonWeakPtr = sdf_python_ptr;
90 MOFEM_LOG("EP", Sev::inform) << "SdfPython initialized";
91 } else {
92 MOFEM_LOG("EP", Sev::warning) << sdf_file_name << " file NOT found";
93 }
94
95#endif // ENABLE_PYTHON_BINDING
96
97 return sdf_python_ptr;
98};
99
101 : public PostProcBrokenMeshInMoabBase<FaceElementForcesAndSourcesCore> {
102
105
107 boost::shared_ptr<moab::Core> core_mesh_ptr, int max_order,
108 std::map<int, Range> &&body_map);
111 inline boost::shared_ptr<OrientedBoxTreeTool> &getTreeSurfPtr() {
112 return treeSurfPtr;
113 }
114 inline auto getRootSetSurf() { return rootSetSurf; }
115 int getMaxLevel() const { return refElementsMap.at(MBTRI)->defMaxLevel; }
116
117 friend struct OpMoveNode;
118 friend struct OpTreeSearch;
119
120 struct FaceData {
121 int gaussPtNb; //< integration points number
122 std::array<double, 3> slavePoint;
123 std::array<double, 3> masterPoint;
124 std::array<double, 3> rayPoint;
125 std::array<double, 3> unitRay;
126 double eleRadius;
127
128 // std::vector<int> dofsSlaveIds;
129 // std::vector<double> dofsSlaveCoeff;
130 // std::vector<double> baseSlaveFuncs;
131
132 std::array<double, 9> masterPointNodes;
133 std::array<double, 9> masterTractionNodes;
134 std::array<double, 9> slavePointNodes;
135 std::array<double, 9> slaveTractionNodes;
136
137 FaceData() = default;
138 };
139
140 using MapFaceData = std::map<EntityHandle, std::vector<FaceData>>;
141
142 inline auto findFaceDataVecPtr(EntityHandle fe_ent) {
143 auto &map_face_data = shadowDataMap;
144 auto it = map_face_data.find(fe_ent);
145 if (it == map_face_data.end()) {
146 return (std::vector<FaceData> *)nullptr;
147 }
148 return &(it->second);
149 }
150
151 inline auto getFaceDataPtr(std::vector<FaceData>::iterator &it, int gg,
152 std::vector<FaceData> *vec_ptr) {
153 FaceData *face_data_ptr = nullptr;
154 if (it != vec_ptr->end()) {
155 if (it->gaussPtNb == gg) {
156 face_data_ptr = &(*it);
157 ++it;
158 }
159 }
160 return face_data_ptr;
161 }
162
163protected:
165 boost::shared_ptr<OrientedBoxTreeTool> treeSurfPtr;
167
169 // Tag thCoeff;
170 // Tag thIds;
171 // Tag thBases;
176
178
179 std::map<int, Range> bodyMap;
180
181 const int maxOrder;
182};
183
184auto checkSdf(EntityHandle fe_ent, std::map<int, Range> &sdf_map_range) {
185 for (auto &m_sdf : sdf_map_range) {
186 if (m_sdf.second.find(fe_ent) != m_sdf.second.end())
187 return m_sdf.first;
188 }
189 return -1;
190}
191
192template <typename OP_PTR>
193auto getSdf(OP_PTR op_ptr, MatrixDouble &contact_disp, int block_id,
194 bool eval_hessian) {
195
196 auto nb_gauss_pts = op_ptr->getGaussPts().size2();
197
198 auto ts_time = op_ptr->getTStime();
199 auto ts_time_step = op_ptr->getTStimeStep();
202 ts_time_step = EshelbianCore::dynamicStep;
203 }
204
205 auto m_spatial_coords = ContactOps::get_spatial_coords(
206 op_ptr->getFTensor1CoordsAtGaussPts(),
207 getFTensor1FromMat<3>(contact_disp), nb_gauss_pts);
208 auto m_normals_at_pts = ContactOps::get_normalize_normals(
209 op_ptr->getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
211 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
212 block_id);
214 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
215 block_id);
216
217#ifndef NDEBUG
218 if (v_sdf.size() != nb_gauss_pts)
220 "Wrong number of integration pts");
221 if (m_grad_sdf.size2() != nb_gauss_pts)
223 "Wrong number of integration pts");
224 if (m_grad_sdf.size1() != 3)
225 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Should be size of 3");
226#endif // NDEBUG
227
228 if (eval_hessian) {
230 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
231 block_id);
232 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
233 m_hess_sdf);
234 } else {
235
236 return std::make_tuple(block_id, m_normals_at_pts, v_sdf, m_grad_sdf,
237 MatrixDouble(6, nb_gauss_pts, 0.));
238 }
239};
240
241/**
242 * @brief Calculate points data on contact surfaces
243 *
244 * @tparam T1
245 * @param unit_ray
246 * @param point
247 * @param elem_point_nodes
248 * @param elem_traction_nodes
249 * @param t_spatial_coords
250 * @return auto
251 */
252template <typename T1>
254
255 std::array<double, 3> &unit_ray, std::array<double, 3> &point,
256
257 std::array<double, 9> &elem_point_nodes,
258 std::array<double, 9> &elem_traction_nodes,
259
260 FTensor::Tensor1<T1, 3> &t_spatial_coords) {
261 FTensor::Index<'i', 3> i;
262
263 auto t_unit_ray = getFTensor1FromPtr<3>(unit_ray.data());
264 auto t_point = getFTensor1FromPtr<3>(point.data());
265
266 auto get_normal = [](auto &ele_coords) {
268 Tools::getTriNormal(ele_coords.data(), &t_normal(0));
269 return t_normal;
270 };
271
272 auto t_normal = get_normal(elem_point_nodes);
273 t_normal(i) /= t_normal.l2();
274
275 auto sn = t_normal(i) * t_point(i);
276 auto nm = t_normal(i) * t_spatial_coords(i);
277 auto nr = t_normal(i) * t_unit_ray(i);
278
279 auto gamma = (sn - nm) / nr;
280
281 FTensor::Tensor1<T1, 3> t_point_current;
282 t_point_current(i) = t_spatial_coords(i) + gamma * t_unit_ray(i);
283
284 auto get_local_point_shape_functions = [&](auto &&t_elem_coords,
285 auto &t_point) {
286 std::array<T1, 2> loc_coords;
289 &t_elem_coords(0, 0), &t_point(0), 1, loc_coords.data()),
290 "get local coords");
291 return FTensor::Tensor1<T1, 3>{N_MBTRI0(loc_coords[0], loc_coords[1]),
292 N_MBTRI1(loc_coords[0], loc_coords[1]),
293 N_MBTRI2(loc_coords[0], loc_coords[1])};
294 };
295
296 auto eval_position = [&](auto &&t_field, auto &t_shape_fun) {
297 FTensor::Index<'i', 3> i;
298 FTensor::Index<'j', 3> j;
299 FTensor::Tensor1<T1, 3> t_point_out;
300 t_point_out(i) = t_shape_fun(j) * t_field(j, i);
301 return t_point_out;
302 };
303
304 auto t_shape_fun = get_local_point_shape_functions(
305 getFTensor2FromPtr<3, 3>(elem_point_nodes.data()), t_point_current);
306 auto t_slave_point_updated = eval_position(
307 getFTensor2FromPtr<3, 3>(elem_point_nodes.data()), t_shape_fun);
308 auto t_traction_updated = eval_position(
309 getFTensor2FromPtr<3, 3>(elem_traction_nodes.data()), t_shape_fun);
310
311 return std::make_tuple(t_slave_point_updated, t_traction_updated, t_normal);
312};
313
314template <typename T1>
316
317 ContactTree::FaceData *face_data_ptr,
318 FTensor::Tensor1<T1, 3> &t_spatial_coords
319
320) {
321
322 return multiPoint(face_data_ptr->unitRay, face_data_ptr->masterPoint,
323 face_data_ptr->masterPointNodes,
324 face_data_ptr->masterTractionNodes, t_spatial_coords);
325};
326
327template <typename T1>
329
330 ContactTree::FaceData *face_data_ptr,
331 FTensor::Tensor1<T1, 3> &t_spatial_coords
332
333) {
334
335 return multiPoint(face_data_ptr->unitRay, face_data_ptr->slavePoint,
336 face_data_ptr->slavePointNodes,
337 face_data_ptr->slaveTractionNodes, t_spatial_coords);
338};
339
340/**
341 * Evaluate gap and tractions between master and slave points
342 */
343template <typename T1>
345 FTensor::Tensor1<T1, 3> &t_spatial_coords) {
346
347 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
348 multiSlavePoint(face_data_ptr, t_spatial_coords);
349 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
350 multiMasterPoint(face_data_ptr, t_spatial_coords);
351
352 FTensor::Index<'i', 3> i;
353
355 t_normal(i) = t_master_normal(i) - t_slave_normal(i);
356 t_normal.normalize();
357
358 auto gap = t_normal(i) * (t_slave_point_current(i) - t_spatial_coords(i));
359 auto tn_master = t_master_traction_current(i) * t_normal(i);
360 auto tn_slave = t_slave_traction_current(i) * t_normal(i);
361 auto tn = std::max(-tn_master, tn_slave);
362
363 return std::make_tuple(gap, tn_master, tn_slave,
364 ContactOps::constrain(gap, tn),
365 t_master_traction_current, t_slave_traction_current);
366};
367
369
370/** Caluclate rhs term for contact
371 */
372template <typename T1, typename T2, typename T3>
374
375 ContactTree::FaceData *face_data_ptr, FTensor::Tensor1<T1, 3> &t_coords,
376 FTensor::Tensor1<T2, 3> &t_spatial_coords,
377 FTensor::Tensor1<T3, 3> &t_master_traction, MultiPointRhsType type,
378 bool debug = false
379
380) {
381 FTensor::Index<'i', 3> i;
382 FTensor::Index<'j', 3> j;
383
385 t_u(i) = t_spatial_coords(i) - t_coords(i);
386
388
389 switch (type) {
391
392 auto [t_slave_point_current, t_slave_traction_current, t_slave_normal] =
393 multiSlavePoint(face_data_ptr, t_spatial_coords);
394 auto [t_master_point_current, t_master_traction_current, t_master_normal] =
395 multiMasterPoint(face_data_ptr, t_spatial_coords);
396
397 // average normal surface
399 t_normal(i) = t_master_normal(i) - t_slave_normal(i);
400 t_normal.normalize();
401
402 // get projection operators
404 t_P(i, j) = t_normal(i) * t_normal(j);
406 t_Q(i, j) = kronecker_delta(i, j) - t_P(i, j);
407
408 constexpr double beta = 0.5;
409
410 auto zeta = 1e-8 * face_data_ptr->eleRadius;
413
414 // this is regularised std::min(d,s)
415 auto f_min_gap = [zeta](auto d, auto s) {
416 return 0.5 * (d + s - std::sqrt((d - s) * (d - s) + zeta));
417 };
418 // this derivative of regularised std::min(d,s)
419 auto f_diff_min_gap = [zeta](auto d, auto s) {
420 return 0.5 * (1 - (d - s) / std::sqrt((d - s) * (d - s) + zeta));
421 };
422
423 // add penalty on penetration side
424 auto f_barrier = [alpha1, alpha2, f_min_gap, f_diff_min_gap](auto g,
425 auto tn) {
426 auto d = alpha1 * g;
427 auto b1 =
428 0.5 * (tn + f_min_gap(d, tn) + f_diff_min_gap(d, tn) * (d - tn));
429 auto b2 = alpha2 * f_min_gap(g, 0) * g;
430 return b1 - b2;
431 };
432
433 FTensor::Tensor1<T2, 3> t_gap_vec;
434 t_gap_vec(i) = beta * t_spatial_coords(i) +
435 (1 - beta) * t_slave_point_current(i) - t_spatial_coords(i);
437 t_traction_vec(i) =
438 -beta * t_master_traction(i) + (beta - 1) * t_slave_traction_current(i);
439
440 auto t_gap = t_normal(i) * t_gap_vec(i);
441 auto t_tn = t_normal(i) * t_traction_vec(i);
442 auto barrier = f_barrier(t_gap, t_tn);
443
445 // add penalty on penetration side
446 t_traction_bar(i) = t_normal(i) * barrier;
447
448 t_rhs(i) =
449
450 t_Q(i, j) * t_master_traction(j) +
451
452 t_P(i, j) * (t_master_traction(j) - t_traction_bar(j));
453
454 if (debug) {
455 auto is_nan_or_inf = [](double value) -> bool {
456 return std::isnan(value) || std::isinf(value);
457 };
458
459 double v = std::complex<double>(t_rhs(i) * t_rhs(i)).real();
460 if (is_nan_or_inf(v)) {
461 MOFEM_LOG_CHANNEL("SELF");
462 MOFEM_LOG("SELF", Sev::error) << "t_rhs " << t_rhs;
463 CHK_MOAB_THROW(MOFEM_DATA_INCONSISTENCY, "rhs is nan or inf");
464 }
465 }
466
467 } break;
470 break;
471 }
472
473 return t_rhs;
474};
475
478 const std::string row_field_name,
479 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
480 boost::shared_ptr<ContactTree> contact_tree_ptr,
481 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr = nullptr);
482
484
485private:
486 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
487 boost::shared_ptr<ContactTree> contactTreePtr;
488 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
489};
490
492 const std::string row_field_name,
493 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
494 boost::shared_ptr<ContactTree> contact_tree_ptr,
495 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
496 : ContactOps::AssemblyBoundaryEleOp(row_field_name, row_field_name,
497 ContactOps::BoundaryEleOp::OPROW),
498 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
499 sdfMapRangePtr(sdf_map_range_ptr) {
500 CHK_THROW_MESSAGE(PetscOptionsGetScalar(PETSC_NULLPTR, "", "-cn",
502 PETSC_NULLPTR),
503 "get cn failed");
506 PetscOptionsGetScalar(PETSC_NULLPTR, "", "-alpha_contact_const",
507 &ContactOps::alpha_contact_const, PETSC_NULLPTR),
508 "get alpha contact failed");
510 PETSC_NULLPTR, "", "-alpha_contact_quadratic",
512 "get alpha contact failed");
514 PetscOptionsGetScalar(PETSC_NULLPTR, "", "-airplane_ray_distance",
515 &ContactOps::airplane_ray_distance, PETSC_NULLPTR),
516 "get alpha contact failed");
517
518 MOFEM_LOG("EP", Sev::inform) << "cn " << ContactOps::cn_contact;
519 MOFEM_LOG("EP", Sev::inform)
520 << "alpha_contact_const " << ContactOps::alpha_contact_const;
521 MOFEM_LOG("EP", Sev::inform)
522 << "alpha_contact_quadratic " << ContactOps::alpha_contact_quadratic;
523 MOFEM_LOG("EP", Sev::inform)
524 << "airplane_ray_distance " << ContactOps::airplane_ray_distance;
525}
526
530
531 FTensor::Index<'i', 3> i;
532 FTensor::Index<'j', 3> j;
533 FTensor::Index<'k', 3> k;
534 FTensor::Index<'l', 3> l;
535
536 const size_t nb_gauss_pts = getGaussPts().size2();
537
538#ifndef NDEBUG
539 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
540 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
541 "Wrong number of integration pts %ld != %ld",
542 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
543 }
544#endif // !NDEBUG
545
546 auto &nf = locF;
547 locF.clear();
548
549 auto t_w = getFTensor0IntegrationWeight();
550 auto t_coords = getFTensor1CoordsAtGaussPts();
551 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
552 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
553
554 // placeholder to pass boundary block id to python. default SDF is set on
555 // block = -1, one can choose different block by making block "CONTACT_SDF",
556 // then specific SDF can be set to that block.
557 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
558 getSdf(this, commonDataPtr->contactDisp,
559 checkSdf(getFEEntityHandle(), *sdfMapRangePtr), false);
560
561 auto t_sdf_v = getFTensor0FromVec(v_sdf);
562 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
563 auto t_normalize_normal = getFTensor1FromMat<3>(m_normals_at_pts);
564
565 auto next = [&]() {
566 ++t_w;
567 ++t_coords;
568 ++t_disp;
569 ++t_traction;
570 ++t_normalize_normal;
571 ++t_sdf_v;
572 ++t_grad_sdf_v;
573 };
574
575 auto face_data_vec_ptr =
576 contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
577 auto face_gauss_pts_it = face_data_vec_ptr->begin();
578
579 auto nb_base_functions = data.getN().size2();
580 auto t_base = data.getFTensor0N();
581 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
582
584 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
585 face_data_vec_ptr);
586
587 auto check_face_contact = [&]() {
588 if (checkSdf(getFEEntityHandle(), *sdfMapRangePtr) != -1)
589 return false;
590
591 if (face_data_ptr) {
592 return true;
593 }
594 return false;
595 };
596
597#ifdef ENABLE_PYTHON_BINDING
598 double c = 0.;
599 if (ContactOps::sdfPythonWeakPtr.lock()) {
600 auto tn = t_traction(i) * t_grad_sdf_v(i);
601 c = ContactOps::constrain(t_sdf_v, tn);
602 }
603#else
604 constexpr double c = 0;
605#endif
606
607 if (!c && check_face_contact()) {
608 FTensor::Tensor1<double, 3> t_spatial_coords;
609 t_spatial_coords(i) = t_coords(i) + t_disp(i);
610 auto t_rhs_tmp = multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
611 t_traction, MultiPointRhsType::U, true);
612 t_rhs(i) = t_rhs_tmp(i);
613
614 } else {
615
616#ifdef ENABLE_PYTHON_BINDING
617 auto inv_cn = 1. / ContactOps::cn_contact;
618
619 if (ContactOps::sdfPythonWeakPtr.lock()) {
621 t_cP(i, j) = (c * t_grad_sdf_v(i)) * t_grad_sdf_v(j);
622 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
623 t_rhs(i) = t_cQ(i, j) * t_traction(j) +
624 (c * inv_cn * t_sdf_v) * t_grad_sdf_v(i);
625 } else {
626 t_rhs(i) = t_traction(i);
627 }
628#else
629 t_rhs(i) = t_traction(i);
630#endif
631 }
632
633 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
634 const double alpha = t_w * getMeasure();
635
636 size_t bb = 0;
637 for (; bb != nbRows / 3; ++bb) {
638 const double beta = alpha * t_base;
639 t_nf(i) -= beta * t_rhs(i);
640 ++t_nf;
641 ++t_base;
642 }
643 for (; bb < nb_base_functions; ++bb)
644 ++t_base;
645
646 next();
647 }
648
650}
651
652template <AssemblyType A, IntegrationType I> struct OpConstrainBoundaryHDivRhs;
653
654template <AssemblyType A>
656 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
657
660
662 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
663 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
664 boost::shared_ptr<ContactTree> contact_tree_ptr);
665
666 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data);
667
668private:
669 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
670 boost::shared_ptr<ContactTree> contactTreePtr;
671};
672
673template <AssemblyType A>
676 boost::shared_ptr<std::vector<BrokenBaseSideData>>
677 broken_base_side_data,
678 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
679 boost::shared_ptr<ContactTree> contact_tree_ptr)
680 : OP(broken_base_side_data), commonDataPtr(common_data_ptr),
681 contactTreePtr(contact_tree_ptr) {}
682
683template <AssemblyType A>
687
692
693 const size_t nb_gauss_pts = OP::getGaussPts().size2();
694
695#ifndef NDEBUG
696 if (commonDataPtr->contactDisp.size2() != nb_gauss_pts) {
697 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
698 "Wrong number of integration pts %ld != %ld",
699 commonDataPtr->contactDisp.size2(), nb_gauss_pts);
700 }
701#endif // !NDEBUG
702
703 auto &nf = OP::locF;
704 OP::locF.clear();
705
706 auto t_w = OP::getFTensor0IntegrationWeight();
707 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
708 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
709
710 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
711 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
712
713 auto next = [&]() {
714 ++t_w;
715 ++t_disp;
716 ++t_traction;
717 ++t_coords;
718 ++t_material_normal;
719 };
720
721 auto face_data_vec_ptr =
722 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
723 auto face_gauss_pts_it = face_data_vec_ptr->begin();
724
725 auto nb_base_functions = data.getN().size2() / 3;
726 auto t_base = data.getFTensor1N<3>();
727 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
728
729 auto t_nf = getFTensor1FromPtr<3>(&nf[0]);
730 const double alpha = t_w / 2.;
731
732 size_t bb = 0;
733 for (; bb != OP::nbRows / 3; ++bb) {
734 const double beta = alpha * t_base(i) * t_material_normal(i);
735 t_nf(i) += beta * t_disp(i);
736 ++t_nf;
737 ++t_base;
738 }
739 for (; bb < nb_base_functions; ++bb)
740 ++t_base;
741
742 next();
743 }
744
746}
747
749
751 const std::string row_field_name, const std::string col_field_name,
752 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
753 boost::shared_ptr<ContactTree> contact_tree_ptr,
754 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr = nullptr);
755
758
759private:
760 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
761 boost::shared_ptr<ContactTree> contactTreePtr;
762 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
763};
764
766 const std::string row_field_name, const std::string col_field_name,
767 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
768 boost::shared_ptr<ContactTree> contact_tree_ptr,
769 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
770 : ContactOps::AssemblyBoundaryEleOp(row_field_name, col_field_name,
771 ContactOps::BoundaryEleOp::OPROWCOL),
772 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
773 sdfMapRangePtr(sdf_map_range_ptr) {
774
775 sYmm = false;
776}
777
780 EntitiesFieldData::EntData &col_data) {
782
783 using namespace ContactOps;
784
785 FTensor::Index<'i', 3> i;
786 FTensor::Index<'j', 3> j;
787 FTensor::Index<'k', 3> k;
788
789 auto nb_rows = row_data.getIndices().size();
790 auto nb_cols = col_data.getIndices().size();
791
792 auto &locMat = AssemblyBoundaryEleOp::locMat;
793 locMat.resize(nb_rows, nb_cols, false);
794 locMat.clear();
795
796 if (nb_cols && nb_rows) {
797
798 auto nb_gauss_pts = getGaussPts().size2();
799 auto t_w = getFTensor0IntegrationWeight();
800 auto t_coords = getFTensor1CoordsAtGaussPts();
801 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
802 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
803
804 // placeholder to pass boundary block id to python
805 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
806 getSdf(this, commonDataPtr->contactDisp,
807 checkSdf(getFEEntityHandle(), *sdfMapRangePtr), true);
808
809 auto t_sdf_v = getFTensor0FromVec(v_sdf);
810 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
811 auto t_hess_sdf_v = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
812 auto t_normalized_normal = getFTensor1FromMat<3>(m_normals_at_pts);
813
814 auto next = [&]() {
815 ++t_w;
816 ++t_coords;
817 ++t_disp;
818 ++t_traction;
819 ++t_sdf_v;
820 ++t_grad_sdf_v;
821 ++t_hess_sdf_v;
822 ++t_normalized_normal;
823 };
824
825 auto face_data_vec_ptr =
826 contactTreePtr->findFaceDataVecPtr(getFEEntityHandle());
827 auto face_gauss_pts_it = face_data_vec_ptr->begin();
828
829 auto t_row_base = row_data.getFTensor0N();
830 auto nb_face_functions = row_data.getN().size2() / 3;
831 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
832
833 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
834
835 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
836 face_data_vec_ptr);
837
838 auto check_face_contact = [&]() {
839 if (checkSdf(getFEEntityHandle(), *sdfMapRangePtr) != -1)
840 return false;
841
842 if (face_data_ptr) {
843 return true;
844 }
845 return false;
846 };
847
849
850#ifdef ENABLE_PYTHON_BINDING
851 double c = 0.;
852 if (ContactOps::sdfPythonWeakPtr.lock()) {
853 auto tn = t_traction(i) * t_grad_sdf_v(i);
854 c = ContactOps::constrain(t_sdf_v, tn);
855 }
856#else
857 constexpr double c = 0;
858#endif
859
860 if (!c && check_face_contact()) {
861 FTensor::Tensor1<double, 3> t_spatial_coords;
862 t_spatial_coords(i) = t_coords(i) + t_disp(i);
863 constexpr double eps = std::numeric_limits<float>::epsilon();
864 for (auto ii = 0; ii < 3; ++ii) {
865 FTensor::Tensor1<std::complex<double>, 3> t_spatial_coords_cx{
866 t_spatial_coords(0), t_spatial_coords(1), t_spatial_coords(2)};
867 t_spatial_coords_cx(ii) += eps * 1i;
868 auto t_rhs_tmp =
869 multiPointRhs(face_data_ptr, t_coords, t_spatial_coords_cx,
870 t_traction, MultiPointRhsType::U);
871 for (int jj = 0; jj != 3; ++jj) {
872 auto v = t_rhs_tmp(jj).imag();
873 t_res_dU(jj, ii) = v / eps;
874 }
875 }
876
877 } else {
878
879#ifdef ENABLE_PYTHON_BINDING
880
881 if (ContactOps::sdfPythonWeakPtr.lock()) {
882 auto inv_cn = 1. / ContactOps::cn_contact;
883 t_res_dU(i, j) =
884
885 (-c) * (t_hess_sdf_v(i, j) * t_grad_sdf_v(k) * t_traction(k) +
886 t_grad_sdf_v(i) * t_hess_sdf_v(k, j) * t_traction(k))
887
888 + (c * inv_cn) * (t_sdf_v * t_hess_sdf_v(i, j) +
889
890 t_grad_sdf_v(j) * t_grad_sdf_v(i));
891 } else {
892 t_res_dU(i, j) = 0;
893 }
894#else
895 t_res_dU(i, j) = 0;
896#endif
897 }
898
899 auto alpha = t_w * getMeasure();
900
901 size_t rr = 0;
902 for (; rr != nb_rows / 3; ++rr) {
903
904 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
905 auto t_col_base = col_data.getFTensor0N(gg, 0);
906
907 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
908 auto beta = alpha * t_row_base * t_col_base;
909 t_mat(i, j) -= beta * t_res_dU(i, j);
910 ++t_col_base;
911 ++t_mat;
912 }
913
914 ++t_row_base;
915 }
916 for (; rr < nb_face_functions; ++rr)
917 ++t_row_base;
918
919 next();
920 }
921 }
922
924}
925
926template <AssemblyType A, IntegrationType I> struct OpConstrainBoundaryL2Lhs_dP;
927
928template <AssemblyType A>
930 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
931
934
936 std::string row_field_name,
937 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
938 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
939 boost::shared_ptr<ContactTree> contact_tree_ptr,
940 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr = nullptr);
941
942 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
944
945private:
946 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
947 boost::shared_ptr<ContactTree> contactTreePtr;
948 boost::shared_ptr<std::map<int, Range>> sdfMapRangePtr;
949};
950
951template <AssemblyType A>
954 std::string row_field_name,
955 boost::shared_ptr<std::vector<BrokenBaseSideData>>
956 broken_base_side_data,
957 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
958 boost::shared_ptr<ContactTree> contact_tree_ptr,
959 boost::shared_ptr<std::map<int, Range>> sdf_map_range_ptr)
960 : OP(row_field_name, broken_base_side_data, false, false),
961 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr),
962 sdfMapRangePtr(sdf_map_range_ptr) {
963 OP::sYmm = false;
964}
965
966template <AssemblyType A>
970 EntitiesFieldData::EntData &col_data) {
972
973 using namespace ContactOps;
974
975 FTensor::Index<'i', 3> i;
976 FTensor::Index<'j', 3> j;
977 FTensor::Index<'k', 3> k;
978
979 auto nb_rows = row_data.getIndices().size();
980 auto nb_cols = col_data.getIndices().size();
981
982 auto &locMat = AssemblyBoundaryEleOp::locMat;
983 locMat.resize(nb_rows, nb_cols, false);
984 locMat.clear();
985
986 if (nb_cols && nb_rows) {
987
988 const size_t nb_gauss_pts = OP::getGaussPts().size2();
989
990 auto t_w = OP::getFTensor0IntegrationWeight();
991 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
992 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
993 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
994 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
995
996 // placeholder to pass boundary block id to python
997 auto [block_id, m_normals_at_pts, v_sdf, m_grad_sdf, m_hess_sdf] =
998 getSdf(this, commonDataPtr->contactDisp,
999 checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr), false);
1000
1001 auto t_sdf_v = getFTensor0FromVec(v_sdf);
1002 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
1003
1004 auto next = [&]() {
1005 ++t_w;
1006 ++t_disp;
1007 ++t_traction;
1008 ++t_coords;
1009 ++t_material_normal;
1010 ++t_sdf_v;
1011 ++t_grad_sdf_v;
1012 };
1013
1014 auto face_data_vec_ptr =
1015 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1016 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1017
1018 auto t_row_base = row_data.getFTensor0N();
1019 auto nb_face_functions = row_data.getN().size2();
1020
1021 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
1022
1023 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1024
1025 auto face_data_ptr = contactTreePtr->getFaceDataPtr(face_gauss_pts_it, gg,
1026 face_data_vec_ptr);
1027
1028 auto check_face_contact = [&]() {
1029 if (checkSdf(OP::getFEEntityHandle(), *sdfMapRangePtr) != -1)
1030 return false;
1031
1032 if (face_data_ptr) {
1033 return true;
1034 }
1035 return false;
1036 };
1037
1039
1040#ifdef ENABLE_PYTHON_BINDING
1041 double c = 0.;
1042 if (ContactOps::sdfPythonWeakPtr.lock()) {
1043 auto tn = t_traction(i) * t_grad_sdf_v(i);
1044 c = ContactOps::constrain(t_sdf_v, tn);
1045 }
1046#else
1047 constexpr double c = 0;
1048#endif
1049
1050 if (!c && check_face_contact()) {
1051 FTensor::Tensor1<double, 3> t_spatial_coords;
1052 t_spatial_coords(i) = t_coords(i) + t_disp(i);
1053 constexpr double eps = std::numeric_limits<float>::epsilon();
1054 for (auto ii = 0; ii != 3; ++ii) {
1055 FTensor::Tensor1<std::complex<double>, 3> t_traction_cx{
1056 t_traction(0), t_traction(1), t_traction(2)};
1057 t_traction_cx(ii) += eps * 1i;
1058 auto t_rhs_tmp =
1059 multiPointRhs(face_data_ptr, t_coords, t_spatial_coords,
1060 t_traction_cx, MultiPointRhsType::U);
1061 for (int jj = 0; jj != 3; ++jj) {
1062 auto v = t_rhs_tmp(jj).imag();
1063 t_res_dP(jj, ii) = v / eps;
1064 }
1065 }
1066 } else {
1067
1068#ifdef ENABLE_PYTHON_BINDING
1069 if (ContactOps::sdfPythonWeakPtr.lock()) {
1071 t_cP(i, j) = (c * t_grad_sdf_v(i)) * t_grad_sdf_v(j);
1072 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
1073 t_res_dP(i, j) = t_cQ(i, j);
1074 } else {
1075 t_res_dP(i, j) = t_kd(i, j);
1076 }
1077#else
1078 t_res_dP(i, j) = t_kd(i, j);
1079#endif
1080 }
1081
1082 const double alpha = t_w / 2.;
1083 size_t rr = 0;
1084 for (; rr != nb_rows / 3; ++rr) {
1085
1086 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1087 auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
1088
1089 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
1090 auto col_base = t_col_base(i) * t_material_normal(i);
1091 const double beta = alpha * t_row_base * col_base;
1092 t_mat(i, j) -= beta * t_res_dP(i, j);
1093 ++t_col_base;
1094 ++t_mat;
1095 }
1096
1097 ++t_row_base;
1098 }
1099 for (; rr < nb_face_functions; ++rr)
1100 ++t_row_base;
1101
1102 next();
1103 }
1104 }
1105
1107}
1108
1109template <AssemblyType A, IntegrationType I>
1111
1112template <AssemblyType A>
1114 : public FormsIntegrators<FaceUserDataOperator>::Assembly<A>::OpBrokenBase {
1115
1118
1120 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
1121 std::string col_field_name,
1122 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1123 boost::shared_ptr<ContactTree> contact_tree_ptr);
1124
1125 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
1126 EntitiesFieldData::EntData &col_data);
1127
1128private:
1129 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
1130 boost::shared_ptr<ContactTree> contactTreePtr;
1131};
1132
1133template <AssemblyType A>
1136 boost::shared_ptr<std::vector<BrokenBaseSideData>>
1137 broken_base_side_data,
1138 std::string col_field_name,
1139 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1140 boost::shared_ptr<ContactTree> contact_tree_ptr)
1141 : OP(col_field_name, broken_base_side_data, true, true),
1142 commonDataPtr(common_data_ptr), contactTreePtr(contact_tree_ptr) {
1143 OP::sYmm = false;
1144}
1145
1146template <AssemblyType A>
1150 EntitiesFieldData::EntData &row_data) {
1152
1153 // Note: col_data and row_data are swapped in this function, we going to
1154 // transpose locMat at the end
1155
1156 using namespace ContactOps;
1157
1161
1162 auto nb_rows = row_data.getIndices().size();
1163 auto nb_cols = col_data.getIndices().size();
1164
1165 auto &locMat = AssemblyBoundaryEleOp::locMat;
1166 locMat.resize(nb_rows, nb_cols, false);
1167 locMat.clear();
1168
1169 if (nb_cols && nb_rows) {
1170
1171 const size_t nb_gauss_pts = OP::getGaussPts().size2();
1172
1173 auto t_w = OP::getFTensor0IntegrationWeight();
1174 auto t_disp = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
1175 auto t_traction = getFTensor1FromMat<3>(commonDataPtr->contactTraction);
1176 auto t_coords = OP::getFTensor1CoordsAtGaussPts();
1177 auto t_material_normal = OP::getFTensor1NormalsAtGaussPts();
1178
1179 auto next = [&]() {
1180 ++t_w;
1181 ++t_disp;
1182 ++t_traction;
1183 ++t_coords;
1184 ++t_material_normal;
1185 };
1186
1187 constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
1188
1189 auto face_data_vec_ptr =
1190 contactTreePtr->findFaceDataVecPtr(OP::getFEEntityHandle());
1191 auto face_gauss_pts_it = face_data_vec_ptr->begin();
1192
1193 auto t_row_base = row_data.getFTensor1N<3>();
1194 auto nb_face_functions = row_data.getN().size2() / 3;
1195 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1196
1197 const auto alpha = t_w / 2.;
1198
1199 size_t rr = 0;
1200 for (; rr != nb_rows / 3; ++rr) {
1201
1202 auto row_base = alpha * (t_row_base(i) * t_material_normal(i));
1203
1204 auto t_mat = getFTensor2FromArray<3, 3, 3>(locMat, 3 * rr);
1205 auto t_col_base = col_data.getFTensor0N(gg, 0);
1206
1207 for (size_t cc = 0; cc != nb_cols / 3; ++cc) {
1208 const auto beta = row_base * t_col_base;
1209 t_mat(i, j) += beta * t_kd(i, j);
1210 ++t_col_base;
1211 ++t_mat;
1212 }
1213
1214 ++t_row_base;
1215 }
1216 for (; rr < nb_face_functions; ++rr)
1217 ++t_row_base;
1218
1219 next();
1220 }
1221 }
1222
1223 locMat = trans(locMat);
1224
1226}
1227
1229 boost::shared_ptr<moab::Core> core_mesh_ptr,
1230 int max_order, std::map<int, Range> &&body_map)
1231 : Base(m_field, core_mesh_ptr, "contact"), maxOrder(max_order),
1232 bodyMap(body_map) {
1233
1234 auto ref_ele_ptr = boost::make_shared<PostProcGenerateRefMesh<MBTRI>>();
1235 ref_ele_ptr->hoNodes =
1236 PETSC_FALSE; ///< So far only linear geometry is implemented for contact
1237
1238 CHK_THROW_MESSAGE(ref_ele_ptr->getOptions(optionsPrefix), "getOptions");
1239 CHK_THROW_MESSAGE(ref_ele_ptr->generateReferenceElementMesh(),
1240 "Error when generating reference element");
1241
1242 MOFEM_LOG("EP", Sev::inform) << "Contact hoNodes " << ref_ele_ptr->hoNodes
1243 ? "true"
1244 : "false";
1245 MOFEM_LOG("EP", Sev::inform)
1246 << "Contact maxOrder " << ref_ele_ptr->defMaxLevel;
1247
1248 refElementsMap[MBTRI] = ref_ele_ptr;
1249
1250 int def_ele_id = -1;
1251 CHKERR getPostProcMesh().tag_get_handle("ELE_ID", 1, MB_TYPE_INTEGER, thEleId,
1252 MB_TAG_CREAT | MB_TAG_DENSE,
1253 &def_ele_id);
1254 CHKERR getPostProcMesh().tag_get_handle("BODY_ID", 1, MB_TYPE_INTEGER,
1255 thBodyId, MB_TAG_CREAT | MB_TAG_DENSE,
1256 &def_ele_id);
1257
1258 // std::vector<int> def_ids(3 * nb_bases, 0);
1259 // CHKERR getPostProcMesh().tag_get_handle("IDS", 3 * nb_bases,
1260 // MB_TYPE_INTEGER,
1261 // thIds, MB_TAG_CREAT | MB_TAG_DENSE,
1262 // &*def_ids.begin());
1263 // std::vector<double> def_coeffs(3 * nb_bases, 0);
1264 // CHKERR getPostProcMesh().tag_get_handle("COEFF", 3 * nb_bases,
1265 // MB_TYPE_DOUBLE,
1266 // thCoeff, MB_TAG_CREAT |
1267 // MB_TAG_DENSE,
1268 // &*def_coeffs.begin());
1269 // std::vector<double> def_basses(nb_bases, 0);
1270 // CHKERR getPostProcMesh().tag_get_handle("BASES", nb_bases, MB_TYPE_DOUBLE,
1271 // thBases, MB_TAG_CREAT |
1272 // MB_TAG_DENSE,
1273 // &*def_basses.begin());
1274
1275 std::array<double, 3> def_small_x{0., 0., 0.};
1276 CHKERR getPostProcMesh().tag_get_handle("x", 3, MB_TYPE_DOUBLE, thSmallX,
1277 MB_TAG_CREAT | MB_TAG_DENSE,
1278 def_small_x.data());
1279 CHKERR getPostProcMesh().tag_get_handle("X", 3, MB_TYPE_DOUBLE, thLargeX,
1280 MB_TAG_CREAT | MB_TAG_DENSE,
1281 def_small_x.data());
1282
1283 std::array<double, 3> def_tractions{0., 0., 0.};
1284 CHKERR getPostProcMesh().tag_get_handle(
1285 "TRACTION", 3, MB_TYPE_DOUBLE, thTraction, MB_TAG_CREAT | MB_TAG_DENSE,
1286 &*def_tractions.begin());
1287}
1288
1294
1297
1299 shadowDataMap.clear();
1300
1301 PetscBarrier(nullptr);
1302
1303 ParallelComm *pcomm_post_proc_mesh =
1304 ParallelComm::get_pcomm(&getPostProcMesh(), MYPCOMM_INDEX);
1305 if (pcomm_post_proc_mesh == nullptr)
1306 SETERRQ(PETSC_COMM_WORLD, MOFEM_DATA_INCONSISTENCY, "PComm not allocated");
1307
1308 // MOFEM_LOG("EP", Sev::verbose) << "ContactTree broadcast entities";
1309
1310 auto brodacts = [&](auto &brodacts_ents) {
1312 Range recived_ents;
1313 for (auto rr = 0; rr != mField.get_comm_size(); ++rr) {
1314 pcomm_post_proc_mesh->broadcast_entities(
1315 rr, rr == mField.get_comm_rank() ? brodacts_ents : recived_ents,
1316 false, true);
1317 }
1319 };
1320
1321 Range brodacts_ents;
1322 CHKERR getPostProcMesh().get_entities_by_dimension(0, 2, brodacts_ents, true);
1323 CHKERR brodacts(brodacts_ents);
1324
1325 // MOFEM_LOG("EP", Sev::verbose) << "ContactTree build tree";
1326
1327 Range ents;
1328 CHKERR getPostProcMesh().get_entities_by_dimension(0, 2, ents, true);
1329 CHKERR buildTree(ents);
1330
1331#ifndef NDEBUG
1332 CHKERR writeFile("debug_tree.h5m");
1333#endif // NDEBUG
1334
1336}
1337
1340 // MOFEM_LOG("SELF", Sev::inform) << ents << endl;
1341
1342 treeSurfPtr = boost::shared_ptr<OrientedBoxTreeTool>(
1343 new OrientedBoxTreeTool(&getPostProcMesh(), "ROOTSETSURF", true));
1344 CHKERR treeSurfPtr->build(ents, rootSetSurf);
1345
1346 // CHKERR treeSurfPtr->stats(rootSetSurf, std::cerr);
1347
1349}
1350
1352
1354
1355 OpMoveNode(boost::shared_ptr<ContactTree> contact_tree_ptr,
1356 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1357 boost::shared_ptr<MatrixDouble> u_h1_ptr);
1358 MoFEMErrorCode doWork(int side, EntityType type,
1360
1361protected:
1362 boost::shared_ptr<ContactTree> contactTreePtr;
1363 boost::shared_ptr<MatrixDouble> uH1Ptr;
1364 boost::shared_ptr<ContactOps::CommonData> commonDataPtr;
1365};
1366
1368 boost::shared_ptr<ContactTree> contact_tree_ptr,
1369 boost::shared_ptr<ContactOps::CommonData> common_data_ptr,
1370 boost::shared_ptr<MatrixDouble> u_h1_ptr)
1371 : UOP(NOSPACE, UOP::OPSPACE), contactTreePtr(contact_tree_ptr),
1372 uH1Ptr(u_h1_ptr), commonDataPtr(common_data_ptr) {}
1373
1374MoFEMErrorCode OpMoveNode::doWork(int side, EntityType type,
1377
1378 auto get_body_id = [&](auto fe_ent) {
1379 for (auto &m : contactTreePtr->bodyMap) {
1380 if (m.second.find(fe_ent) != m.second.end()) {
1381 return m.first;
1382 }
1383 }
1384 return -1;
1385 };
1386
1387 auto &moab_post_proc_mesh = contactTreePtr->getPostProcMesh();
1388 auto &post_proc_ents = contactTreePtr->getPostProcElements();
1389
1390 auto fe_ent = getNumeredEntFiniteElementPtr()->getEnt();
1391 auto fe_id = id_from_handle(fe_ent);
1392 auto body_id = get_body_id(fe_ent);
1393 auto &map_gauss_pts = contactTreePtr->getMapGaussPts();
1394
1395 CHKERR moab_post_proc_mesh.tag_clear_data(contactTreePtr->thEleId,
1396 post_proc_ents, &fe_id);
1397 CHKERR moab_post_proc_mesh.tag_clear_data(contactTreePtr->thBodyId,
1398 post_proc_ents, &body_id);
1399
1400 auto nb_gauss_pts = getGaussPts().size2();
1401 auto t_u_h1 = getFTensor1FromMat<3>(*uH1Ptr);
1402 auto t_u_l2 = getFTensor1FromMat<3>(commonDataPtr->contactDisp);
1403 auto t_coords = getFTensor1CoordsAtGaussPts();
1404
1405 MatrixDouble x_h1(nb_gauss_pts, 3);
1406 auto t_x_h1 = getFTensor1FromPtr<3>(&x_h1(0, 0));
1407 MatrixDouble x_l2(nb_gauss_pts, 3);
1408 auto t_x_l2 = getFTensor1FromPtr<3>(&x_l2(0, 0));
1409 MatrixDouble tractions = trans(commonDataPtr->contactTraction);
1410 MatrixDouble coords = getCoordsAtGaussPts();
1411
1412 FTensor::Index<'i', 3> i;
1413
1414 // VectorDouble bases(nb_bases, 0);
1415 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1416 t_x_h1(i) = t_coords(i) + t_u_h1(i);
1417 t_x_l2(i) = t_coords(i) + t_u_l2(i);
1418
1419 ++t_coords;
1420 ++t_u_h1;
1421 ++t_u_l2;
1422 ++t_x_h1;
1423 ++t_x_l2;
1424 }
1425
1426 CHKERR moab_post_proc_mesh.set_coords(
1427 &*map_gauss_pts.begin(), map_gauss_pts.size(), &*x_h1.data().begin());
1428 CHKERR moab_post_proc_mesh.tag_set_data(
1429 contactTreePtr->thSmallX, &*map_gauss_pts.begin(), map_gauss_pts.size(),
1430 &*x_h1.data().begin());
1431 CHKERR moab_post_proc_mesh.tag_set_data(
1432 contactTreePtr->thLargeX, &*map_gauss_pts.begin(), map_gauss_pts.size(),
1433 &*coords.data().begin());
1434 CHKERR moab_post_proc_mesh.tag_set_data(
1435 contactTreePtr->thTraction, &*map_gauss_pts.begin(), map_gauss_pts.size(),
1436 &*tractions.data().begin());
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 ParallelComm *pcomm_post_proc_mesh = ParallelComm::get_pcomm(
1667 &(contactTreePtr->getPostProcMesh()), MYPCOMM_INDEX);
1668 int part;
1669 if (contactTreePtr->getPostProcMesh().tag_get_data(
1670 pcomm_post_proc_mesh->part_tag(), &face, 1, &part) == MB_SUCCESS) {
1671 return part;
1672 }
1673 return -1;
1674 };
1675
1676 auto check_face = [&](auto face, auto fe_id, auto part) {
1677 auto face_id = get_face_id(face);
1678 auto face_part = get_face_part(face);
1679 if (face_id == fe_id && face_part == part)
1680 return true;
1681 return false;
1682 };
1683
1684 // vectors
1685 VectorDouble3 v(3);
1686 FTensor::Tensor1<FTensor::PackPtr<double *, 0>, 3> t_v(&v[0], &v[1], &v[2]);
1687 auto save_vec_tag = [&](auto &th, auto &t_d, const int gg) {
1689 if (postProcMeshPtr) {
1690 t_v(i) = t_d(i);
1691 for (auto &a : v.data())
1692 a = set_float_precision(a);
1693 CHKERR postProcMeshPtr->tag_set_data(th, &(*mapGaussPtsPtr)[gg], 1,
1694 &*v.data().begin());
1695 }
1697 };
1698
1699 Tag th_mark = create_tag("contact_mark", 1);
1700 Tag th_mark_slave = create_tag("contact_mark_slave", 1);
1701 Tag th_body_id = create_tag("contact_body_id", 1);
1702 Tag th_gap = create_tag("contact_gap", 1);
1703 Tag th_tn_master = create_tag("contact_tn_master", 1);
1704 Tag th_tn_slave = create_tag("contact_tn_slave", 1);
1705 Tag th_contact_traction = create_tag("contact_traction", 3);
1706 Tag th_contact_traction_master = create_tag("contact_traction_master", 3);
1707 Tag th_contact_traction_slave = create_tag("contact_traction_slave", 3);
1708 Tag th_c = create_tag("contact_c", 1);
1709 Tag th_normal = create_tag("contact_normal", 3);
1710 Tag th_dist = create_tag("contact_dip", 3);
1711
1712 auto t_ele_centre = get_ele_centre(getFTensor1Coords());
1713 auto ele_radius = get_ele_radius(t_ele_centre, getFTensor1Coords());
1714
1715 contactTreePtr->shadowDataMap.clear();
1716 auto &shadow_vec = contactTreePtr->shadowDataMap[fe_ent];
1717 shadow_vec.clear();
1718
1719 auto adj_fe_ids = get_fe_adjacencies(fe_ent);
1720
1721 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1722
1723 FTensor::Tensor1<double, 3> t_spatial_coords;
1724 t_spatial_coords(i) = t_coords(i) + t_disp_h1(i);
1725
1726 if (postProcMeshPtr) {
1727 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1728 }
1729
1730 auto faces_close = get_closet_face(&t_spatial_coords(0), ele_radius);
1731 for (auto face_close : faces_close) {
1732 if (check_face(face_close, fe_id, m_field.get_comm_rank())) {
1733
1734 auto body_id = get_body_id(face_close);
1735
1736 auto master_face_conn = get_face_conn(face_close);
1737 std::array<double, 9> master_coords;
1738 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1739 contactTreePtr->thSmallX, master_face_conn, 3,
1740 master_coords.data());
1741 std::array<double, 9> master_traction;
1742 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1743 contactTreePtr->thTraction, master_face_conn, 3,
1744 master_traction.data());
1745 auto t_normal_face_close = get_normal(master_coords);
1746 t_normal_face_close.normalize();
1747
1748 if (postProcMeshPtr) {
1749 double m = 1;
1750 CHKERR save_scal_tag(th_mark, m, gg);
1751 CHKERR save_scal_tag(th_body_id, static_cast<double>(body_id), gg);
1752 CHKERR save_vec_tag(th_normal, t_normal_face_close, gg);
1753 CHKERR save_vec_tag(th_contact_traction, t_traction, gg);
1754 }
1755
1756 FTensor::Tensor1<double, 3> t_unit_ray;
1757 t_unit_ray(i) = -t_normal_face_close(i);
1758 FTensor::Tensor1<double, 3> t_ray_point;
1759 t_ray_point(i) =
1760 t_spatial_coords(i) -
1761 t_unit_ray(i) * ContactOps::airplane_ray_distance * ele_radius;
1762
1763 constexpr double eps = 1e-3;
1764 auto [faces_out, faces_dist] =
1765 get_faces_out(&t_ray_point(0), &t_unit_ray(0),
1766 2 * ContactOps::airplane_ray_distance * ele_radius,
1767 eps * ele_radius);
1768
1769 auto m = make_map(faces_out, faces_dist, t_ray_point, t_unit_ray,
1770 t_spatial_coords);
1771 for (auto m_it = m.begin(); m_it != m.end(); ++m_it) {
1772 auto face = m_it->second;
1773 if (face != face_close) {
1774
1775 if (
1776
1777 (adj_fe_ids.find(get_face_id(face)) == adj_fe_ids.end() ||
1778 get_face_part(face) != m_field.get_comm_rank())
1779
1780 ) {
1781
1782 shadow_vec.push_back(ContactTree::FaceData());
1783 shadow_vec.back().gaussPtNb = gg;
1784
1785 auto slave_face_conn = get_face_conn(face);
1786 std::array<double, 9> slave_coords;
1787 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1788 contactTreePtr->thSmallX, slave_face_conn, 3,
1789 slave_coords.data());
1790 auto t_normal_face = get_normal(slave_coords);
1791 std::array<double, 9> slave_tractions;
1792 CHKERR contactTreePtr->getPostProcMesh().tag_get_data(
1793 contactTreePtr->thTraction, slave_face_conn, 3,
1794 slave_tractions.data());
1795
1796 auto t_master_point =
1797 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data());
1798 auto t_slave_point =
1799 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data());
1800 auto t_ray_point_data =
1801 getFTensor1FromPtr<3>(shadow_vec.back().rayPoint.data());
1802 auto t_unit_ray_data =
1803 getFTensor1FromPtr<3>(shadow_vec.back().unitRay.data());
1804
1805 t_slave_point(i) = t_ray_point(i) + m_it->first * t_unit_ray(i);
1806
1807 auto eval_position = [&](auto &&t_elem_coords, auto &&t_point) {
1808 std::array<double, 2> loc_coords;
1811 &t_elem_coords(0, 0), &t_point(0), 1,
1812 loc_coords.data()),
1813 "get local coords");
1814 FTensor::Tensor1<double, 3> t_shape_fun;
1815 CHK_THROW_MESSAGE(Tools::shapeFunMBTRI<0>(&t_shape_fun(0),
1816 &loc_coords[0],
1817 &loc_coords[1], 1),
1818 "calc shape fun");
1819 FTensor::Index<'i', 3> i;
1820 FTensor::Index<'j', 3> j;
1821 FTensor::Tensor1<double, 3> t_point_out;
1822 t_point_out(i) = t_shape_fun(j) * t_elem_coords(j, i);
1823 return t_point_out;
1824 };
1825
1826 auto t_master_point_updated = eval_position(
1827 getFTensor2FromPtr<3, 3>(master_coords.data()),
1828 getFTensor1FromPtr<3>(shadow_vec.back().masterPoint.data()));
1829 t_master_point(i) = t_master_point_updated(i);
1830
1831 auto t_slave_point_updated = eval_position(
1832 getFTensor2FromPtr<3, 3>(slave_coords.data()),
1833 getFTensor1FromPtr<3>(shadow_vec.back().slavePoint.data()));
1834 t_slave_point(i) = t_slave_point_updated(i);
1835
1836 t_ray_point_data(i) = t_ray_point(i);
1837 t_unit_ray_data(i) = t_unit_ray(i);
1838
1839 std::copy(master_coords.begin(), master_coords.end(),
1840 shadow_vec.back().masterPointNodes.begin());
1841 std::copy(master_traction.begin(), master_traction.end(),
1842 shadow_vec.back().masterTractionNodes.begin());
1843 std::copy(slave_coords.begin(), slave_coords.end(),
1844 shadow_vec.back().slavePointNodes.begin());
1845 std::copy(slave_tractions.begin(), slave_tractions.end(),
1846 shadow_vec.back().slaveTractionNodes.begin());
1847
1848 shadow_vec.back().eleRadius = ele_radius;
1849
1850 // CHKERR get_tag_data(contactTreePtr->thIds, face,
1851 // shadow_vec.back().dofsSlaveIds);
1852 // CHKERR get_tag_data(contactTreePtr->thCoeff, face,
1853 // shadow_vec.back().dofsSlaveCoeff);
1854 // CHKERR get_tag_data(contactTreePtr->thBases, face,
1855 // shadow_vec.back().baseSlaveFuncs);
1856
1857 if (postProcMeshPtr) {
1858 auto [gap, tn_master, tn_slave, c, t_master_traction,
1859 t_slave_traction] =
1860 multiGetGap(&(shadow_vec.back()), t_spatial_coords);
1862 t_gap_vec(i) = t_slave_point(i) - t_spatial_coords(i);
1863 CHKERR save_scal_tag(th_gap, gap, gg);
1864 CHKERR save_scal_tag(th_tn_master, tn_master, gg);
1865 CHKERR save_scal_tag(th_tn_slave, tn_slave, gg);
1866 CHKERR save_scal_tag(th_c, c, gg);
1867 double m = 1;
1868 CHKERR save_scal_tag(th_mark_slave, m, gg);
1869 CHKERR save_vec_tag(th_dist, t_gap_vec, gg);
1870 CHKERR save_vec_tag(th_contact_traction_master,
1871 t_master_traction, gg);
1872 CHKERR save_vec_tag(th_contact_traction_slave, t_slave_traction,
1873 gg);
1874 }
1875
1876 break;
1877 }
1878 }
1879 }
1880 break;
1881 }
1882 }
1883 next();
1884 }
1885
1887}
1888
1890 const std::string block_name, int dim) {
1891 Range r;
1892
1893 auto mesh_mng = m_field.getInterface<MeshsetsManager>();
1894 auto bcs = mesh_mng->getCubitMeshsetPtr(
1895
1896 std::regex((boost::format("%s(.*)") % block_name).str())
1897
1898 );
1899
1900 for (auto bc : bcs) {
1901 Range faces;
1902 CHK_MOAB_THROW(bc->getMeshsetIdEntitiesByDimension(m_field.get_moab(), dim,
1903 faces, true),
1904 "get meshset ents");
1905 r.merge(faces);
1906 }
1907
1908 return r;
1909};
1910
1911boost::shared_ptr<ForcesAndSourcesCore>
1913
1914 auto &m_field = ep.mField;
1915
1916 boost::shared_ptr<ContactTree> fe_contact_tree;
1917
1918 auto impl = [&]() {
1920
1921 /** Contact requires that body is marked */
1922 auto get_body_range = [&](auto name, int dim, auto sev) {
1923 std::map<int, Range> map;
1924
1925 for (auto m_ptr :
1926 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
1927 std::regex(
1928
1929 (boost::format("%s(.*)") % name).str()
1930
1931 ))
1932
1933 ) {
1934 Range ents;
1935 CHK_MOAB_THROW(m_ptr->getMeshsetIdEntitiesByDimension(
1936 m_field.get_moab(), dim, ents, true),
1937 "by dim");
1938 map[m_ptr->getMeshsetId()] = ents;
1939 MOFEM_LOG("EPSYNC", sev) << "Meshset: " << m_ptr->getMeshsetId() << " "
1940 << ents.size() << " entities";
1941 }
1942
1943 MOFEM_LOG_SEVERITY_SYNC(m_field.get_comm(), sev);
1944 return map;
1945 };
1946
1947 auto get_map_skin = [&](auto &&map) {
1948 ParallelComm *pcomm =
1949 ParallelComm::get_pcomm(&m_field.get_moab(), MYPCOMM_INDEX);
1950
1951 Skinner skin(&m_field.get_moab());
1952 for (auto &m : map) {
1953 Range skin_faces;
1954 CHKERR skin.find_skin(0, m.second, false, skin_faces);
1955 CHK_MOAB_THROW(pcomm->filter_pstatus(
1956 skin_faces, PSTATUS_SHARED | PSTATUS_MULTISHARED,
1957 PSTATUS_NOT, -1, nullptr),
1958 "filter");
1959 m.second.swap(skin_faces);
1960 }
1961 return map;
1962 };
1963
1964 /* The above code is written in C++ and it appears to be defining and using
1965 various operations on boundary elements and side elements. */
1967 using BoundaryEleOp = BoundaryEle::UserDataOperator;
1968
1969 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1970
1971 auto calcs_side_traction = [&](auto &pip) {
1973 using EleOnSide =
1975 using SideEleOp = EleOnSide::UserDataOperator;
1976 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1977 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
1978 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1979 boost::make_shared<CGGUserPolynomialBase>();
1980 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
1981 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
1983 op_loop_domain_side->getOpPtrVector().push_back(
1985 ep.piolaStress, contact_common_data_ptr->contactTractionPtr(),
1986 boost::make_shared<double>(1.0)));
1987 pip.push_back(op_loop_domain_side);
1989 };
1990
1991 auto add_contact_three = [&]() {
1993 auto tree_moab_ptr = boost::make_shared<moab::Core>();
1994 fe_contact_tree = boost::make_shared<ContactTree>(
1995 m_field, tree_moab_ptr, ep.spaceOrder,
1996 get_body_range("CONTACT", SPACE_DIM - 1, Sev::inform));
1997 fe_contact_tree->getOpPtrVector().push_back(
1999 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2000 CHKERR calcs_side_traction(fe_contact_tree->getOpPtrVector());
2001 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2002 fe_contact_tree->getOpPtrVector().push_back(
2004 fe_contact_tree->getOpPtrVector().push_back(
2005 new OpMoveNode(fe_contact_tree, contact_common_data_ptr, u_h1_ptr));
2007 };
2008
2009 CHKERR add_contact_three();
2010
2012 };
2013
2014 CHK_THROW_MESSAGE(impl(), "createContactDetectionFiniteElement");
2015
2016 struct exclude_sdf {
2017 exclude_sdf(Range &&r) : map(r) {}
2018 bool operator()(FEMethod *fe_method_ptr) {
2019 auto ent = fe_method_ptr->getFEEntityHandle();
2020 if (map.find(ent) != map.end()) {
2021 return false;
2022 }
2023 return true;
2024 }
2025
2026 private:
2027 Range map;
2028 };
2029
2030 fe_contact_tree->exeTestHook =
2031 exclude_sdf(get_range_from_block(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2032
2033 return fe_contact_tree;
2034};
2035
2036static auto get_body_range(MoFEM::Interface &m_field, const std::string name,
2037 int dim) {
2038 std::map<int, Range> map;
2039
2040 for (auto m_ptr :
2041 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
2042
2043 (boost::format("%s(.*)") % name).str()
2044
2045 ))
2046
2047 ) {
2048 Range ents;
2049 CHK_MOAB_THROW(m_ptr->getMeshsetIdEntitiesByDimension(m_field.get_moab(),
2050 dim, ents, true),
2051 "by dim");
2052 map[m_ptr->getMeshsetId()] = ents;
2053 }
2054
2055 return map;
2056};
2057
2059 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2060 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2062
2063 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2064
2065 auto &m_field = ep.mField;
2066
2067 using BoundaryEle =
2069 using EleOnSide =
2071 using SideEleOp = EleOnSide::UserDataOperator;
2072 using BdyEleOp = BoundaryEle::UserDataOperator;
2073
2074 // First: Iterate over skeleton FEs adjacent to Domain FEs
2075 // Note: BoundaryEle, i.e. uses skeleton interation rule
2076 auto op_loop_skeleton_side = new OpLoopSide<BoundaryEle>(
2077 m_field, ep.contactElement, SPACE_DIM - 1, Sev::noisy);
2078
2079 auto rule_contact = [](int, int, int o) { return -1; };
2081
2082 auto set_rule_contact = [refine](
2083
2084 ForcesAndSourcesCore *fe_raw_ptr, int order_row,
2085 int order_col, int order_data
2086
2087 ) {
2089 auto rule = 2 * order_data;
2090 fe_raw_ptr->gaussPts = Tools::refineTriangleIntegrationPts(rule, refine);
2092 };
2093
2094 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2095 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2096 CHKERR
2097 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2098 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
2099 ep.frontAdjEdges);
2100
2101 // Second: Iterate over domain FEs adjacent to skelton, particularly
2102 // one domain element.
2103 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2104
2105 // Data storing contact fields
2106 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2107
2108 auto add_ops_domain_side = [&](auto &pip) {
2110 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
2111 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
2112 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
2113 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2114 boost::make_shared<CGGUserPolynomialBase>();
2115 CHKERR
2116 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
2117 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
2119 op_loop_domain_side->getOpPtrVector().push_back(
2121 broken_data_ptr));
2122 op_loop_domain_side->getOpPtrVector().push_back(
2124 ep.piolaStress, contact_common_data_ptr->contactTractionPtr()));
2125 pip.push_back(op_loop_domain_side);
2127 };
2128
2129 auto add_ops_contact_rhs = [&](auto &pip) {
2131 // get body id and SDF range
2132 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2133 get_body_range(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2134
2135 pip.push_back(new OpCalculateVectorFieldValues<3>(
2136 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2137 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2138 pip.push_back(
2140 pip.push_back(new OpTreeSearch(
2141 contact_tree_ptr, u_h1_ptr,
2142 contact_common_data_ptr->contactTractionPtr(),
2143 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1), nullptr,
2144 nullptr));
2146 ep.contactDisp, contact_common_data_ptr, contact_tree_ptr,
2147 contact_sfd_map_range_ptr));
2149 broken_data_ptr, contact_common_data_ptr, contact_tree_ptr));
2150
2152 };
2153
2154 // push ops to face/side pipeline
2155 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2156 CHKERR add_ops_contact_rhs(op_loop_skeleton_side->getOpPtrVector());
2157
2158 // Add skeleton to domain pipeline
2159 pip.push_back(op_loop_skeleton_side);
2160
2162};
2163
2165 EshelbianCore &ep, boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2166 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
2168
2169 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2170 auto &m_field = ep.mField;
2171
2172 using BoundaryEle =
2174 using EleOnSide =
2176 using SideEleOp = EleOnSide::UserDataOperator;
2177 using BdyEleOp = BoundaryEle::UserDataOperator;
2178
2179 // First: Iterate over skeleton FEs adjacent to Domain FEs
2180 // Note: BoundaryEle, i.e. uses skeleton interation rule
2181 auto op_loop_skeleton_side = new OpLoopSide<BoundaryEle>(
2182 m_field, ep.contactElement, SPACE_DIM - 1, Sev::noisy);
2183
2184 auto rule_contact = [](int, int, int o) { return -1; };
2186
2187 auto set_rule_contact = [refine](
2188
2189 ForcesAndSourcesCore *fe_raw_ptr, int order_row,
2190 int order_col, int order_data
2191
2192 ) {
2194 auto rule = 2 * order_data;
2195 fe_raw_ptr->gaussPts = Tools::refineTriangleIntegrationPts(rule, refine);
2197 };
2198
2199 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = rule_contact;
2200 op_loop_skeleton_side->getSideFEPtr()->setRuleHook = set_rule_contact;
2201 CHKERR
2202 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
2203 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
2204 ep.frontAdjEdges);
2205
2206 // Second: Iterate over domain FEs adjacent to skelton, particularly
2207 // one domain element.
2208 auto broken_data_ptr = boost::make_shared<std::vector<BrokenBaseSideData>>();
2209
2210 // Data storing contact fields
2211 auto contact_common_data_ptr = boost::make_shared<ContactOps::CommonData>();
2212
2213 auto add_ops_domain_side = [&](auto &pip) {
2215 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
2216 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
2217 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
2218 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
2219 boost::make_shared<CGGUserPolynomialBase>();
2220 CHKERR
2221 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
2222 op_loop_domain_side->getOpPtrVector(), {HDIV, H1, L2},
2224 op_loop_domain_side->getOpPtrVector().push_back(
2226 broken_data_ptr));
2227 op_loop_domain_side->getOpPtrVector().push_back(
2229 ep.piolaStress, contact_common_data_ptr->contactTractionPtr()));
2230 pip.push_back(op_loop_domain_side);
2232 };
2233
2234 auto add_ops_contact_lhs = [&](auto &pip) {
2236 pip.push_back(new OpCalculateVectorFieldValues<3>(
2237 ep.contactDisp, contact_common_data_ptr->contactDispPtr()));
2238 auto u_h1_ptr = boost::make_shared<MatrixDouble>();
2239 pip.push_back(
2241 pip.push_back(new OpTreeSearch(
2242 contact_tree_ptr, u_h1_ptr,
2243 contact_common_data_ptr->contactTractionPtr(),
2244 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1), nullptr,
2245 nullptr));
2246
2247 // get body id and SDF range
2248 auto contact_sfd_map_range_ptr = boost::make_shared<std::map<int, Range>>(
2249 get_body_range(m_field, "CONTACT_SDF", SPACE_DIM - 1));
2250
2252 ep.contactDisp, ep.contactDisp, contact_common_data_ptr,
2253 contact_tree_ptr, contact_sfd_map_range_ptr));
2254 pip.push_back(
2256 ep.contactDisp, broken_data_ptr, contact_common_data_ptr,
2257 contact_tree_ptr, contact_sfd_map_range_ptr));
2258 pip.push_back(
2260 broken_data_ptr, ep.contactDisp, contact_common_data_ptr,
2261 contact_tree_ptr));
2262
2264 };
2265
2266 // push ops to face/side pipeline
2267 CHKERR add_ops_domain_side(op_loop_skeleton_side->getOpPtrVector());
2268 CHKERR add_ops_contact_lhs(op_loop_skeleton_side->getOpPtrVector());
2269
2270 // Add skeleton to domain pipeline
2271 pip.push_back(op_loop_skeleton_side);
2272
2274};
2275
2278 boost::shared_ptr<ForcesAndSourcesCore> fe_ptr,
2279 boost::shared_ptr<MatrixDouble> u_h1_ptr,
2280 boost::shared_ptr<MatrixDouble> contact_traction_ptr,
2281 Range r, moab::Interface *post_proc_mesh_ptr,
2282 std::vector<EntityHandle> *map_gauss_pts_ptr) {
2283
2284 auto &m_field = ep.mField;
2285 auto contact_tree_ptr = boost::dynamic_pointer_cast<ContactTree>(fe_ptr);
2286 return new OpTreeSearch(
2287 contact_tree_ptr, u_h1_ptr, contact_traction_ptr,
2288 get_range_from_block(m_field, "CONTACT", SPACE_DIM - 1),
2289 post_proc_mesh_ptr, map_gauss_pts_ptr);
2290}
2291
2292} // 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
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
OpBaseImpl< PETSC, EdgeEleOp > OpBase
Definition radiation.cpp:29
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)
Integrate grad-grad operator.
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)
Class dedicated to integrate operator.
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.
VectorDouble locF
local entity vector
int nbRows
number of dofs on rows
MatrixDouble locMat
local entity block matrix
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