v0.15.5
Loading...
Searching...
No Matches
ContactOps.hpp
Go to the documentation of this file.
1
2
3/**
4 * \file ContactOps.hpp
5 * \example mofem/tutorials/adv-1_contact/src/ContactOps.hpp
6 */
7
8#ifndef __CONTACTOPS_HPP__
9#define __CONTACTOPS_HPP__
10
11namespace ContactOps {
12
13//! [Common data]
14struct CommonData : public boost::enable_shared_from_this<CommonData> {
15 // MatrixDouble contactStress;
16 MatrixDouble contactTraction;
17 MatrixDouble contactDisp;
18 MatrixDouble contactDispGrad;
19
20 VectorDouble sdfVals; ///< size is equal to number of gauss points on element
21 MatrixDouble gradsSdf; ///< nb of rows is equals to dimension, and nb of cols
22 ///< is equals to number of gauss points on element
23 MatrixDouble hessSdf; ///< nb of rows is equals to nb of element of symmetric
24 ///< matrix, and nb of cols is equals to number of gauss
25 ///< points on element
26 VectorDouble constraintVals;
27
28 static SmartPetscObj<Vec>
29 totalTraction; // User have to release and create vector when appropiate.
30
31 static auto createTotalTraction(MoFEM::Interface &m_field) {
32 constexpr int ghosts[] = {0, 1, 2, 3, 4};
34 createGhostVector(m_field.get_comm(),
35
36 (m_field.get_comm_rank() == 0) ? 5 : 0, 5,
37
38 (m_field.get_comm_rank() == 0) ? 0 : 5, ghosts);
39 return totalTraction;
40 }
41
44 const double *t_ptr;
45 CHK_THROW_MESSAGE(VecGetArrayRead(CommonData::totalTraction, &t_ptr),
46 "get array");
47 FTensor::Tensor1<double, 5> t{t_ptr[0], t_ptr[1], t_ptr[2], t_ptr[3],
48 t_ptr[4]};
49 CHK_THROW_MESSAGE(VecRestoreArrayRead(CommonData::totalTraction, &t_ptr),
50 "restore array");
51 return t;
52 } else {
53 return FTensor::Tensor1<double, 5>{0., 0., 0., 0., 0.};
54 }
55 }
56
57 inline auto contactTractionPtr() {
58 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
60 }
61
62 inline auto contactDispPtr() {
63 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &contactDisp);
64 }
65
66 inline auto contactDispGradPtr() {
67 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
69 }
70
71 inline auto sdfPtr() {
72 return boost::shared_ptr<VectorDouble>(shared_from_this(), &sdfVals);
73 }
74
75 inline auto gradSdfPtr() {
76 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &gradsSdf);
77 }
78
79 inline auto hessSdfPtr() {
80 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &hessSdf);
81 }
82
83 inline auto constraintPtr() {
84 return boost::shared_ptr<VectorDouble>(shared_from_this(), &constraintVals);
85 }
86};
87
88SmartPetscObj<Vec> CommonData::totalTraction;
89
90//! [Common data]
91
92//! [Surface distance function from python]
93#ifdef ENABLE_PYTHON_BINDING
94struct SDFPython {
95 SDFPython() = default;
96 virtual ~SDFPython() = default;
97
98 MoFEMErrorCode sdfInit(const std::string py_file) {
100 try {
101
102 // create main module
103 auto main_module = bp::import("__main__");
104 mainNamespace = main_module.attr("__dict__");
105 bp::exec_file(py_file.c_str(), mainNamespace, mainNamespace);
106 // create a reference to python function
107 sdfFun = mainNamespace["sdf"];
108 sdfGradFun = mainNamespace["grad_sdf"];
109 sdfHessFun = mainNamespace["hess_sdf"];
110
111 } catch (bp::error_already_set const &) {
112 // print all other errors to stderr
113 PyErr_Print();
115 }
117 };
118
119 template <typename T>
120 inline std::vector<T>
121 py_list_to_std_vector(const boost::python::object &iterable) {
122 return std::vector<T>(boost::python::stl_input_iterator<T>(iterable),
123 boost::python::stl_input_iterator<T>());
124 }
125
126 MoFEMErrorCode evalSdf(
127
128 double delta_t, double t, np::ndarray x, np::ndarray y, np::ndarray z,
129 np::ndarray tx, np::ndarray ty, np::ndarray tz, int block_id,
130 np::ndarray &sdf
131
132 ) {
134 try {
135
136 // call python function
137 sdf = bp::extract<np::ndarray>(
138 sdfFun(delta_t, t, x, y, z, tx, ty, tz, block_id));
139
140 } catch (bp::error_already_set const &) {
141 // print all other errors to stderr
142 PyErr_Print();
144 }
146 }
147
148 MoFEMErrorCode evalGradSdf(
149
150 double delta_t, double t, np::ndarray x, np::ndarray y, np::ndarray z,
151 np::ndarray tx, np::ndarray ty, np::ndarray tz, int block_id,
152 np::ndarray &grad_sdf
153
154 ) {
156 try {
157
158 // call python function
159 grad_sdf = bp::extract<np::ndarray>(
160 sdfGradFun(delta_t, t, x, y, z, tx, ty, tz, block_id));
161
162 } catch (bp::error_already_set const &) {
163 // print all other errors to stderr
164 PyErr_Print();
166 }
168 }
169
170 MoFEMErrorCode evalHessSdf(
171
172 double delta_t, double t, np::ndarray x, np::ndarray y, np::ndarray z,
173 np::ndarray tx, np::ndarray ty, np::ndarray tz, int block_id,
174 np::ndarray &hess_sdf
175
176 ) {
178 try {
179
180 // call python function
181 hess_sdf = bp::extract<np::ndarray>(
182 sdfHessFun(delta_t, t, x, y, z, tx, ty, tz, block_id));
183
184 } catch (bp::error_already_set const &) {
185 // print all other errors to stderr
186 PyErr_Print();
188 }
190 }
191
192private:
193 bp::object mainNamespace;
194 bp::object sdfFun;
195 bp::object sdfGradFun;
196 bp::object sdfHessFun;
197};
198
199static boost::weak_ptr<SDFPython> sdfPythonWeakPtr;
200
201inline np::ndarray convert_to_numpy(VectorDouble &data, int nb_gauss_pts,
202 int id) {
203 auto dtype = np::dtype::get_builtin<double>();
204 auto size = bp::make_tuple(nb_gauss_pts);
205 auto stride = bp::make_tuple(3 * sizeof(double));
206 return (np::from_data(&data[id], dtype, size, stride, bp::object()));
207};
208#endif
209//! [Surface distance function from python]
210
211using SurfaceDistanceFunction = boost::function<VectorDouble(
212 double delta_t, double t, int nb_gauss_pts, MatrixDouble &spatial_coords,
213 MatrixDouble &normals_at_pts, int block_id)>;
214
215using GradSurfaceDistanceFunction = boost::function<MatrixDouble(
216 double delta_t, double t, int nb_gauss_pts, MatrixDouble &spatial_coords,
217 MatrixDouble &normals_at_pts, int block_id)>;
218
219using HessSurfaceDistanceFunction = boost::function<MatrixDouble(
220 double delta_t, double t, int nb_gauss_pts, MatrixDouble &spatial_coords,
221 MatrixDouble &normals_at_pts, int block_id)>;
222
223inline VectorDouble surface_distance_function(double delta_t, double t,
224 int nb_gauss_pts,
225 MatrixDouble &m_spatial_coords,
226 MatrixDouble &m_normals_at_pts,
227 int block_id) {
228
229#ifdef ENABLE_PYTHON_BINDING
230 if (auto sdf_ptr = sdfPythonWeakPtr.lock()) {
231
232 VectorDouble v_spatial_coords = m_spatial_coords.data();
233 VectorDouble v_normal_at_pts = m_normals_at_pts.data();
234
235 bp::list python_coords;
236 bp::list python_normals;
237
238 for (int idx = 0; idx < 3; ++idx) {
239 python_coords.append(
240 convert_to_numpy(v_spatial_coords, nb_gauss_pts, idx));
241 python_normals.append(
242 convert_to_numpy(v_normal_at_pts, nb_gauss_pts, idx));
243 }
244
245 np::ndarray np_sdf = np::empty(bp::make_tuple(nb_gauss_pts),
246 np::dtype::get_builtin<double>());
247 CHK_MOAB_THROW(sdf_ptr->evalSdf(delta_t, t,
248 bp::extract<np::ndarray>(python_coords[0]),
249 bp::extract<np::ndarray>(python_coords[1]),
250 bp::extract<np::ndarray>(python_coords[2]),
251 bp::extract<np::ndarray>(python_normals[0]),
252 bp::extract<np::ndarray>(python_normals[1]),
253 bp::extract<np::ndarray>(python_normals[2]),
254 block_id, np_sdf),
255 "Failed python call");
256
257 // check the shape of returned array
258 if (np_sdf.get_nd() != 1 || np_sdf.get_shape()[0] != nb_gauss_pts) {
261 "Wrong number of dimensions or size of SDF returned from "
262 "python, expected: dim 1, got: " +
263 std::to_string(np_sdf.get_nd()) + ", expected size: (" +
264 std::to_string(nb_gauss_pts) + ")");
265 }
266
267 double *sdf_val_ptr = reinterpret_cast<double *>(np_sdf.get_data());
268
269 VectorDouble v_sdf;
270 v_sdf.resize(nb_gauss_pts, false);
271
272 for (size_t gg = 0; gg < nb_gauss_pts; ++gg)
273 v_sdf[gg] = *(sdf_val_ptr + gg);
274
275 return v_sdf;
276 }
277#endif
278 VectorDouble v_sdf;
279 v_sdf.resize(nb_gauss_pts, false);
280 auto t_coords = getFTensor1FromPtr<3>(&m_spatial_coords(0, 0));
281
282 for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
283 v_sdf[gg] = -t_coords(2) - 0.1;
284 ++t_coords;
285 }
286
287 return v_sdf;
288}
289
290inline MatrixDouble
291grad_surface_distance_function(double delta_t, double t, int nb_gauss_pts,
292 MatrixDouble &m_spatial_coords,
293 MatrixDouble &m_normals_at_pts, int block_id) {
294#ifdef ENABLE_PYTHON_BINDING
295 if (auto sdf_ptr = sdfPythonWeakPtr.lock()) {
296
297 VectorDouble v_spatial_coords = m_spatial_coords.data();
298 VectorDouble v_normal_at_pts = m_normals_at_pts.data();
299
300 bp::list python_coords;
301 bp::list python_normals;
302
303 for (int idx = 0; idx < 3; ++idx) {
304 python_coords.append(
305 convert_to_numpy(v_spatial_coords, nb_gauss_pts, idx));
306 python_normals.append(
307 convert_to_numpy(v_normal_at_pts, nb_gauss_pts, idx));
308 }
309
310 np::ndarray np_grad_sdf = np::empty(bp::make_tuple(nb_gauss_pts, 3),
311 np::dtype::get_builtin<double>());
312 CHK_MOAB_THROW(sdf_ptr->evalGradSdf(
313 delta_t, t, bp::extract<np::ndarray>(python_coords[0]),
314 bp::extract<np::ndarray>(python_coords[1]),
315 bp::extract<np::ndarray>(python_coords[2]),
316 bp::extract<np::ndarray>(python_normals[0]),
317 bp::extract<np::ndarray>(python_normals[1]),
318 bp::extract<np::ndarray>(python_normals[2]), block_id,
319 np_grad_sdf),
320 "Failed python call");
321
322 // check the shape of returned array
323 if (np_grad_sdf.get_shape()[0] != nb_gauss_pts ||
324 np_grad_sdf.get_shape()[1] != 3) {
326 "Wrong shape of gradient of SDF returned from "
327 "python, expected: (" +
328 std::to_string(nb_gauss_pts) + ", 3), got: (" +
329 std::to_string(np_grad_sdf.get_shape()[0]) + ", " +
330 std::to_string(np_grad_sdf.get_shape()[1]) + ")");
331 }
332
333 double *grad_ptr = reinterpret_cast<double *>(np_grad_sdf.get_data());
334
335 MatrixDouble m_grad_sdf;
336 m_grad_sdf.resize(3, nb_gauss_pts, false);
337 for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
338 for (int idx = 0; idx < 3; ++idx)
339 m_grad_sdf(idx, gg) = *(grad_ptr + (3 * gg + idx));
340 }
341 return m_grad_sdf;
342 }
343#endif
344 MatrixDouble m_grad_sdf;
345 m_grad_sdf.resize(3, nb_gauss_pts, false);
346 FTensor::Index<'i', 3> i;
347 FTensor::Tensor1<double, 3> t_grad_sdf_set{0.0, 0.0, -1.0};
348 auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
349
350 for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
351 t_grad_sdf(i) = t_grad_sdf_set(i);
352 ++t_grad_sdf;
353 }
354
355 return m_grad_sdf;
356}
357
358inline MatrixDouble
359hess_surface_distance_function(double delta_t, double t, int nb_gauss_pts,
360 MatrixDouble &m_spatial_coords,
361 MatrixDouble &m_normals_at_pts, int block_id) {
362#ifdef ENABLE_PYTHON_BINDING
363 if (auto sdf_ptr = sdfPythonWeakPtr.lock()) {
364
365 VectorDouble v_spatial_coords = m_spatial_coords.data();
366 VectorDouble v_normal_at_pts = m_normals_at_pts.data();
367
368 bp::list python_coords;
369 bp::list python_normals;
370
371 for (int idx = 0; idx < 3; ++idx) {
372 python_coords.append(
373 convert_to_numpy(v_spatial_coords, nb_gauss_pts, idx));
374 python_normals.append(
375 convert_to_numpy(v_normal_at_pts, nb_gauss_pts, idx));
376 };
377
378 np::ndarray np_hess_sdf = np::empty(bp::make_tuple(nb_gauss_pts, 6),
379 np::dtype::get_builtin<double>());
380 CHK_MOAB_THROW(sdf_ptr->evalHessSdf(
381 delta_t, t, bp::extract<np::ndarray>(python_coords[0]),
382 bp::extract<np::ndarray>(python_coords[1]),
383 bp::extract<np::ndarray>(python_coords[2]),
384 bp::extract<np::ndarray>(python_normals[0]),
385 bp::extract<np::ndarray>(python_normals[1]),
386 bp::extract<np::ndarray>(python_normals[2]), block_id,
387 np_hess_sdf),
388 "Failed python call");
389
390 // check the shape of returned array
391 if (np_hess_sdf.get_shape()[0] != nb_gauss_pts ||
392 np_hess_sdf.get_shape()[1] != 6) {
394 "Wrong shape of Hessian of SDF returned from "
395 "python, expected: (" +
396 std::to_string(nb_gauss_pts) + ", 6), got: (" +
397 std::to_string(np_hess_sdf.get_shape()[0]) + ", " +
398 std::to_string(np_hess_sdf.get_shape()[1]) + ")");
399 }
400
401 double *hess_ptr = reinterpret_cast<double *>(np_hess_sdf.get_data());
402
403 MatrixDouble m_hess_sdf;
404 m_hess_sdf.resize(6, nb_gauss_pts, false);
405 for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
406 for (int idx = 0; idx < 6; ++idx)
407 m_hess_sdf(idx, gg) =
408 *(hess_ptr + (6 * gg + idx));
409 }
410 return m_hess_sdf;
411 }
412#endif
413 MatrixDouble m_hess_sdf;
414 m_hess_sdf.resize(6, nb_gauss_pts, false);
415 FTensor::Index<'i', 3> i;
416 FTensor::Index<'j', 3> j;
417 FTensor::Tensor2_symmetric<double, 3> t_hess_sdf_set{0., 0., 0., 0., 0., 0.};
418 auto t_hess_sdf = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
419
420 for (size_t gg = 0; gg < nb_gauss_pts; ++gg) {
421 t_hess_sdf(i, j) = t_hess_sdf_set(i, j);
422 ++t_hess_sdf;
423 }
424 return m_hess_sdf;
425}
426
427template <int DIM, IntegrationType I, typename BoundaryEleOp>
429
430template <int DIM, IntegrationType I, typename BoundaryEleOp>
432
433template <int DIM, IntegrationType I, typename BoundaryEleOp>
435
436template <int DIM, IntegrationType I, typename AssemblyBoundaryEleOp>
438
439template <int DIM, IntegrationType I, typename AssemblyBoundaryEleOp>
441
442template <int DIM, IntegrationType I, typename AssemblyBoundaryEleOp>
444
445template <typename T1, typename T2, int DIM1, int DIM2>
448 size_t nb_gauss_pts) {
449 MatrixDouble m_spatial_coords(nb_gauss_pts, 3);
450 m_spatial_coords.clear();
451 auto t_spatial_coords = getFTensor1FromPtr<3>(&m_spatial_coords(0, 0));
452 FTensor::Index<'i', DIM2> i;
453 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
454 t_spatial_coords(i) = t_coords(i) + t_disp(i);
455 ++t_spatial_coords;
456 ++t_coords;
457 ++t_disp;
458 }
459 return m_spatial_coords;
460}
461
462template <typename T1, int DIM1>
464 size_t nb_gauss_pts) {
465 MatrixDouble m_normals_at_pts(3, nb_gauss_pts);
466 m_normals_at_pts.clear();
467 FTensor::Index<'i', DIM1> i;
468 auto t_set_normal = getFTensor1FromMat<3>(m_normals_at_pts);
469 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
470 t_set_normal(i) = t_normal_at_pts(i) / t_normal_at_pts.l2();
471 ++t_set_normal;
472 ++t_normal_at_pts;
473 }
474 return m_normals_at_pts;
475}
476
477template <int DIM, typename BoundaryEleOp>
479 : public BoundaryEleOp {
481 boost::shared_ptr<CommonData> common_data_ptr, double scale = 1,
482 bool is_axisymmetric = false);
483 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
484
485private:
486 boost::shared_ptr<CommonData> commonDataPtr;
487 const double scaleTraction;
489};
490
491template <int DIM, typename BoundaryEleOp>
493 : public BoundaryEleOp {
495 boost::shared_ptr<CommonData> common_data_ptr,
496 bool is_axisymmetric = false,
497 boost::shared_ptr<Range> contact_range_ptr = nullptr);
498 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
499
501 GradSurfaceDistanceFunction gradSurfaceDistanceFunction =
503
504private:
505 boost::shared_ptr<CommonData> commonDataPtr;
507 boost::shared_ptr<Range> contactRange;
508};
509
510template <int DIM, typename BoundaryEleOp>
511struct OpEvaluateSDFImpl<DIM, GAUSS, BoundaryEleOp> : public BoundaryEleOp {
512 OpEvaluateSDFImpl(boost::shared_ptr<CommonData> common_data_ptr);
513 MoFEMErrorCode doWork(int side, EntityType type, EntData &data);
514
515private:
516 boost::shared_ptr<CommonData> commonDataPtr;
517
519 GradSurfaceDistanceFunction gradSurfaceDistanceFunction =
521 HessSurfaceDistanceFunction hessSurfaceDistanceFunction =
523};
524
525template <int DIM, typename AssemblyBoundaryEleOp>
527 : public AssemblyBoundaryEleOp {
528 OpConstrainBoundaryRhsImpl(const std::string field_name,
529 boost::shared_ptr<CommonData> common_data_ptr,
530 bool is_axisymmetric = false);
531 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &data);
532
534 GradSurfaceDistanceFunction gradSurfaceDistanceFunction =
536
537private:
538 boost::shared_ptr<CommonData> commonDataPtr;
540};
541
542template <int DIM, typename AssemblyBoundaryEleOp>
544 : public AssemblyBoundaryEleOp {
545 OpConstrainBoundaryLhs_dUImpl(const std::string row_field_name,
546 const std::string col_field_name,
547 boost::shared_ptr<CommonData> common_data_ptr,
548 bool is_axisymmetric = false);
549 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
550 EntitiesFieldData::EntData &col_data);
551
553 GradSurfaceDistanceFunction gradSurfaceDistanceFunction =
555 HessSurfaceDistanceFunction hessSurfaceDistanceFunction =
557
558 boost::shared_ptr<CommonData> commonDataPtr;
560};
561
562template <int DIM, typename AssemblyBoundaryEleOp>
564 : public AssemblyBoundaryEleOp {
566 const std::string row_field_name, const std::string col_field_name,
567 boost::shared_ptr<CommonData> common_data_ptr,
568 bool is_axisymmetric = false);
569 MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
570 EntitiesFieldData::EntData &col_data);
571
573 GradSurfaceDistanceFunction gradSurfaceDistanceFunction =
575
576private:
577 boost::shared_ptr<CommonData> commonDataPtr;
579};
580
581template <typename BoundaryEleOp> struct ContactIntegrators {
582 template <int DIM, IntegrationType I>
585
586 template <int DIM, IntegrationType I>
589
590 template <int DIM, IntegrationType I>
592
593 template <AssemblyType A> struct Assembly {
594
596 typename FormsIntegrators<BoundaryEleOp>::template Assembly<A>::OpBase;
597
598 template <int DIM, IntegrationType I>
601
602 template <int DIM, IntegrationType I>
605
606 template <int DIM, IntegrationType I>
609 };
610};
611
612inline double sign(double x) {
613 constexpr auto eps = std::numeric_limits<float>::epsilon();
614 if (std::abs(x) < eps)
615 return 0;
616 else if (x > eps)
617 return 1;
618 else
619 return -1;
620};
621
622inline double w(const double sdf, const double tn) {
623 return sdf - cn_contact * tn;
624}
625
626/**
627 * @brief constrain function
628 *
629 * return 1 if negative sdf or positive tn
630 *
631 * @param sdf signed distance
632 * @param tn traction
633 * @return double
634 */
635inline double constrain(double sdf, double tn) {
636 const auto s = sign(w(sdf, tn));
637 return (1 - s) / 2;
638}
639
640template <int DIM, typename BoundaryEleOp>
643 boost::shared_ptr<CommonData> common_data_ptr, double scale,
644 bool is_axisymmetric)
646 commonDataPtr(common_data_ptr), scaleTraction(scale),
647 isAxisymmetric(is_axisymmetric) {}
648
649template <int DIM, typename BoundaryEleOp>
650MoFEMErrorCode
652 int side, EntityType type, EntData &data) {
654
655 FTensor::Index<'i', DIM> i;
656 FTensor::Tensor1<double, 3> t_sum_t{0., 0., 0.};
657
658 auto t_w = BoundaryEleOp::getFTensor0IntegrationWeight();
659 auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
660 auto t_coords = BoundaryEleOp::getFTensor1CoordsAtGaussPts();
661
662 const auto nb_gauss_pts = BoundaryEleOp::getGaussPts().size2();
663 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
664 double jacobian = 1.;
665 if (isAxisymmetric) {
666 jacobian = 2. * M_PI * t_coords(0);
667 }
668 const double alpha = t_w * jacobian * BoundaryEleOp::getMeasure();
669 t_sum_t(i) += alpha * t_traction(i);
670 ++t_w;
671 ++t_traction;
672 ++t_coords;
673 }
674
675 t_sum_t(i) *= scaleTraction;
676
677 constexpr int ind[] = {0, 1, 2};
678 CHKERR VecSetValues(commonDataPtr->totalTraction, 3, ind, &t_sum_t(0),
679 ADD_VALUES);
680
682}
683template <int DIM, typename BoundaryEleOp>
686 boost::shared_ptr<CommonData> common_data_ptr, bool is_axisymmetric,
687 boost::shared_ptr<Range> contact_range_ptr)
689 commonDataPtr(common_data_ptr), isAxisymmetric(is_axisymmetric),
690 contactRange(contact_range_ptr) {}
691
692template <int DIM, typename BoundaryEleOp>
693MoFEMErrorCode
695 int side, EntityType type, EntData &data) {
697
698 auto fe_type = BoundaryEleOp::getFEType();
699
700 const auto fe_ent = BoundaryEleOp::getFEEntityHandle();
701
702 if (contactRange->find(fe_ent) != contactRange->end()) {
703 FTensor::Index<'i', DIM> i;
704 FTensor::Index<'j', DIM> j;
705 FTensor::Tensor1<double, 2> t_sum_a{0., 0.};
706
707 auto t_w = BoundaryEleOp::getFTensor0IntegrationWeight();
708 auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
709 auto t_coords = BoundaryEleOp::getFTensor1CoordsAtGaussPts();
710
711 auto t_grad = getFTensor2FromMat<DIM, DIM>(commonDataPtr->contactDispGrad);
712 auto t_normal_at_pts = BoundaryEleOp::getFTensor1NormalsAtGaussPts();
713
714 const auto nb_gauss_pts = BoundaryEleOp::getGaussPts().size2();
715 auto m_spatial_coords = get_spatial_coords(
716 BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
717 getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
718 auto m_normals_at_pts = get_normalize_normals(
719 BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
720
721 auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
722 auto ts_time = BoundaryEleOp::getTStime();
723 auto ts_time_step = BoundaryEleOp::getTStimeStep();
724 int block_id = 0;
725 auto v_sdf =
726 surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
727 m_spatial_coords, m_normals_at_pts, block_id);
728 auto m_grad_sdf = gradSurfaceDistanceFunction(
729 ts_time_step, ts_time, nb_gauss_pts, m_spatial_coords, m_normals_at_pts,
730 block_id);
731 auto t_sdf = getFTensor0FromVec(v_sdf);
732 auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
733 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
734 double jacobian = 1.;
735 if (isAxisymmetric) {
736 jacobian = 2. * M_PI * t_coords(0); // Axisymmetric Jacobian
737 }
738 auto tn = -t_traction(i) * t_grad_sdf(i);
739 auto c = constrain(t_sdf, tn);
740 double alpha = t_w * jacobian;
741
744 FTensor::Tensor1<double, DIM> t_normal_current;
745
746 F(i, j) = t_grad(i, j) + kronecker_delta(i, j);
747 auto det = determinantTensor(F);
748 CHKERR invertTensor(F, det, invF);
749 t_normal_current(i) = det * (invF(j, i) * t_normal_at_pts(j));
750
751 alpha *= sqrt(t_normal_current(i) * t_normal_current(i));
752
753 if (fe_type == MBTRI) {
754 alpha /= 2;
755 }
756 if (c > 1e-12) {
757 t_sum_a(0) += alpha; // real area
758 }
759 t_sum_a(1) += alpha; // Potential area
760 ++t_w;
761 ++t_traction;
762 ++t_coords;
763 ++t_sdf;
764 ++t_grad_sdf;
765
766 ++t_grad;
767 ++t_normal_at_pts;
768 }
769 constexpr int ind[] = {3, 4};
770 CHKERR VecSetValues(commonDataPtr->totalTraction, 2, ind, &t_sum_a(0),
771 ADD_VALUES);
772 }
774}
775
776template <int DIM, typename BoundaryEleOp>
778 boost::shared_ptr<CommonData> common_data_ptr)
780 commonDataPtr(common_data_ptr) {}
781
782template <int DIM, typename BoundaryEleOp>
783MoFEMErrorCode
785 EntData &data) {
787
788 const auto nb_gauss_pts = BoundaryEleOp::getGaussPts().size2();
789 auto &sdf_vec = commonDataPtr->sdfVals;
790 auto &grad_mat = commonDataPtr->gradsSdf;
791 auto &hess_mat = commonDataPtr->hessSdf;
792 auto &constraint_vec = commonDataPtr->constraintVals;
793 auto &contactTraction_mat = commonDataPtr->contactTraction;
794
795 sdf_vec.resize(nb_gauss_pts, false);
796 grad_mat.resize(DIM, nb_gauss_pts, false);
797 hess_mat.resize((DIM * (DIM + 1)) / 2, nb_gauss_pts, false);
798 constraint_vec.resize(nb_gauss_pts, false);
799
800 auto t_traction = getFTensor1FromMat<DIM>(contactTraction_mat);
801
802 auto t_sdf = getFTensor0FromVec(sdf_vec);
803 auto t_grad_sdf = getFTensor1FromMat<DIM>(grad_mat);
804 auto t_hess_sdf = getFTensor2SymmetricFromMat<DIM>(hess_mat);
805 auto t_constraint = getFTensor0FromVec(constraint_vec);
806
807 auto t_disp = getFTensor1FromMat<DIM>(commonDataPtr->contactDisp);
808 auto t_coords = BoundaryEleOp::getFTensor1CoordsAtGaussPts();
809 auto t_normal_at_pts = BoundaryEleOp::getFTensor1NormalsAtGaussPts();
810
811 FTensor::Index<'i', DIM> i;
812 FTensor::Index<'j', DIM> j;
813
814 auto ts_time = BoundaryEleOp::getTStime();
815 auto ts_time_step = BoundaryEleOp::getTStimeStep();
816
817 auto m_spatial_coords = get_spatial_coords(
818 BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
819 getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
820 auto m_normals_at_pts = get_normalize_normals(
821 BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
822
823 // placeholder to pass boundary block id to python
824 int block_id = 0;
825
826 auto v_sdf =
827 surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
828 m_spatial_coords, m_normals_at_pts, block_id);
829
830 auto m_grad_sdf =
831 gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
832 m_spatial_coords, m_normals_at_pts, block_id);
833
834 auto m_hess_sdf =
835 hessSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
836 m_spatial_coords, m_normals_at_pts, block_id);
837
838 auto t_sdf_v = getFTensor0FromVec(v_sdf);
839 auto t_grad_sdf_v = getFTensor1FromMat<3>(m_grad_sdf);
840 auto t_hess_sdf_v = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
841
842 auto next = [&]() {
843 ++t_sdf;
844 ++t_sdf_v;
845 ++t_grad_sdf;
846 ++t_grad_sdf_v;
847 ++t_hess_sdf;
848 ++t_hess_sdf_v;
849 ++t_disp;
850 ++t_traction;
851 ++t_constraint;
852 };
853
854 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
855
856 auto tn = -t_traction(i) * t_grad_sdf_v(i);
857 auto c = constrain(t_sdf_v, tn);
858
859 t_sdf = t_sdf_v;
860 t_grad_sdf(i) = t_grad_sdf_v(i);
861 t_hess_sdf(i, j) = t_hess_sdf_v(i, j);
862 t_constraint = c;
863
864 next();
865 }
866
868}
869
870template <int DIM, typename AssemblyBoundaryEleOp>
873 boost::shared_ptr<CommonData> common_data_ptr,
874 bool is_axisymmetric)
876 AssemblyBoundaryEleOp::OPROW),
877 commonDataPtr(common_data_ptr), isAxisymmetric(is_axisymmetric) {}
878
879template <int DIM, typename AssemblyBoundaryEleOp>
880MoFEMErrorCode
882 EntitiesFieldData::EntData &data) {
884
885 FTensor::Index<'i', DIM> i;
886 FTensor::Index<'j', DIM> j;
887 FTensor::Index<'k', DIM> k;
888 FTensor::Index<'l', DIM> l;
889
890 const size_t nb_gauss_pts = AssemblyBoundaryEleOp::getGaussPts().size2();
891
892 auto &nf = AssemblyBoundaryEleOp::locF;
893
894 auto t_normal_at_pts = AssemblyBoundaryEleOp::getFTensor1NormalsAtGaussPts();
895
896 auto t_w = AssemblyBoundaryEleOp::getFTensor0IntegrationWeight();
897 auto t_disp = getFTensor1FromMat<DIM>(commonDataPtr->contactDisp);
898 auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
899 auto t_coords = AssemblyBoundaryEleOp::getFTensor1CoordsAtGaussPts();
900
901 size_t nb_base_functions = data.getN().size2() / 3;
902 auto t_base = data.getFTensor1N<3>();
903
904 auto m_spatial_coords = get_spatial_coords(
905 BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
906 getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
907 auto m_normals_at_pts = get_normalize_normals(
908 BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
909
910 auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
911
912 auto ts_time = AssemblyBoundaryEleOp::getTStime();
913 auto ts_time_step = AssemblyBoundaryEleOp::getTStimeStep();
914
915 // placeholder to pass boundary block id to python
916 int block_id = 0;
917
918 auto v_sdf =
919 surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
920 m_spatial_coords, m_normals_at_pts, block_id);
921
922 auto m_grad_sdf =
923 gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
924 m_spatial_coords, m_normals_at_pts, block_id);
925
926 auto t_sdf = getFTensor0FromVec(v_sdf);
927 auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
928
929 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
930
931 auto t_nf = getFTensor1FromPtr<DIM>(&nf[0]);
932
933 double jacobian = 1.;
934 if (isAxisymmetric) {
935 jacobian = 2. * M_PI * t_coords(0);
936 }
937 const double alpha = t_w * jacobian * AssemblyBoundaryEleOp::getMeasure();
938
939 auto tn = -t_traction(i) * t_grad_sdf(i);
940 auto c = constrain(t_sdf, tn);
941
943 t_cP(i, j) = (c * t_grad_sdf(i)) * t_grad_sdf(j);
945 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
946
948 t_rhs(i) =
949
950 t_cQ(i, j) * (t_disp(j) - cn_contact * t_traction(j))
951
952 +
953
954 t_cP(i, j) * t_disp(j) +
955 c * (t_sdf * t_grad_sdf(i)); // add gap0 displacements
956
957 size_t bb = 0;
958 for (; bb != AssemblyBoundaryEleOp::nbRows / DIM; ++bb) {
959 const double beta = alpha * (t_base(i) * t_normal(i));
960 t_nf(i) -= beta * t_rhs(i);
961
962 ++t_nf;
963 ++t_base;
964 }
965 for (; bb < nb_base_functions; ++bb)
966 ++t_base;
967
968 ++t_disp;
969 ++t_traction;
970 ++t_coords;
971 ++t_w;
972 ++t_normal;
973 ++t_sdf;
974 ++t_grad_sdf;
975 }
976
978}
979
980template <int DIM, typename AssemblyBoundaryEleOp>
982 OpConstrainBoundaryLhs_dUImpl(const std::string row_field_name,
983 const std::string col_field_name,
984 boost::shared_ptr<CommonData> common_data_ptr,
985 bool is_axisymmetric)
986 : AssemblyBoundaryEleOp(row_field_name, col_field_name,
987 AssemblyBoundaryEleOp::OPROWCOL),
988 commonDataPtr(common_data_ptr), isAxisymmetric(is_axisymmetric) {
989 AssemblyBoundaryEleOp::sYmm = false;
990}
991
992template <int DIM, typename AssemblyBoundaryEleOp>
993MoFEMErrorCode
995 EntitiesFieldData::EntData &row_data,
996 EntitiesFieldData::EntData &col_data) {
998
999 FTensor::Index<'i', DIM> i;
1000 FTensor::Index<'j', DIM> j;
1001 FTensor::Index<'k', DIM> k;
1002
1003 const size_t nb_gauss_pts = AssemblyBoundaryEleOp::getGaussPts().size2();
1004 auto &locMat = AssemblyBoundaryEleOp::locMat;
1005
1006 auto t_normal_at_pts = AssemblyBoundaryEleOp::getFTensor1NormalsAtGaussPts();
1007 auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
1008 auto t_coords = AssemblyBoundaryEleOp::getFTensor1CoordsAtGaussPts();
1009
1010 auto t_w = AssemblyBoundaryEleOp::getFTensor0IntegrationWeight();
1011 auto t_row_base = row_data.getFTensor1N<3>();
1012 size_t nb_face_functions = row_data.getN().size2() / 3;
1013
1014 auto m_spatial_coords = get_spatial_coords(
1015 BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
1016 getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
1017 auto m_normals_at_pts = get_normalize_normals(
1018 BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
1019
1020 auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
1021
1022 auto ts_time = AssemblyBoundaryEleOp::getTStime();
1023 auto ts_time_step = AssemblyBoundaryEleOp::getTStimeStep();
1024
1025 // placeholder to pass boundary block id to python
1026 int block_id = 0;
1027
1028 auto v_sdf =
1029 surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
1030 m_spatial_coords, m_normals_at_pts, block_id);
1031
1032 auto m_grad_sdf =
1033 gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
1034 m_spatial_coords, m_normals_at_pts, block_id);
1035
1036 auto m_hess_sdf =
1037 hessSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
1038 m_spatial_coords, m_normals_at_pts, block_id);
1039
1040 auto t_sdf = getFTensor0FromVec(v_sdf);
1041 auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
1042 auto t_hess_sdf = getFTensor2SymmetricFromMat<3>(m_hess_sdf);
1043
1044 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1045
1046 double jacobian = 1.;
1047 if (isAxisymmetric) {
1048 jacobian = 2. * M_PI * t_coords(0);
1049 }
1050 const double alpha = t_w * jacobian * AssemblyBoundaryEleOp::getMeasure();
1051
1052 auto tn = -t_traction(i) * t_grad_sdf(i);
1053 auto c = constrain(t_sdf, tn);
1054
1056 t_cP(i, j) = (c * t_grad_sdf(i)) * t_grad_sdf(j);
1058 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
1059
1061 t_res_dU(i, j) = kronecker_delta(i, j) + t_cP(i, j);
1062
1063 if (c > 0) {
1064 t_res_dU(i, j) +=
1065 (c * cn_contact) *
1066 (t_hess_sdf(i, j) * (t_grad_sdf(k) * t_traction(k)) +
1067 t_grad_sdf(i) * t_hess_sdf(k, j) * t_traction(k)) +
1068 c * t_sdf * t_hess_sdf(i, j);
1069 }
1070
1071 size_t rr = 0;
1072 for (; rr != AssemblyBoundaryEleOp::nbRows / DIM; ++rr) {
1073
1074 auto t_mat = getFTensor2FromArray<DIM, DIM, DIM>(locMat, DIM * rr);
1075
1076 const double row_base = t_row_base(i) * t_normal(i);
1077
1078 auto t_col_base = col_data.getFTensor0N(gg, 0);
1079 for (size_t cc = 0; cc != AssemblyBoundaryEleOp::nbCols / DIM; ++cc) {
1080 const double beta = alpha * row_base * t_col_base;
1081
1082 t_mat(i, j) -= beta * t_res_dU(i, j);
1083
1084 ++t_col_base;
1085 ++t_mat;
1086 }
1087
1088 ++t_row_base;
1089 }
1090 for (; rr < nb_face_functions; ++rr)
1091 ++t_row_base;
1092
1093 ++t_traction;
1094 ++t_coords;
1095 ++t_w;
1096 ++t_normal;
1097 ++t_sdf;
1098 ++t_grad_sdf;
1099 ++t_hess_sdf;
1100 }
1101
1103}
1104
1105template <int DIM, typename AssemblyBoundaryEleOp>
1108 const std::string row_field_name, const std::string col_field_name,
1109 boost::shared_ptr<CommonData> common_data_ptr, bool is_axisymmetric)
1110 : AssemblyBoundaryEleOp(row_field_name, col_field_name,
1111 AssemblyBoundaryEleOp::OPROWCOL),
1112 commonDataPtr(common_data_ptr), isAxisymmetric(is_axisymmetric) {
1113 AssemblyBoundaryEleOp::sYmm = false;
1114}
1115
1116template <int DIM, typename AssemblyBoundaryEleOp>
1117MoFEMErrorCode
1119 iNtegrate(EntitiesFieldData::EntData &row_data,
1120 EntitiesFieldData::EntData &col_data) {
1122
1123 FTensor::Index<'i', DIM> i;
1124 FTensor::Index<'j', DIM> j;
1125 FTensor::Index<'k', DIM> k;
1126
1127 const size_t nb_gauss_pts = AssemblyBoundaryEleOp::getGaussPts().size2();
1128 auto &locMat = AssemblyBoundaryEleOp::locMat;
1129
1130 auto t_normal_at_pts = AssemblyBoundaryEleOp::getFTensor1NormalsAtGaussPts();
1131 auto t_traction = getFTensor1FromMat<DIM>(commonDataPtr->contactTraction);
1132 auto t_coords = AssemblyBoundaryEleOp::getFTensor1CoordsAtGaussPts();
1133
1134 auto t_w = AssemblyBoundaryEleOp::getFTensor0IntegrationWeight();
1135 auto t_row_base = row_data.getFTensor1N<3>();
1136 size_t nb_face_functions = row_data.getN().size2() / 3;
1137
1138 auto m_spatial_coords = get_spatial_coords(
1139 BoundaryEleOp::getFTensor1CoordsAtGaussPts(),
1140 getFTensor1FromMat<DIM>(commonDataPtr->contactDisp), nb_gauss_pts);
1141 auto m_normals_at_pts = get_normalize_normals(
1142 BoundaryEleOp::getFTensor1NormalsAtGaussPts(), nb_gauss_pts);
1143
1144 auto t_normal = getFTensor1FromMat<3>(m_normals_at_pts);
1145
1146 auto ts_time = AssemblyBoundaryEleOp::getTStime();
1147 auto ts_time_step = AssemblyBoundaryEleOp::getTStimeStep();
1148
1149 // placeholder to pass boundary block id to python
1150 int block_id = 0;
1151
1152 auto v_sdf =
1153 surfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
1154 m_spatial_coords, m_normals_at_pts, block_id);
1155
1156 auto m_grad_sdf =
1157 gradSurfaceDistanceFunction(ts_time_step, ts_time, nb_gauss_pts,
1158 m_spatial_coords, m_normals_at_pts, block_id);
1159
1160 auto t_sdf = getFTensor0FromVec(v_sdf);
1161 auto t_grad_sdf = getFTensor1FromMat<3>(m_grad_sdf);
1162
1163 for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
1164
1165 double jacobian = 1.;
1166 if (isAxisymmetric) {
1167 jacobian = 2. * M_PI * t_coords(0);
1168 }
1169 const double alpha = t_w * jacobian * AssemblyBoundaryEleOp::getMeasure();
1170
1171 auto tn = -t_traction(i) * t_grad_sdf(i);
1172 auto c = constrain(t_sdf, tn);
1173
1175 t_cP(i, j) = (c * t_grad_sdf(i)) * t_grad_sdf(j);
1177 t_cQ(i, j) = kronecker_delta(i, j) - t_cP(i, j);
1178
1180 t_res_dt(i, j) = -cn_contact * t_cQ(i, j);
1181
1182 size_t rr = 0;
1183 for (; rr != AssemblyBoundaryEleOp::nbRows / DIM; ++rr) {
1184
1185 auto t_mat = getFTensor2FromArray<DIM, DIM, DIM>(locMat, DIM * rr);
1186 const double row_base = t_row_base(i) * t_normal(i);
1187
1188 auto t_col_base = col_data.getFTensor1N<3>(gg, 0);
1189 for (size_t cc = 0; cc != AssemblyBoundaryEleOp::nbCols / DIM; ++cc) {
1190 const double col_base = t_col_base(i) * t_normal(i);
1191 const double beta = alpha * row_base * col_base;
1192
1193 t_mat(i, j) -= beta * t_res_dt(i, j);
1194
1195 ++t_col_base;
1196 ++t_mat;
1197 }
1198
1199 ++t_row_base;
1200 }
1201 for (; rr < nb_face_functions; ++rr)
1202 ++t_row_base;
1203
1204 ++t_traction;
1205 ++t_coords;
1206 ++t_w;
1207 ++t_normal;
1208 ++t_sdf;
1209 ++t_grad_sdf;
1210 }
1211
1213}
1214
1215template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
1216MoFEMErrorCode opFactoryDomainRhs(
1217 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1218 std::string sigma, std::string u, bool is_axisymmetric = false) {
1220
1221 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
1222 A>::template LinearForm<I>;
1223 using OpMixDivURhs = typename B::template OpMixDivTimesU<3, DIM, DIM>;
1224 using OpMixDivUCylRhs =
1225 typename B::template OpMixDivTimesU<3, DIM, DIM, CYLINDRICAL>;
1226
1227 using OpMixLambdaGradURhs = typename B::template OpMixTensorTimesGradU<DIM>;
1228 using OpMixUTimesDivLambdaRhs =
1229 typename B::template OpMixVecTimesDivLambda<SPACE_DIM>;
1230 using OpMixUTimesLambdaRhs =
1231 typename B::template OpGradTimesTensor<1, DIM, DIM>;
1232
1233 auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1234 auto mat_grad_ptr = boost::make_shared<MatrixDouble>();
1235 auto div_stress_ptr = boost::make_shared<MatrixDouble>();
1236 auto contact_stress_ptr = boost::make_shared<MatrixDouble>();
1237
1238 auto jacobian = [is_axisymmetric](const double r, const double,
1239 const double) {
1240 if (is_axisymmetric)
1241 return 2. * M_PI * r;
1242 else
1243 return 1.;
1244 };
1245
1246 pip.push_back(new OpCalculateVectorFieldValues<DIM>(
1247 u, common_data_ptr->contactDispPtr()));
1248 pip.push_back(
1249 new OpCalculateHVecTensorField<DIM, DIM>(sigma, contact_stress_ptr));
1250
1251 if (!is_axisymmetric) {
1252 pip.push_back(
1253 new OpCalculateHVecTensorDivergence<DIM, DIM>(sigma, div_stress_ptr));
1254 } else {
1255 pip.push_back(new OpCalculateHVecTensorDivergence<DIM, DIM, CYLINDRICAL>(
1256 sigma, div_stress_ptr));
1257 }
1258
1259 pip.push_back(new OpCalculateVectorFieldGradient<DIM, DIM>(u, mat_grad_ptr));
1260
1261 if (!is_axisymmetric) {
1262 pip.push_back(
1263 new OpMixDivURhs(sigma, common_data_ptr->contactDispPtr(), jacobian));
1264 } else {
1265 pip.push_back(new OpMixDivUCylRhs(sigma, common_data_ptr->contactDispPtr(),
1266 jacobian));
1267 }
1268
1269 pip.push_back(new OpMixLambdaGradURhs(sigma, mat_grad_ptr, jacobian));
1270 pip.push_back(new OpMixUTimesDivLambdaRhs(u, div_stress_ptr, jacobian));
1271 pip.push_back(new OpMixUTimesLambdaRhs(u, contact_stress_ptr, jacobian));
1272
1274}
1275
1276template <typename OpMixLhs> struct OpMixLhsSide : public OpMixLhs {
1277 using OpMixLhs::OpMixLhs;
1278 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
1279 EntityType col_type,
1280 EntitiesFieldData::EntData &row_data,
1281 EntitiesFieldData::EntData &col_data) {
1283 auto side_fe_entity = OpMixLhs::getSidePtrFE()->getFEEntityHandle();
1284 auto side_fe_data = OpMixLhs::getSideEntity(row_side, row_type);
1285 // Only assemble side which correspond to edge entity on boundary
1286 if (side_fe_entity == side_fe_data) {
1287 CHKERR OpMixLhs::doWork(row_side, col_side, row_type, col_type, row_data,
1288 col_data);
1289 }
1291 }
1292};
1293
1294template <int DIM, AssemblyType A, IntegrationType I, typename DomainEle>
1296 MoFEM::Interface &m_field,
1297 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1298 std::string fe_domain_name, std::string sigma, std::string u,
1299 std::string geom, ForcesAndSourcesCore::RuleHookFun rule,
1300 bool is_axisymmetric = false) {
1302
1303 using DomainEleOp = typename DomainEle::UserDataOperator;
1304
1305 auto op_loop_side = new OpLoopSide<DomainEle>(
1306 m_field, fe_domain_name, DIM, Sev::noisy,
1307 boost::make_shared<ForcesAndSourcesCore::UserDataOperator::AdjCache>());
1308 pip.push_back(op_loop_side);
1309
1310 CHKERR AddHOOps<DIM, DIM, DIM>::add(op_loop_side->getOpPtrVector(),
1311 {H1, HDIV}, geom);
1312
1313 using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
1314 A>::template BiLinearForm<I>;
1315
1316 using OpMixDivULhs = typename B::template OpMixDivTimesVec<DIM>;
1317 using OpMixDivUCylLhs =
1318 typename B::template OpMixDivTimesVec<DIM, CYLINDRICAL>;
1319 using OpLambdaGraULhs = typename B::template OpMixTensorTimesGrad<DIM>;
1320
1321 using OpMixDivULhsSide = OpMixLhsSide<OpMixDivULhs>;
1322 using OpMixDivUCylLhsSide = OpMixLhsSide<OpMixDivUCylLhs>;
1323 using OpLambdaGraULhsSide = OpMixLhsSide<OpLambdaGraULhs>;
1324
1325 auto unity = []() { return 1; };
1326 auto jacobian = [is_axisymmetric](const double r, const double,
1327 const double) {
1328 if (is_axisymmetric)
1329 return 2. * M_PI * r;
1330 else
1331 return 1.;
1332 };
1333
1334 if (!is_axisymmetric) {
1335 op_loop_side->getOpPtrVector().push_back(
1336 new OpMixDivULhsSide(sigma, u, unity, jacobian, true));
1337 } else {
1338 op_loop_side->getOpPtrVector().push_back(
1339 new OpMixDivUCylLhsSide(sigma, u, unity, jacobian, true));
1340 }
1341 op_loop_side->getOpPtrVector().push_back(
1342 new OpLambdaGraULhsSide(sigma, u, unity, jacobian, true));
1343
1344 op_loop_side->getSideFEPtr()->getRuleHook = rule;
1346}
1347
1348template <int DIM, AssemblyType A, IntegrationType I, typename BoundaryEleOp>
1350 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1351 std::string sigma, std::string u, bool is_axisymmetric = false) {
1353
1355
1356 auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1357
1358 pip.push_back(new OpCalculateVectorFieldValues<DIM>(
1359 u, common_data_ptr->contactDispPtr()));
1360 pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
1361 sigma, common_data_ptr->contactTractionPtr()));
1362 pip.push_back(
1363 new typename C::template Assembly<A>::template OpConstrainBoundaryLhs_dU<
1364 DIM, GAUSS>(sigma, u, common_data_ptr, is_axisymmetric));
1365 pip.push_back(new typename C::template Assembly<A>::
1366 template OpConstrainBoundaryLhs_dTraction<DIM, GAUSS>(
1367 sigma, sigma, common_data_ptr, is_axisymmetric));
1368
1370}
1371
1372template <int DIM, AssemblyType A, IntegrationType I, typename BoundaryEleOp>
1374 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1375 std::string sigma, std::string u, bool is_axisymmetric = false) {
1377
1379
1380 auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1381
1382 pip.push_back(new OpCalculateVectorFieldValues<DIM>(
1383 u, common_data_ptr->contactDispPtr()));
1384 pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
1385 sigma, common_data_ptr->contactTractionPtr()));
1386 pip.push_back(
1387 new typename C::template Assembly<A>::template OpConstrainBoundaryRhs<
1388 DIM, GAUSS>(sigma, common_data_ptr, is_axisymmetric));
1389
1391}
1392
1393template <int DIM, IntegrationType I, typename BoundaryEleOp>
1395 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1396 std::string sigma, bool is_axisymmetric = false) {
1398
1400
1401 auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1402 pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
1403 sigma, common_data_ptr->contactTractionPtr()));
1404 pip.push_back(new typename C::template OpAssembleTotalContactTraction<DIM, I>(
1405 common_data_ptr, 1. / scale, is_axisymmetric));
1406
1408}
1409
1410template <int DIM, IntegrationType I, typename BoundaryEleOp>
1412 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1413 OpLoopSide<SideEle> *op_loop_side, std::string sigma, std::string u,
1414 bool is_axisymmetric = false,
1415 boost::shared_ptr<Range> contact_range_ptr = nullptr) {
1418
1419 auto common_data_ptr = boost::make_shared<ContactOps::CommonData>();
1420
1421 op_loop_side->getOpPtrVector().push_back(
1422 new OpCalculateVectorFieldGradient<SPACE_DIM, SPACE_DIM>(
1423 "U", common_data_ptr->contactDispGradPtr()));
1424
1425 if (contact_range_ptr) {
1426 pip.push_back(new OpCalculateVectorFieldValues<DIM>(
1427 u, common_data_ptr->contactDispPtr()));
1428 pip.push_back(new OpCalculateHVecTensorTrace<DIM, BoundaryEleOp>(
1429 sigma, common_data_ptr->contactTractionPtr()));
1430 pip.push_back(op_loop_side);
1431 pip.push_back(new typename C::template OpAssembleTotalContactArea<DIM, I>(
1432 common_data_ptr, is_axisymmetric, contact_range_ptr));
1433 }
1435}
1436
1437}; // namespace ContactOps
1438
1439#endif // __CONTACTOPS_HPP__
static const double eps
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
@ NOSPACE
Definition definitions.h:83
#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_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
@ 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.
@ F
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
constexpr int DIM2
Definition level_set.cpp:22
constexpr int DIM1
Definition level_set.cpp:21
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
boost::function< VectorDouble(double delta_t, double t, int nb_gauss_pts, MatrixDouble &spatial_coords, MatrixDouble &normals_at_pts, int block_id)> SurfaceDistanceFunction
[Common data]
MoFEMErrorCode opFactoryDomainRhs(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string sigma, std::string u, bool is_axisymmetric=false)
auto get_normalize_normals(FTensor::Tensor1< T1, DIM1 > &&t_normal_at_pts, size_t nb_gauss_pts)
double cn_contact
Definition contact.cpp:97
MoFEMErrorCode opFactoryCalculateArea(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, OpLoopSide< SideEle > *op_loop_side, std::string sigma, std::string u, bool is_axisymmetric=false, boost::shared_ptr< Range > contact_range_ptr=nullptr)
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)
MoFEMErrorCode opFactoryBoundaryToDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string fe_domain_name, std::string sigma, std::string u, std::string geom, ForcesAndSourcesCore::RuleHookFun rule, bool is_axisymmetric=false)
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)
MoFEMErrorCode opFactoryBoundaryRhs(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string sigma, std::string u, bool is_axisymmetric=false)
double w(const double sdf, const double tn)
boost::function< MatrixDouble(double delta_t, double t, int nb_gauss_pts, MatrixDouble &spatial_coords, MatrixDouble &normals_at_pts, int block_id)> GradSurfaceDistanceFunction
EntitiesFieldData::EntData EntData
boost::function< MatrixDouble(double delta_t, double t, int nb_gauss_pts, MatrixDouble &spatial_coords, MatrixDouble &normals_at_pts, int block_id)> HessSurfaceDistanceFunction
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)
MoFEMErrorCode opFactoryBoundaryLhs(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string sigma, std::string u, bool is_axisymmetric=false)
MoFEMErrorCode opFactoryCalculateTraction(boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string sigma, bool is_axisymmetric=false)
double sign(double x)
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition sdf.py:1
hess_sdf(delta_t, t, x, y, z, tx, ty, tz, block_id)
Definition sdf.py:20
grad_sdf(delta_t, t, x, y, z, tx, ty, tz, block_id)
Definition sdf.py:16
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
constexpr AssemblyType A
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr auto field_name
VectorDouble sdfVals
size is equal to number of gauss points on element
VectorDouble constraintVals
MatrixDouble contactTraction
static auto getFTensor1TotalTraction()
MatrixDouble contactDisp
static SmartPetscObj< Vec > totalTraction
MatrixDouble contactDispGrad
static auto createTotalTraction(MoFEM::Interface &m_field)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
virtual MPI_Comm & get_comm() const =0
virtual int get_comm_rank() const =0
Deprecated interface functions.
PetscBool is_axisymmetric
Definition contact.cpp:91