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