v0.15.5
Loading...
Searching...
No Matches
SurfaceSlidingConstrains.hpp
Go to the documentation of this file.
1/** \file SurfaceSlidingConstrains.hpp
2 * \brief Implementing surface sliding constrains
3 */
4
5
6
7#ifndef __SURFACE_SLIDING_CONSTRAINS_HPP__
8#define __SURFACE_SLIDING_CONSTRAINS_HPP__
9
10#ifdef WITH_ADOL_C
11
12/**
13 * @brief Implementation of surface sliding constrains
14 *
15 */
16struct GenericSliding {
17
18 GenericSliding() = default;
19 virtual ~GenericSliding() = default;
20
23 boost::shared_ptr<VectorDouble> activeVariablesPtr;
25 boost::shared_ptr<VectorDouble> &active_variables_ptr)
28 activeVariablesPtr(active_variables_ptr) {}
29 MoFEMErrorCode doWork(int side, EntityType type,
30 EntitiesFieldData::EntData &data) {
32 if (type == MBVERTEX) {
33 for (unsigned int dd = 0; dd != data.getFieldData().size(); ++dd)
34 (*activeVariablesPtr)[dd] = data.getFieldData()[dd];
35 }
37 }
38 };
39
40 template <int SizeLambda>
41 struct OpGetActiveDofsPositions
43 boost::shared_ptr<VectorDouble> activeVariablesPtr;
45 const std::string field_name,
46 boost::shared_ptr<VectorDouble> &active_variables_ptr)
49 activeVariablesPtr(active_variables_ptr) {}
50 MoFEMErrorCode doWork(int side, EntityType type,
51 EntitiesFieldData::EntData &data) {
53 if (type == MBVERTEX) {
54 for (unsigned int dd = 0; dd != data.getFieldData().size(); ++dd)
55 (*activeVariablesPtr)[SizeLambda + dd] = data.getFieldData()[dd];
56 }
58 }
59 };
60
61 template <int SizeLambda, int SizePositions>
62 struct OpAssembleRhs : public MoFEM::ForcesAndSourcesCore::UserDataOperator {
63
64 boost::shared_ptr<VectorDouble> resultsPtr;
65
66 OpAssembleRhs(const std::string field_name,
67 boost::shared_ptr<VectorDouble> &results_ptr)
70 resultsPtr(results_ptr) {}
71
72 MoFEMErrorCode doWork(int side, EntityType type,
73 EntitiesFieldData::EntData &data) {
75 if (type != MBVERTEX)
77 VectorInt &indices = data.getIndices();
78 int shift = 0;
79 if (indices.empty())
81 else if (indices.size() == SizeLambda)
82 shift = 0;
83 else if (indices.size() == SizePositions)
84 shift = SizeLambda;
85 else
86 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
87 "Element %s: Data inconsistency nb of indices %ld",
88 getFEName().c_str(), indices.size());
89
90 CHKERR VecSetOption(getSNESf(), VEC_IGNORE_NEGATIVE_INDICES, PETSC_TRUE);
91 CHKERR VecSetValues(getSNESf(), indices.size(), &indices[0],
92 &(*resultsPtr)[shift], ADD_VALUES);
94 }
95 };
96 template <int SizeLambda, int SizePositions>
97 struct OpAssembleLhs : public MoFEM::ForcesAndSourcesCore::UserDataOperator {
98
99 boost::shared_ptr<MatrixDouble> jacobianPtr;
100
101 OpAssembleLhs(const std::string field_name_row,
102 const std::string field_name_col,
103 boost::shared_ptr<MatrixDouble> &jacobian_ptr)
105 field_name_row, field_name_col, UserDataOperator::OPROWCOL),
106 jacobianPtr(jacobian_ptr) {
107 sYmm = false;
108 }
109
110 MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type,
111 EntityType col_type,
112 EntitiesFieldData::EntData &row_data,
113 EntitiesFieldData::EntData &col_data) {
115 if (row_type != MBVERTEX)
117 if (col_type != MBVERTEX)
119 VectorInt &row_indices = row_data.getIndices();
120 VectorInt &col_indices = col_data.getIndices();
121 if (row_indices.empty() || col_indices.empty())
123
124 int shift_row = 0;
125 if (row_indices.size() == SizeLambda)
126 shift_row = 0;
127 else if (row_indices.size() == SizePositions)
128 shift_row = SizeLambda;
129 else
130 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
131 "Data inconsistency");
132
133 int shift_col = 0;
134 if (col_indices.size() == SizeLambda)
135 shift_col = 0;
136 else if (col_indices.size() == SizePositions)
137 shift_col = SizeLambda;
138 else
139 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
140 "Data inconsistency");
141
142 MatrixDouble jac(row_indices.size(), col_indices.size());
143 for (unsigned int rr = 0; rr != row_indices.size(); ++rr) {
144 for (unsigned int cc = 0; cc != col_indices.size(); ++cc) {
145 jac(rr, cc) = (*jacobianPtr)(shift_row + rr, shift_col + cc);
146 if (jac(rr, cc) != jac(rr, cc))
147 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
148 "Jacobian assemble not a number jac(rr,cc) = %3.4f",
149 jac(rr, cc));
150 }
151 }
152 CHKERR MatSetValues(getSNESB(), row_indices.size(), &row_indices[0],
153 col_indices.size(), &col_indices[0], &jac(0, 0),
154 ADD_VALUES);
156 }
157 };
158};
159
160/** \brief Shape preserving constrains, i.e. nodes sliding on body surface.
161
162 Derivation and implementation of constrains preserving body surface,
163 i.e. body shape and volume.
164
165 The idea starts form observation that body shape can be globally characterized
166 by constant calculated as volume over its area
167 \f[
168 \frac{V}{A} = C
169 \f]
170 Above equation expressed in integral form is
171 \f[
172 \int_\Omega \textrm{d}V = C \int_\Gamma \textrm{d}S
173 \f]
174 where notting that,
175 \f[
176 \frac{1}{3}
177 \int_\Omega \textrm{div}[\mathbf{X}] \textrm{d}\Omega
178 =
179 C \int_\Gamma \textrm{d}S
180 \f]
181 and applying Gauss theorem we get
182 \f[
183 \int_\Gamma
184 \mathbf{X}\cdot \frac{\mathbf{N}}{\|\mathbf{N}\|}
185 \textrm{d}\Gamma
186 =
187 3C \int_\Gamma \textrm{d}S.
188 \f]
189 Drooping integrals on both sides, and linearizing equation, we get
190 \f[
191 \frac{\mathbf{N}}{\|\mathbf{N}\|} \cdot \delta \mathbf{X}
192 =
193 3C - \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot \mathbf{X}
194 \f]
195 where \f$\delta \mathbf{X}\f$ is displacement sub-inctrement. Above equation
196 is a constrain if satisfied in body shape and volume is conserved. Final form
197 of constrain equation is \f[ \mathcal{r} =
198 \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot \mathbf{X}
199 -
200 \frac{\mathbf{N_0}}{\|\mathbf{N_0}\|}\cdot \mathbf{X}_0 =
201 \frac{\mathbf{N}}{\|\mathbf{N}\|}\cdot (\mathbf{X}-\mathbf{X}_0)
202 \f]
203
204 In the spirit of finite element method the constrain equation is multiplied
205 by shape functions and enforce using Lagrange multiplier method
206 \f[
207 \int_\Gamma \mathbf{N}^\mathsf{T}_\lambda
208 \left(
209 \frac{\mathbf{N}}{\|\mathbf{N}\|}\mathbf{N}_\mathbf{X}\cdot
210 (\overline{\mathbf{X}}-\overline{\mathbf{X}}_0)
211 \right)
212 \|\mathbf{N}\|
213 \textrm{d}\Gamma
214 =
215 \mathbf{0}.
216 \f]
217 Above equation is nonlinear, applying to it Taylor expansion, we can get form
218 which can be used with Newton interactive method \f[ \begin{split}
219 &\int_\Gamma \mathbf{N}^\mathsf{T}_\lambda
220 \left\{
221 \mathbf{N}\mathbf{N}_\mathbf{X}
222 +
223 \left(\mathbf{X}-\mathbf{X}_0\right) \cdot
224 \left(
225 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta
226 -
227 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi
228 \right)
229 \right\}
230 \textrm{d}\Gamma
231 \cdot
232 \delta
233 \overline{\mathbf{X}}\\
234 =
235 &\int_\Gamma \mathbf{N}^\mathsf{T}_\lambda
236 \mathbf{N}\cdot(\mathbf{X}-\mathbf{X}_0)
237 \textrm{d}\Gamma
238 \end{split}.
239 \f]
240 Equation expressing forces on shape as result of constrains, as result
241 Lagrange multiplier method have following form \f[ \begin{split}
242 &\int_\Gamma
243 \mathbf{N}^\mathsf{T}_\mathbf{X} \cdot \mathbf{N}
244 \mathbf{N}_\lambda
245 \textrm{d}\Gamma
246 \cdot
247 \delta\overline{\lambda}\\
248 +
249 &\int_\Gamma
250 \lambda
251 \mathbf{N}^\mathsf{T}_\mathbf{X}
252 \left(
253 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta
254 -
255 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi
256 \right)
257 \textrm{d}\Gamma
258 \delta\overline{\mathbf{X}}\\
259 =
260 &\int_\Gamma
261 \lambda
262 \mathbf{N}^\mathsf{T}_\mathbf{X} \cdot \mathbf{N}
263 \textrm{d}\Gamma
264 \end{split}
265 \f]
266
267 Above equations are assembled into global system of equations as following
268 \f[
269 \left[
270 \begin{array}{cc}
271 \mathbf{K} + \mathbf{B} & \mathbf{C}^\mathsf{T} \\
272 \mathbf{C} + \mathbf{A} & 0
273 \end{array}
274 \right]
275 \left\{
276 \begin{array}{c}
277 \delta \overline{\mathbf{X}} \\
278 \delta \overline{\lambda}
279 \end{array}
280 \right\}=
281 \left[
282 \begin{array}{c}
283 \mathbf{f} - \mathbf{C}^\mathsf{T}\overline{\lambda} \\
284 \overline{\mathbf{r}}
285 \end{array}
286 \right]
287 \f]
288 where
289 \f[
290 \mathbf{C}=
291 \int_\Gamma
292 \mathbf{N}_\lambda^\mathsf{T}
293 \mathbf{N} \cdot
294 \mathbf{N}_\mathbf{X}
295 \textrm{d}\Gamma,
296 \f]
297 \f[
298 \mathbf{B}=
299 \int_\Gamma
300 \lambda
301 (\mathbf{X}-\mathbf{X}_0)\cdot
302 \left(
303 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta
304 -
305 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi
306 \right)
307 \textrm{d}\Gamma
308 \f]
309 and
310 \f[
311 \mathbf{A}=
312 \int_\Gamma
313 \mathbf{N}^\mathsf{T}_\lambda
314 \left(\mathbf{X}-\mathbf{X}_0\right) \cdot
315 \left(
316 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\xi}\right]\cdot\mathbf{B}_\eta
317 -
318 \textrm{Spin}\left[\frac{\partial\mathbf{X}}{\partial\eta}\right]\cdot\mathbf{B}_\xi
319 \right).
320 \f]
321
322*/
324
325 /** \brief Class implemented by user to detect face orientation
326
327 If mesh generated is with surface mesher, usually you don't have to do
328 nothing, all elements on the surface have consistent orientation. In case of
329 internal faces or if you do something with mesh connectivity which breaks
330 orientation on the face, you have to implement method which will set
331 orientation to face.
332
333 */
334 struct DriverElementOrientation {
335
337
338 virtual MoFEMErrorCode
345 virtual MoFEMErrorCode getElementOrientation(MoFEM::Interface &m_field,
346 const EntityHandle face) {
350 }
351 };
352
354
355 struct MyTriangleFE : public MoFEM::FaceElementForcesAndSourcesCore {
356
357 Mat B;
358 Vec F;
359
361 : MoFEM::FaceElementForcesAndSourcesCore(m_field), B(PETSC_NULLPTR),
362 F(PETSC_NULLPTR) {}
363 int getRule(int order) { return 2 * order; };
364
365 MoFEMErrorCode preProcess() {
367
369
370 if (B != PETSC_NULLPTR) {
371 snes_B = B;
372 }
373
374 if (F != PETSC_NULLPTR) {
375 snes_f = F;
376 }
377
379 }
380 };
381
382 boost::shared_ptr<MyTriangleFE> feRhsPtr, feLhsPtr;
383
384 MyTriangleFE &feRhs;
386 MyTriangleFE &feLhs;
388
389 DriverElementOrientation &crackFrontOrientation;
390
391 double aLpha;
392 MoFEMErrorCode getOptions() {
394 PetscOptionsBegin(PETSC_COMM_WORLD, "",
395 "Get surface sliding constrains element scaling",
396 "none");
397 CHKERR PetscOptionsScalar("-surface_sliding_alpha", "scaling parameter", "",
398 aLpha, &aLpha, PETSC_NULLPTR);
399 PetscOptionsEnd();
401 }
402
404 DriverElementOrientation &orientation)
405 : mField(m_field), feRhsPtr(new MyTriangleFE(m_field)),
406 feLhsPtr(new MyTriangleFE(m_field)), feRhs(*feRhsPtr), feLhs(*feLhsPtr),
407 crackFrontOrientation(orientation), aLpha(1) {
408
409 ierr = getOptions();
410 CHKERRABORT(PETSC_COMM_WORLD, ierr);
411 }
412
413 virtual ~SurfaceSlidingConstrains() = default;
414
415 struct OpJacobian
417
418 const int tAg;
419 boost::shared_ptr<VectorDouble> activeVariablesPtr;
420 boost::shared_ptr<VectorDouble> resultsPtr;
421 boost::shared_ptr<MatrixDouble> jacobianPtr;
423 bool evaluateJacobian;
424
425 const double &aLpha;
426
427 OpJacobian(int tag, const std::string field_name,
428 boost::shared_ptr<VectorDouble> &active_variables_ptr,
429 boost::shared_ptr<VectorDouble> &results_ptr,
430 boost::shared_ptr<MatrixDouble> &jacobian_ptr,
431 DriverElementOrientation &orientation, bool evaluate_jacobian,
432 const double &alpha)
435 tAg(tag), activeVariablesPtr(active_variables_ptr),
436 resultsPtr(results_ptr), jacobianPtr(jacobian_ptr),
437 oRientation(orientation), evaluateJacobian(evaluate_jacobian),
438 aLpha(alpha) {}
439
440 MoFEMErrorCode doWork(int side, EntityType type,
441 EntitiesFieldData::EntData &data) {
443 if (type != MBVERTEX)
445
446
447
448 FTensor::Index<'i', 3> i;
449 FTensor::Index<'j', 3> j;
450 FTensor::Index<'k', 3> k;
453
455 getFEMethod());
457
458 int nb_gauss_pts = data.getN().size1();
459 int nb_base_functions = data.getN().size2();
460 auto t_base1 = data.getFTensor0N();
461 auto t_base2 = data.getFTensor0N();
462 auto t_diff_base = data.getFTensor1DiffN<2>();
464 FTensor::Tensor1<adouble, 3> t_position_ksi;
465 FTensor::Tensor1<adouble, 3> t_position_eta;
469
470 ublas::vector<adouble> c_vec(3);
471 ublas::vector<adouble> f_vec(9);
472 c_vec.clear();
473 f_vec.clear();
474
475 trace_on(tAg);
476
477 ublas::vector<adouble> lambda_dofs(3);
478 for (int dd = 0; dd != 3; ++dd) {
479 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
480 }
481 ublas::vector<adouble> position_dofs(9);
482 for (int dd = 0; dd != 9; ++dd) {
483 position_dofs[dd] <<= (*activeVariablesPtr)[3 + dd];
484 }
485
486 auto t_w = getFTensor0IntegrationWeight();
487 auto t_coord_ref = getFTensor1CoordsAtGaussPts();
488 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
489
491 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
493 &lambda_dofs[0]);
494
495 t_position(i) = 0;
496 t_position_ksi(i) = 0;
497 t_position_eta(i) = 0;
498 lambda = 0;
499
500 for (int bb = 0; bb != nb_base_functions; ++bb) {
501 t_position(i) += t_base1 * t_position_dofs(i);
502 t_position_ksi(i) += t_diff_base(N0) * t_position_dofs(i);
503 t_position_eta(i) += t_diff_base(N1) * t_position_dofs(i);
504 lambda += t_base1 * t_lambda_dof;
505 ++t_base1;
506 ++t_position_dofs;
507 ++t_lambda_dof;
508 ++t_diff_base;
509 }
510
511 t_delta(i) = t_position(i) - t_coord_ref(i);
512 t_normal(k) = FTensor::cross(t_position_ksi(i), t_position_eta(j), k);
513
516 &f_vec[0], &f_vec[1], &f_vec[2]);
517
518 double w = t_w * (0.5 * aLpha * eo);
519
520 for (int bb = 0; bb != nb_base_functions; ++bb) {
521 double alpha = w * t_base2;
522 t_c += alpha * t_normal(i) * t_delta(i);
523 t_f(i) += alpha * (lambda * t_normal(i));
524 ++t_c;
525 ++t_f;
526 ++t_base2;
527 }
528
529 ++t_w;
530 ++t_coord_ref;
531 }
532
533 for (int rr = 0; rr != 3; ++rr)
534 c_vec[rr] >>= (*resultsPtr)[rr];
535
536 for (int rr = 0; rr != 9; ++rr)
537 f_vec(rr) >>= (*resultsPtr)[3 + rr];
538
539 trace_off();
540
541 if (evaluateJacobian) {
542 double *jac_ptr[3 + 9];
543 for (int rr = 0; rr != 3 + 9; ++rr)
544 jac_ptr[rr] = &(*jacobianPtr)(rr, 0);
545
546 // play recorder for jacobians
547 int r =
548 ::jacobian(tAg, 3 + 9, 3 + 9, &(*activeVariablesPtr)[0], jac_ptr);
549 if (r < 0)
550 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
551 "ADOL-C function evaluation with error");
552 }
553
555 }
556 };
557
558 MoFEMErrorCode setOperators(int tag,
559 const std::string lagrange_multipliers_field_name,
560 const std::string material_field_name,
561 const double *alpha = nullptr) {
563
564 if (alpha != nullptr) {
565 aLpha = *alpha;
566 }
567
568 boost::shared_ptr<VectorDouble> active_variables_ptr(
569 new VectorDouble(3 + 9));
570 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(3 + 9));
571 boost::shared_ptr<MatrixDouble> jacobian_ptr(
572 new MatrixDouble(3 + 9, 3 + 9));
573
574 feRhs.getOpPtrVector().clear();
576 lagrange_multipliers_field_name, active_variables_ptr));
578 material_field_name, active_variables_ptr));
579 feRhs.getOpPtrVector().push_back(new OpJacobian(
580 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
581 jacobian_ptr, crackFrontOrientation, false, aLpha));
582 feRhs.getOpPtrVector().push_back(
583 new OpAssembleRhs<3, 9>(lagrange_multipliers_field_name, results_ptr));
584 feRhs.getOpPtrVector().push_back(
585 new OpAssembleRhs<3, 9>(material_field_name, results_ptr));
586
587 // Adding operators to calculate the left hand side
588 feLhs.getOpPtrVector().clear();
590 lagrange_multipliers_field_name, active_variables_ptr));
592 material_field_name, active_variables_ptr));
593 feLhs.getOpPtrVector().push_back(new OpJacobian(
594 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
595 jacobian_ptr, crackFrontOrientation, true, aLpha));
597 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
599 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
601 material_field_name, material_field_name, jacobian_ptr));
602
604 }
605
606 MoFEMErrorCode
608 const std::string lagrange_multipliers_field_name,
609 const std::string material_field_name) {
611
612 boost::shared_ptr<VectorDouble> active_variables_ptr(
613 new VectorDouble(3 + 9));
614 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(3 + 9));
615 boost::shared_ptr<MatrixDouble> jacobian_ptr(
616 new MatrixDouble(3 + 9, 3 + 9));
617
618 // Adding operators to calculate the left hand side
619 feLhs.getOpPtrVector().clear();
621 lagrange_multipliers_field_name, active_variables_ptr));
623 material_field_name, active_variables_ptr));
624 feLhs.getOpPtrVector().push_back(new OpJacobian(
625 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
626 jacobian_ptr, crackFrontOrientation, true, aLpha));
628 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
629
631 }
632};
633
635
636 struct CalculateEdgeBase {
637
638 static MoFEMErrorCode createTag(moab::Interface &moab, Tag &th0, Tag &th1,
639 Tag &th2, Tag &th3) {
641 double def_val[] = {0, 0, 0};
642 CHKERR moab.tag_get_handle("EDGE_BASE0", 3, MB_TYPE_DOUBLE, th0,
643 MB_TAG_CREAT | MB_TAG_SPARSE, def_val);
644 CHKERR moab.tag_get_handle("EDGE_BASE1", 3, MB_TYPE_DOUBLE, th1,
645 MB_TAG_CREAT | MB_TAG_SPARSE, def_val);
646 int def_numb_val[] = {-1};
647 CHKERR moab.tag_get_handle("PATCH_NUMBER", 1, MB_TYPE_INTEGER, th2,
648 MB_TAG_CREAT | MB_TAG_SPARSE, def_numb_val);
649 int def_orientation_val[] = {1};
650 CHKERR moab.tag_get_handle("PATCH_ORIENTATION", 1, MB_TYPE_INTEGER, th3,
651 MB_TAG_CREAT | MB_TAG_SPARSE,
652 def_orientation_val);
653
655 }
656
657 static MoFEMErrorCode numberSurfaces(moab::Interface &moab, Range edges,
658 Range tris) {
660
661 auto get_edges = [&](const Range &ents) {
662 Range edges;
663 CHKERR moab.get_adjacencies(ents, 1, false, edges,
664 moab::Interface::UNION);
665 return edges;
666 };
667
668 auto get_face_adj = [&, get_edges](const Range &faces) {
669 Range adj_faces;
670 CHKERR moab.get_adjacencies(subtract(get_edges(faces), edges), 2, false,
671 adj_faces, moab::Interface::UNION);
672 return intersect(adj_faces, tris);
673 };
674
675 auto get_patch = [&, get_face_adj](const EntityHandle face) {
676 Range patch_ents;
677 patch_ents.insert(face);
678 unsigned int nb0;
679 do {
680 nb0 = patch_ents.size();
681 patch_ents.merge(get_face_adj(patch_ents));
682 } while (nb0 != patch_ents.size());
683 return patch_ents;
684 };
685
686 auto get_patches = [&]() {
687 std::vector<Range> patches;
688 while (!tris.empty()) {
689 patches.push_back(get_patch(tris[0]));
690 tris = subtract(tris, patches.back());
691 }
692 return patches;
693 };
694
695 Tag th0, th1, th2, th3;
696 CHKERR createTag(moab, th0, th1, th2, th3);
697
698 auto patches = get_patches();
699 int pp = 0;
700 for (auto patch : patches) {
701 std::vector<int> tags_vals(patch.size(), pp);
702 CHKERR moab.tag_set_data(th2, patch, &*tags_vals.begin());
703 ++pp;
704 }
705
707 }
708
709 static MoFEMErrorCode setTags(
710 moab::Interface &moab, Range edges, Range tris,
711 bool number_pathes = true,
712 boost::shared_ptr<SurfaceSlidingConstrains::DriverElementOrientation>
713 surface_orientation = nullptr,
714 MoFEM::Interface *m_field_ptr = nullptr) {
716 Tag th0, th1, th2, th3;
717 CHKERR createTag(moab, th0, th1, th2, th3);
718 if (number_pathes) {
719 CHKERR numberSurfaces(moab, edges, tris);
720 }
721 for (auto edge : edges) {
722 Range adj_faces;
723 CHKERR moab.get_adjacencies(&edge, 1, 2, false, adj_faces);
724 adj_faces = intersect(adj_faces, tris);
725 if (adj_faces.size() != 2) {
726 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
727 "Should be 2 faces adjacent to edge but is %ld",
728 adj_faces.size());
729 }
730 VectorDouble3 v[2] = {VectorDouble3(3), VectorDouble3(3)};
731 auto get_tensor_from_vec = [](VectorDouble3 &v) {
733 &v[0], &v[1], &v[2]);
734 };
735
736 FTensor::Index<'i', 3> i;
737 FTensor::Index<'j', 3> j;
738 FTensor::Index<'k', 3> k;
739
740 std::array<double, 3> areas;
741 auto calculate_normals = [&, get_tensor_from_vec]() {
742 int ff = 0;
743 for (auto face : adj_faces) {
744 double &x = (v[ff])[0];
745 double &y = (v[ff])[1];
746 double &z = (v[ff])[2];
747 moab::Util::normal(&moab, face, x, y, z);
748 auto t_n = get_tensor_from_vec(v[ff]);
749 areas[ff] = sqrt(t_n(i) * t_n(i));
750 t_n(i) /= areas[ff];
751 int orientation;
752 // FIXME: this tag is empty
753 CHKERR moab.tag_get_data(th3, &face, 1, &orientation);
754 if (orientation == -1) {
755 t_n(i) *= -1;
756 }
757 if (surface_orientation && m_field_ptr) {
758 CHKERR surface_orientation->getElementOrientation(*m_field_ptr,
759 face);
760 int eo = surface_orientation->elementOrientation;
761 if (eo == -1) {
762 t_n(i) *= -1;
763 }
764 }
765 ++ff;
766 }
767 };
768 calculate_normals();
769
770 auto get_patch_number = [&]() {
771 std::vector<int> p = {0, 0};
772 CHKERR moab.tag_get_data(th2, adj_faces, &*p.begin());
773 return p;
774 };
775
776 std::array<EntityHandle, 2> faces_handles = {adj_faces[0],
777 adj_faces[1]};
778 auto order_normals = [&, get_patch_number]() {
779 auto p = get_patch_number();
780 if (p[0] < p[1]) {
781 v[0].swap(v[1]);
782 faces_handles[0] = adj_faces[1];
783 faces_handles[1] = adj_faces[0];
784 }
785 };
786 order_normals();
787
789 auto t_n0 = get_tensor_from_vec(v[0]);
790 auto t_n1 = get_tensor_from_vec(v[1]);
791 t_cross(k) = FTensor::cross(t_n0(i), t_n1(j), k);
792
793 auto get_base_for_coplanar_faces = [&]() {
795
796 std::vector<EntityHandle> face_conn;
797 CHKERR moab.get_connectivity(&faces_handles[1], 1, face_conn, true);
798 std::vector<EntityHandle> edge_conn;
799 CHKERR moab.get_connectivity(&edge, 1, edge_conn, true);
800 std::sort(face_conn.begin(), face_conn.end());
801 std::sort(edge_conn.begin(), edge_conn.end());
802 int n = 0;
803 for (; n != 2; ++n)
804 if (face_conn[n] != edge_conn[n])
805 break;
806 VectorDouble3 coords_face(3);
807 CHKERR moab.get_coords(&face_conn[n], 1,
808 &*coords_face.data().begin());
809 VectorDouble6 coords_edge(6);
810 CHKERR moab.get_coords(&edge_conn[0], 2,
811 &*coords_edge.data().begin());
812 auto t_edge_n0 = FTensor::Tensor1<double, 3>(
813 coords_edge[0], coords_edge[1], coords_edge[2]);
814 auto t_edge_n1 = FTensor::Tensor1<double, 3>(
815 coords_edge[3], coords_edge[4], coords_edge[5]);
816 auto t_face_n = get_tensor_from_vec(coords_face);
817 t_face_n(i) -= t_edge_n0(i);
819 t_delta(i) = t_edge_n1(i) - t_edge_n0(i);
820 t_delta(i) /= sqrt(t_delta(i) * t_delta(i));
821 t_n1(k) = FTensor::cross(t_n0(i), t_delta(j), k);
822 if (t_n1(i) * t_face_n(i) < 0)
823 t_n1(i) *= -1;
824
826 };
827
828 auto get_base_for_not_planar_faces = [&]() {
830 t_n1(k) = FTensor::cross(t_n0(i), t_cross(j), k);
831 t_n1(i) /= sqrt(t_n1(i) * t_n1(i));
833 };
834
835 // This a case when faces adjacent to edge are coplanar !!!
836 constexpr double tol = 1e-6;
837 if ((t_cross(i) * t_cross(i)) < tol)
838 CHKERR get_base_for_coplanar_faces();
839 else
840 CHKERR get_base_for_not_planar_faces();
841
842 VectorDouble3 &v0 = v[0];
843 VectorDouble3 &v1 = v[1];
844 CHKERR moab.tag_set_data(th0, &edge, 1, &v0[0]);
845 CHKERR moab.tag_set_data(th1, &edge, 1, &v1[0]);
846 }
848 }
849
850 static MoFEMErrorCode saveEdges(moab::Interface &moab, std::string name,
851 Range edges, Range *faces = nullptr) {
853 EntityHandle meshset;
854 Tag ths[4];
855 CHKERR createTag(moab, ths[0], ths[1], ths[2], ths[3]);
856 CHKERR moab.create_meshset(MESHSET_SET, meshset);
857 CHKERR moab.add_entities(meshset, edges);
858 if (faces != nullptr) {
859 CHKERR moab.add_entities(meshset, *faces);
860 }
861 CHKERR moab.write_file(name.c_str(), "VTK", "", &meshset, 1, ths, 4);
862 CHKERR moab.delete_entities(&meshset, 1);
864 }
865 };
866
868
869 struct MyEdgeFE : public MoFEM::EdgeElementForcesAndSourcesCore {
870
871 Mat B;
872 Vec F;
873
875 : MoFEM::EdgeElementForcesAndSourcesCore(m_field), B(PETSC_NULLPTR),
876 F(PETSC_NULLPTR) {}
877 int getRule(int order) { return 2 * order; };
878
879 MoFEMErrorCode preProcess() {
881
883
884 if (B != PETSC_NULLPTR) {
885 snes_B = B;
886 }
887
888 if (F != PETSC_NULLPTR) {
889 snes_f = F;
890 }
891
893 }
894 };
895
896 boost::shared_ptr<MyEdgeFE> feRhsPtr, feLhsPtr;
897
898 MyEdgeFE &feRhs;
900 MyEdgeFE &feLhs;
902
903 double aLpha;
904
905 MoFEMErrorCode getOptions() {
907 PetscOptionsBegin(PETSC_COMM_WORLD, "",
908 "Get edge sliding constrains element scaling",
909 "none");
910 CHKERR PetscOptionsScalar("-edge_sliding_alpha", "scaling parameter", "",
911 aLpha, &aLpha, PETSC_NULLPTR);
912 PetscOptionsEnd();
914 }
915
917 : mField(m_field), feRhsPtr(new MyEdgeFE(m_field)),
918 feLhsPtr(new MyEdgeFE(m_field)), feRhs(*feRhsPtr), feLhs(*feLhsPtr),
919 aLpha(1) {
920 ierr = getOptions();
921 CHKERRABORT(PETSC_COMM_WORLD, ierr);
922 }
923
924 struct OpJacobian
926
927 const int tAg;
928 boost::shared_ptr<VectorDouble> activeVariablesPtr;
929 boost::shared_ptr<VectorDouble> resultsPtr;
930 boost::shared_ptr<MatrixDouble> jacobianPtr;
931 bool evaluateJacobian;
932
933 const double &aLpha;
934
935 OpJacobian(int tag, const std::string field_name,
936 boost::shared_ptr<VectorDouble> &active_variables_ptr,
937 boost::shared_ptr<VectorDouble> &results_ptr,
938 boost::shared_ptr<MatrixDouble> &jacobian_ptr,
939 bool evaluate_jacobian, const double &alpha)
942 tAg(tag), activeVariablesPtr(active_variables_ptr),
943 resultsPtr(results_ptr), jacobianPtr(jacobian_ptr),
944 evaluateJacobian(evaluate_jacobian), aLpha(alpha) {}
945
946 MoFEMErrorCode doWork(int side, EntityType type,
947 EntitiesFieldData::EntData &data) {
949 if (type != MBVERTEX)
951
952 FTensor::Index<'i', 3> i;
953 FTensor::Index<'j', 2> j;
956
957 Tag th0, th1, th2, th3;
959 th1, th2, th3);
960 FTensor::Tensor1<double, 3> t_edge_base0, t_edge_base1;
962 CHKERR getEdgeFE()->mField.get_moab().tag_get_data(th0, &fe_ent, 1,
963 &t_edge_base0(0));
964 CHKERR getEdgeFE()->mField.get_moab().tag_get_data(th1, &fe_ent, 1,
965 &t_edge_base1(0));
966
967 VectorInt &indices = data.getIndices();
968
969 trace_on(tAg);
970
971 ublas::vector<adouble> lambda_dofs(4);
972 for (int dd = 0; dd != 4; ++dd) {
973 lambda_dofs[dd] <<= (*activeVariablesPtr)[dd];
974 }
975 ublas::vector<adouble> position_dofs(6);
976 for (int dd = 0; dd != 6; ++dd) {
977 position_dofs[dd] <<= (*activeVariablesPtr)[4 + dd];
978 }
979
981 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
983 &position_dofs[3], &position_dofs[4], &position_dofs[5]);
984
986 t_tangent(i) = t_node1(i) - t_node0(i);
987 adouble l = sqrt(t_tangent(i) * t_tangent(i));
988 t_tangent(i) /= l;
989
990 adouble t_dot0, t_dot1;
991 t_dot0 = t_edge_base0(i) * t_tangent(i);
992 t_dot1 = t_edge_base1(i) * t_tangent(i);
993
994 FTensor::Tensor1<adouble, 3> t_base0, t_base1;
995 t_base0(i) = t_edge_base0(i) - t_dot0 * t_tangent(i);
996 t_base1(i) = t_edge_base1(i) - t_dot1 * t_tangent(i);
997 t_base0(i) /= sqrt(t_base0(i) * t_base0(i));
998 t_base1(i) /= sqrt(t_base1(i) * t_base1(i));
999
1000 auto t_base_fun1 = data.getFTensor0N();
1001 auto t_base_fun2 = data.getFTensor0N();
1005 auto t_coord_ref = getFTensor1CoordsAtGaussPts();
1006
1007 ublas::vector<adouble> c_vec(4);
1008 ublas::vector<adouble> f_vec(6);
1009 c_vec.clear();
1010 f_vec.clear();
1011
1012 auto t_w = getFTensor0IntegrationWeight();
1013 int nb_gauss_pts = data.getN().size1();
1014 int nb_base_functions = data.getN().size2();
1015 for (int gg = 0; gg != nb_gauss_pts; ++gg) {
1016
1018 &position_dofs[0], &position_dofs[1], &position_dofs[2]);
1020 &lambda_dofs[0], &lambda_dofs[1]);
1021
1022 t_position(i) = 0;
1023 t_lambda(j) = 0;
1024 for (int bb = 0; bb != nb_base_functions; ++bb) {
1025 t_position(i) += t_base_fun1 * t_position_dofs(i);
1026 t_lambda(j) += t_base_fun1 * t_lambda_dof(j);
1027 ++t_base_fun1;
1028 ++t_position_dofs;
1029 ++t_lambda_dof;
1030 }
1031
1032 t_delta(i) = t_position(i) - t_coord_ref(i);
1033 adouble dot0 = t_base0(i) * t_delta(i);
1034 adouble dot1 = t_base1(i) * t_delta(i);
1035
1036 double w = t_w * aLpha;
1037 adouble val, val1, val2;
1039 &c_vec[1]);
1041 &f_vec[0], &f_vec[1], &f_vec[2]);
1042 for (int bb = 0; bb != nb_base_functions; ++bb) {
1043 if (indices[2 * bb] != -1) {
1044 val = w * t_base_fun2 * l;
1045 t_c(N0) += val * dot0;
1046 t_c(N1) += val * dot1;
1047 val1 = val * t_lambda(N0);
1048 val2 = val * t_lambda(N1);
1049 t_f(i) += val1 * t_base0(i) + val2 * t_base1(i);
1050 }
1051 ++t_c;
1052 ++t_f;
1053 ++t_base_fun2;
1054 }
1055
1056 ++t_w;
1057 ++t_coord_ref;
1058 }
1059
1060 for (int rr = 0; rr != 4; ++rr) {
1061 c_vec[rr] >>= (*resultsPtr)[rr];
1062 }
1063 for (int rr = 0; rr != 6; ++rr) {
1064 f_vec(rr) >>= (*resultsPtr)[4 + rr];
1065 }
1066
1067 trace_off();
1068
1069 if (evaluateJacobian) {
1070 double *jac_ptr[4 + 6];
1071 for (int rr = 0; rr != 4 + 6; ++rr) {
1072 jac_ptr[rr] = &(*jacobianPtr)(rr, 0);
1073 }
1074 // play recorder for jacobians
1075 int r =
1076 ::jacobian(tAg, 4 + 6, 4 + 6, &(*activeVariablesPtr)[0], jac_ptr);
1077 if (r < 0) {
1078 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1079 "ADOL-C function evaluation with error");
1080 }
1081 }
1082
1084 }
1085 };
1086
1087 MoFEMErrorCode setOperators(int tag, Range edges, Range faces,
1088 const std::string lagrange_multipliers_field_name,
1089 const std::string material_field_name) {
1092 edges, faces);
1093 CHKERR setOperators(tag, lagrange_multipliers_field_name,
1094 material_field_name);
1096 }
1097
1098 MoFEMErrorCode setOperators(int tag,
1099 const std::string lagrange_multipliers_field_name,
1100 const std::string material_field_name,
1101 const double *alpha = nullptr) {
1103
1104 if (alpha) {
1105 aLpha = *alpha;
1106 }
1107
1108 boost::shared_ptr<VectorDouble> active_variables_ptr(
1109 new VectorDouble(4 + 6));
1110 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(4 + 6));
1111 boost::shared_ptr<MatrixDouble> jacobian_ptr(
1112 new MatrixDouble(4 + 6, 4 + 6));
1113
1114 feRhs.getOpPtrVector().clear();
1116 lagrange_multipliers_field_name, active_variables_ptr));
1118 material_field_name, active_variables_ptr));
1119 feRhs.getOpPtrVector().push_back(new OpJacobian(
1120 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1121 jacobian_ptr, false, aLpha));
1122 feRhs.getOpPtrVector().push_back(
1123 new OpAssembleRhs<4, 6>(lagrange_multipliers_field_name, results_ptr));
1124 feRhs.getOpPtrVector().push_back(
1125 new OpAssembleRhs<4, 6>(material_field_name, results_ptr));
1126
1127 // Adding operators to calculate the left hand side
1128 feLhs.getOpPtrVector().clear();
1130 lagrange_multipliers_field_name, active_variables_ptr));
1132 material_field_name, active_variables_ptr));
1133 feLhs.getOpPtrVector().push_back(new OpJacobian(
1134 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1135 jacobian_ptr, true, aLpha));
1137 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
1139 material_field_name, lagrange_multipliers_field_name, jacobian_ptr));
1141 material_field_name, material_field_name, jacobian_ptr));
1142
1144 }
1145
1146 MoFEMErrorCode
1148 const std::string lagrange_multipliers_field_name,
1149 const std::string material_field_name) {
1151
1153 edges, faces);
1154 CHKERR setOperatorsConstrainOnly(tag, lagrange_multipliers_field_name,
1155 material_field_name);
1157 }
1158
1159 MoFEMErrorCode
1161 const std::string lagrange_multipliers_field_name,
1162 const std::string material_field_name) {
1164
1165 boost::shared_ptr<VectorDouble> active_variables_ptr(
1166 new VectorDouble(4 + 6));
1167 boost::shared_ptr<VectorDouble> results_ptr(new VectorDouble(4 + 6));
1168 boost::shared_ptr<MatrixDouble> jacobian_ptr(
1169 new MatrixDouble(4 + 6, 4 + 6));
1170
1171 // Adding operators to calculate the left hand side
1172 feLhs.getOpPtrVector().clear();
1174 lagrange_multipliers_field_name, active_variables_ptr));
1176 material_field_name, active_variables_ptr));
1177 feLhs.getOpPtrVector().push_back(new OpJacobian(
1178 tag, lagrange_multipliers_field_name, active_variables_ptr, results_ptr,
1179 jacobian_ptr, true, aLpha));
1181 lagrange_multipliers_field_name, material_field_name, jacobian_ptr));
1182
1184 }
1185};
1186
1187#endif // WITH_ADOL_C
1188
1189#endif // __SURFACE_SLIDING_CONSTRAINS_HPP__
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static PetscErrorCode ierr
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ 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.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
constexpr int order
FTensor::Index< 'i', SPACE_DIM > i
static double lambda
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
auto cross(const Tensor1_Expr< A, T, 3, i > &a, const Tensor1_Expr< B, U, 3, j > &b, const Index< k, 3 > &)
Definition cross.hpp:10
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
constexpr auto field_name
static MoFEMErrorCode numberSurfaces(moab::Interface &moab, Range edges, Range tris)
static MoFEMErrorCode createTag(moab::Interface &moab, Tag &th0, Tag &th1, Tag &th2, Tag &th3)
static MoFEMErrorCode setTags(moab::Interface &moab, Range edges, Range tris, bool number_pathes=true, boost::shared_ptr< SurfaceSlidingConstrains::DriverElementOrientation > surface_orientation=nullptr, MoFEM::Interface *m_field_ptr=nullptr)
static MoFEMErrorCode saveEdges(moab::Interface &moab, std::string name, Range edges, Range *faces=nullptr)
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
boost::shared_ptr< VectorDouble > activeVariablesPtr
boost::shared_ptr< VectorDouble > resultsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpJacobian(int tag, const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr, boost::shared_ptr< VectorDouble > &results_ptr, boost::shared_ptr< MatrixDouble > &jacobian_ptr, bool evaluate_jacobian, const double &alpha)
boost::shared_ptr< MatrixDouble > jacobianPtr
MoFEMErrorCode setOperators(int tag, Range edges, Range faces, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
EdgeSlidingConstrains(MoFEM::Interface &m_field)
boost::shared_ptr< MyEdgeFE > feRhsPtr
MoFEMErrorCode setOperators(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name, const double *alpha=nullptr)
boost::shared_ptr< MyEdgeFE > feLhsPtr
MoFEMErrorCode setOperatorsConstrainOnly(int tag, Range edges, Range faces, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
MoFEMErrorCode setOperatorsConstrainOnly(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
boost::shared_ptr< MatrixDouble > jacobianPtr
OpAssembleLhs(const std::string field_name_row, const std::string field_name_col, boost::shared_ptr< MatrixDouble > &jacobian_ptr)
MoFEMErrorCode doWork(int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
OpAssembleRhs(const std::string field_name, boost::shared_ptr< VectorDouble > &results_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< VectorDouble > resultsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
OpGetActiveDofsLambda(const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr)
boost::shared_ptr< VectorDouble > activeVariablesPtr
OpGetActiveDofsPositions(const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr)
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< VectorDouble > activeVariablesPtr
Implementation of surface sliding constrains.
virtual ~GenericSliding()=default
GenericSliding()=default
virtual moab::Interface & get_moab()=0
bool sYmm
If true assume that matrix is symmetric structure.
Deprecated interface functions.
const EdgeElementForcesAndSourcesCore * getEdgeFE()
get pointer to this finite element
FaceElementForcesAndSourcesCore * getFaceFE()
return pointer to Generic Triangle Finite Element object
EntityHandle getFEEntityHandle() const
Return finite element entity handle.
auto getFTensor1CoordsAtGaussPts()
Get coordinates at integration points assuming linear geometry.
auto getFTensor0IntegrationWeight()
Get integration weights.
std::string getFEName() const
Get name of the element.
@ OPCOL
operator doWork function is executed on FE columns
@ OPROW
operator doWork function is executed on FE rows
@ OPROWCOL
operator doWork is executed on FE rows &columns
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
virtual MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
boost::ptr_deque< UserDataOperator > & getOpPtrVector()
Use to push back operator for row operator.
Vec & snes_f
Reference to residual vector.
Mat & snes_B
Reference to preconditioner of Jacobian matrix.
Class implemented by user to detect face orientation.
virtual MoFEMErrorCode getElementOrientation(MoFEM::Interface &m_field, const EntityHandle face)
virtual MoFEMErrorCode getElementOrientation(MoFEM::Interface &m_field, const FEMethod *fe_method_ptr)
MoFEMErrorCode preProcess()
Pre-processing function executed at loop initialization.
OpJacobian(int tag, const std::string field_name, boost::shared_ptr< VectorDouble > &active_variables_ptr, boost::shared_ptr< VectorDouble > &results_ptr, boost::shared_ptr< MatrixDouble > &jacobian_ptr, DriverElementOrientation &orientation, bool evaluate_jacobian, const double &alpha)
boost::shared_ptr< VectorDouble > resultsPtr
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
boost::shared_ptr< VectorDouble > activeVariablesPtr
boost::shared_ptr< MatrixDouble > jacobianPtr
Shape preserving constrains, i.e. nodes sliding on body surface.
boost::shared_ptr< MyTriangleFE > feRhsPtr
virtual ~SurfaceSlidingConstrains()=default
DriverElementOrientation & crackFrontOrientation
boost::shared_ptr< MyTriangleFE > feLhsPtr
SurfaceSlidingConstrains(MoFEM::Interface &m_field, DriverElementOrientation &orientation)
MoFEMErrorCode setOperators(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name, const double *alpha=nullptr)
MoFEMErrorCode setOperatorsConstrainOnly(int tag, const std::string lagrange_multipliers_field_name, const std::string material_field_name)
double tol