v0.15.4
Loading...
Searching...
No Matches
EshelbianCohesive.cpp
Go to the documentation of this file.
1/**
2 * @file EshelbianCohesive.cpp
3 * @author your name (you@domain.com)
4 * @brief
5 * @version 0.1
6 * @date 2025-12-21
7 *
8 * @copyright Copyright (c) 2025
9 *
10 */
11
12#define SINGULARITY
13#include <MoFEM.hpp>
14using namespace MoFEM;
15
17#include <TSElasticPostStep.hpp>
18
19#include <boost/math/constants/constants.hpp>
20#include <boost/math/special_functions/lambert_w.hpp>
21
22#include <EshlabianCohesive.hpp>
23
24namespace EshelbianPlasticity {
25
27
28 static inline constexpr double eps = 1e-8;
29
30 static inline constexpr double pi = boost::math::constants::pi<double>();
31 static inline constexpr double root_pi =
32 boost::math::constants::root_pi<double>();
33
34 template <typename T> static inline T getTau(const T &k, double gc) {
35 const T c = static_cast<T>(gc) / root_pi;
36 return c * std::exp(-k) / std::sqrt(k);
37 }
38
39 template <typename T>
40 static inline auto getDiffTau(const T &tau, const T &kappa, double gc) {
41 return -tau * (1.0 + (1.0 / (2.0 * kappa)));
42 }
43
44 template <typename T>
45 static T getAlpha(const T &kappa, double gc, double min_stiffness) {
46 auto tau = getTau(kappa, gc);
47 return kappa / (kappa * min_stiffness + tau);
48 }
49
50 template <typename T> static T getDiffAlpha(const T &kappa, double gc) {
51 T tau = getTau(kappa, gc);
52 T diff_tau = getDiffTau(tau, kappa, gc);
53 return (tau - kappa * diff_tau) / (tau * tau);
54 }
55
56 template <typename T> static auto invTau(const T &tau, double Gf) {
57 const T z = (2.0 * Gf * Gf) / (pi * tau * tau);
58 return T(0.5) * boost::math::lambert_w(z, 0); // k=0 == principal branch
59 }
60
61 template <typename T>
62 static auto calculateY(const T t_eff, const T &kappa, double gc) {
63 T diff_alpha = getDiffAlpha(kappa, gc);
64 return 0.5 * t_eff * t_eff * diff_alpha;
65 }
66
67 template <typename T>
68 static auto calculateDissipation(const T &delta_kappa, const T t_eff,
69 const T &kappa, double gc) {
70 T Y = calculateY(t_eff, kappa + delta_kappa, gc);
71 return Y * delta_kappa;
72 }
73
74 template <typename T>
75 static auto calculateDissipationSurplus(const T &delta_kappa, const T t_eff,
76 const T &kappa, double gc) {
77 T M = calculateY(t_eff, kappa + delta_kappa, gc) -
78 getTau(kappa + delta_kappa, gc);
79 return -M * delta_kappa;
80 }
81
82 template <typename T>
83 static auto calculateDissipationSurplusDiffKappa(const T &delta_kappa,
84 const T &t_eff,
85 const T &kappa, double gc) {
86
87 std::complex<T> cpx_delta = delta_kappa;
88 std::complex<T> cpx_t_eff = t_eff;
89 std::complex<T> cpx_kappa = kappa;
90 cpx_delta += eps * 1i;
91 std::complex<T> cpx_M =
92 calculateDissipationSurplus(cpx_delta, cpx_t_eff, cpx_kappa, gc);
93 return cpx_M.imag() / eps;
94 }
95
96 template <typename T>
98 const FTensor::Tensor1<T, 3> &t_traction, const T &delta_kappa,
99 const T &kappa, FTensor::Tensor1<double, 3> &t_n_normalize, double gc,
100 double beta) {
101 FTENSOR_INDEX(3, i);
102 FTensor::Tensor1<double, 3> t_diff_traction;
103 std::complex<T> cpx_delta = delta_kappa;
104 std::complex<T> cpx_kappa = kappa;
105 FTensor::Tensor1<std::complex<T>, 3> t_cpx_traction;
106 for (auto jj = 0; jj != 3; ++jj) {
107 t_cpx_traction(i) = t_traction(i);
108 t_cpx_traction(jj) += eps * 1i;
109 std::complex<T> cpx_teff =
110 calculateEffectiveTraction(t_cpx_traction, t_n_normalize, beta);
111 std::complex<T> cpx_M =
112 calculateDissipationSurplus(cpx_delta, cpx_teff, cpx_kappa, gc);
113 t_diff_traction(jj) = cpx_M.imag() / eps;
114 }
115 return t_diff_traction;
116 }
117
118 template <typename T>
119 static auto calculateGap(const FTensor::Tensor1<T, 3> &t_traction,
120 FTensor::Tensor1<double, 3> &t_n_normalize,
121 double alpha, double beta) {
122 FTENSOR_INDEX(3, i);
123 FTENSOR_INDEX(3, j);
124 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
126 t_P(i, j) = t_n_normalize(i) * t_n_normalize(j);
128 t_Q(i, j) = t_kd(i, j) - t_P(i, j);
129 FTensor::Tensor1<T, 3> t_normal;
130 t_normal(i) = t_P(i, j) * t_traction(j);
131 FTensor::Tensor1<T, 3> t_tangential;
132 t_tangential(i) = t_Q(i, j) * t_traction(j);
133 FTensor::Tensor1<T, 3> t_gap_np1;
134 t_gap_np1(i) = alpha * (t_normal(i) + (t_tangential(i) / beta));
135 return t_gap_np1;
136 }
137
138 template <typename T>
139 static auto
141 FTensor::Tensor1<double, 3> &t_n_normalize,
142 double alpha, double beta) {
143 FTENSOR_INDEX(3, i);
145 FTensor::Tensor1<std::complex<T>, 3> t_cpx_delta;
146 FTensor::Tensor1<std::complex<T>, 3> t_cpx_traction;
147 for (auto jj = 0; jj != 3; ++jj) {
148 t_cpx_traction(i) = t_traction(i);
149 t_cpx_traction(jj) += eps * 1i;
150 auto t_cpx_gap = calculateGap(t_cpx_traction, t_n_normalize, alpha, beta);
151 for (auto ii = 0; ii != 3; ++ii) {
152 auto v = t_cpx_gap(ii).imag();
153 t_dgap(ii, jj) = v / eps;
154 }
155 }
156 return t_dgap;
157 }
158
159 template <typename T>
160 static auto
162 FTensor::Tensor1<double, 3> &t_n_normalize,
163 double beta) {
164 auto t = calculateGap(t_traction, t_n_normalize, 1.0, beta);
165 FTENSOR_INDEX(3, i);
166 T t_eff = std::sqrt(t(i) * t(i));
167 return t_eff;
168 }
169
170};
171
173 OpGetParameters(boost::shared_ptr<double> gc_ptr, Sev severity = Sev::inform)
175 gcPtr(gc_ptr), logSev(severity) {
176 CHK_THROW_MESSAGE(getOptions(), "Failed to get EshelbianCohesive options");
177 }
178
179 MoFEMErrorCode doWork(int row_side, EntityType row_type,
180 EntitiesFieldData::EntData &row_data) {
182 *gcPtr = defaultGc;
184 }
185
186 static inline double strength = -1;
187 static inline double min_kappa = 1e-12;
188 static inline double min_stiffness = 1e-2;
189 static inline double kappa0 = 1;
190 static inline double beta = 1;
191
192private:
195
196 PetscOptionsBegin(PETSC_COMM_WORLD, "interface_", "", "none");
197
198 CHKERR PetscOptionsScalar("-gc", "Griffith energy release rate", "",
199 defaultGc, &defaultGc, PETSC_NULLPTR);
200 CHKERR PetscOptionsScalar("-min_stiffness", "Minimal interface stiffness",
201 "", min_stiffness, &min_stiffness, PETSC_NULLPTR);
202 CHKERR PetscOptionsScalar("-strength", "Strength of interface", "",
203 strength, &strength, PETSC_NULLPTR);
204 CHKERR PetscOptionsScalar("-min_kappa",
205 "Minimal kappa to avoid singularity", "",
206 min_kappa, &min_kappa, PETSC_NULLPTR);
207 CHKERR PetscOptionsScalar("-kappa0", "Characteristic length kappa0", "",
208 kappa0, &kappa0, PETSC_NULLPTR);
209 CHKERR PetscOptionsScalar("-beta", "Cohesive tangential coupling", "",
210 beta, &beta, PETSC_NULLPTR);
211
212 PetscOptionsEnd();
213
214 MOFEM_LOG("EP", logSev)
215 << "Interface Griffith energy release rate Gc -interface_gc = "
216 << defaultGc;
217 MOFEM_LOG("EP", logSev)
218 << "Interface min stiffness -interface_min_stiffness = "
219 << min_stiffness;
220 MOFEM_LOG("EP", logSev)
221 << "Interface strength -interface_strength = " << strength;
222 MOFEM_LOG("EP", logSev)
223 << "Interface minimal kappa -interface_min_kappa = " << min_kappa;
224 MOFEM_LOG("EP", logSev)
225 << "Interface characteristic length kappa0 -interface_kappa0 = "
226 << kappa0;
227
228 if (strength > 0) {
229 double kappa_strength =
230 GriffithCohesiveLaw::invTau(static_cast<double>(strength), defaultGc);
231 min_kappa = std::max(min_kappa, kappa_strength);
232 MOFEM_LOG("EP", logSev)
233 << "Adjusted interface min kappa to capture strength = " << min_kappa;
234 }
235
237 }
238
239 double defaultGc = 1.0;
240 boost::shared_ptr<double> gcPtr;
242};
243
244
245
246static Tag get_variable_tag(moab::Interface &moab, std::string tag_name) {
247 Tag tag;
248 const int def_size = 1;
249 rval =
250 moab.tag_get_handle(tag_name.c_str(), def_size, MB_TYPE_DOUBLE, tag,
251 MB_TAG_CREAT | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
252 if (rval == MB_ALREADY_ALLOCATED)
253 rval = MB_SUCCESS;
254 CHK_MOAB_THROW(rval, "Failed to get tag " + tag_name);
255 return tag;
256}
257
258static Tag get_tag(moab::Interface &moab, std::string tag_name, int size) {
259 Tag tag;
260 rval = moab.tag_get_handle(tag_name.c_str(), size, MB_TYPE_DOUBLE, tag,
261 MB_TAG_CREAT | MB_TAG_SPARSE, NULL);
262 if (rval == MB_ALREADY_ALLOCATED)
263 rval = MB_SUCCESS;
264 CHK_MOAB_THROW(rval, "Failed to get tag " + tag_name);
265 return tag;
266};
267
268Tag get_kappa_tag(moab::Interface &moab) {
269 return get_tag(moab, "KAPPA", 1);
270}
271
272Tag get_delta_kappa_tag(moab::Interface &moab) {
273 return get_tag(moab, "DELTA_KAPPA", 1);
274}
275
277 : public FormsIntegrators<FaceElementForcesAndSourcesCore::
278 UserDataOperator>::Assembly<A>::OpBrokenBase {
279
280 using OP =
282 Assembly<A>::OpBrokenBase;
284 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_flux_data_ptr,
285 boost::shared_ptr<Range> ents_ptr = nullptr)
286 : OP(broken_flux_data_ptr, ents_ptr) {}
287
288 MoFEMErrorCode doWork(int row_side, EntityType row_type,
289 EntitiesFieldData::EntData &row_data) {
290
292
293 if (OP::entsPtr) {
294 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
296 }
297
298#ifndef NDEBUG
299 if (!brokenBaseSideData) {
300 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
301 }
302#endif // NDEBUG
303
304 auto do_work_rhs = [this](int row_side, EntityType row_type,
305 EntitiesFieldData::EntData &row_data) {
307#ifndef NDEBUG
308 auto base = row_data.getBase();
309 if (base < 0 && base >= LASTBASE) {
310 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
311 "row base not set properly");
312 }
313#endif // NDEBUG
314
315 // get number of dofs on row
316 OP::nbRows = row_data.getIndices().size();
317 if (!OP::nbRows)
319 // get number of integration points
320 OP::nbIntegrationPts = OP::getGaussPts().size2();
321 // get row base functions
323 // resize and clear the right hand side vector
324 OP::locF.resize(OP::nbRows, false);
325 OP::locF.clear();
326 OP::locMat.resize(OP::nbRows, OP::nbRows, false);
327 OP::locMat.clear();
328 if (OP::nbRows) {
329 // integrate local vector
330 CHKERR this->iNtegrate(row_data);
331 // assemble local vector
332 CHKERR this->aSsemble(row_data);
333 }
335 };
336
337 switch (OP::opType) {
338 case OP::OPSPACE:
339 for (auto &bd : *brokenBaseSideData) {
340 fluxMatPtr =
341 boost::shared_ptr<MatrixDouble>(brokenBaseSideData, &bd.getFlux());
342 faceSense = bd.getSense();
343 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData());
344 fluxMatPtr.reset();
345 faceSense = 0;
346 }
347 break;
348 default:
350 (std::string("wrong op type ") +
351 OpBaseDerivativesBase::OpTypeNames[OP::opType])
352 .c_str());
353 }
354
356 }
357
358protected:
359 boost::shared_ptr<MatrixDouble> fluxMatPtr;
360 int faceSense = 0;
361};
362
364
367 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
368 boost::shared_ptr<double> gc_ptr,
369 boost::shared_ptr<VectorDouble> kappa_ptr,
370 boost::shared_ptr<VectorDouble> kappa_delta_ptr,
371 boost::shared_ptr<std::array<MatrixDouble, 2>> lambda_ptr = nullptr,
372 Tag dissipation_tags = Tag(), Tag grad_dissipation_tags = Tag(),
374 boost::shared_ptr<Range> ents_ptr = nullptr)
375 : OpBrokenBaseCohesive(broken_base_side_data, ents_ptr), gcPtr(gc_ptr),
376 kappaPtr(kappa_ptr), kappaDeltaPtr(kappa_delta_ptr),
377 lambdaPtr(lambda_ptr), dissipationTag(dissipation_tags),
378 gradDissipationTag(grad_dissipation_tags), vec_dJdu(vec_dJ_dx) {}
379
382
383 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
384
385 auto nb_dofs = data.getIndices().size();
386 FTENSOR_INDEX(3, i);
387 FTENSOR_INDEX(3, J);
388
389 int nb_integration_pts = OP::getGaussPts().size2();
390
391 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
392 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
393 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
394 auto t_face_normal = getFTensor1NormalsAtGaussPts();
395
396 int nb_base_functions = data.getN().size2() / 3;
397 auto t_row_base_fun = data.getFTensor1N<3>();
398 auto t_w = OP::getFTensor0IntegrationWeight();
399
400 auto next = [&]() {
401 ++t_P;
402 ++t_face_normal;
403 ++t_kappa;
404 ++t_delta_kappa;
405 ++t_w;
406 };
407
408 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
410 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
411 FTensor::Tensor1<double, 3> t_normalized_normal;
412 t_normalized_normal(J) = t_normal(J);
413 t_normalized_normal.normalize();
415 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
416
417 auto fracture = [this](auto &t_traction, auto &t_normalized_normal,
418 auto &t_kappa, auto &t_delta_kappa, auto gc) {
419 double kappa = static_cast<double>(t_kappa + t_delta_kappa);
423 t_traction, t_normalized_normal, alpha, OpGetParameters::beta);
424 FTENSOR_INDEX(3, i);
425 FTensor::Tensor1<double, 3> t_gap_double;
426 t_gap_double(i) = -t_gap(i);
427 return t_gap_double;
428 };
429
430 auto assemble = [&](const FTensor::Tensor1<double, 3> &&t_gap) {
431 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
432 int bb = 0;
433 for (; bb != nb_dofs / SPACE_DIM; ++bb) {
434 double row_base = t_w * (t_row_base_fun(J) * t_normal(J));
435 t_nf(i) += row_base * t_gap(i);
436 ++t_nf;
437 ++t_row_base_fun;
438 }
439 for (; bb != nb_base_functions; ++bb)
440 ++t_row_base_fun;
441 };
442
443 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
444 *gcPtr));
445
446 next();
447 }
448
450 };
451
452protected:
453 boost::shared_ptr<MatrixDouble> uGammaPtr;
454 boost::shared_ptr<double> gcPtr;
455 boost::shared_ptr<VectorDouble> kappaPtr;
456 boost::shared_ptr<VectorDouble> kappaDeltaPtr;
457 boost::shared_ptr<std::array<MatrixDouble, 2>> lambdaPtr;
458 boost::shared_ptr<double> totalDissipation;
459 boost::shared_ptr<double> tatalDissipationGrad;
463};
464
466
469
472
473 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
474
475 auto nb_dofs = data.getIndices().size();
476 FTENSOR_INDEX(3, i);
477 FTENSOR_INDEX(3, j);
478 FTENSOR_INDEX(3, J);
479
480 int nb_integration_pts = OP::getGaussPts().size2();
481
482 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
483 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
484 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
485 auto t_face_normal = getFTensor1NormalsAtGaussPts();
486
487 int nb_base_functions = data.getN().size2() / 3;
488 auto t_row_base_fun = data.getFTensor1N<3>();
489 auto t_w = OP::getFTensor0IntegrationWeight();
490
491 auto next = [&]() {
492 ++t_P;
493 ++t_kappa;
494 ++t_delta_kappa;
495 ++t_face_normal;
496 ++t_w;
497 };
498
499 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
501 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
502 FTensor::Tensor1<double, 3> t_normalized_normal;
503 t_normalized_normal(J) = t_normal(J);
504 t_normalized_normal.normalize();
506 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
507
508 auto fracture = [this](auto &t_traction, auto &t_normalized_normal,
509 auto &t_kappa, auto &t_delta_kappa, auto gc) {
510 double kappa = static_cast<double>(t_kappa + t_delta_kappa);
514 t_traction, t_normalized_normal, alpha, OpGetParameters::beta);
515 FTENSOR_INDEX(3, i);
516 FTENSOR_INDEX(3, j);
517 FTensor::Tensor2<double, 3, 3> t_dgap_double;
518 t_dgap_double(i, j) = -t_dgap(i, j);
519 return t_dgap_double;
520 };
521
522 auto assemble = [&](const FTensor::Tensor2<double, 3, 3> &&t_dgap) {
523 int rr = 0;
524 for (; rr != nb_dofs / SPACE_DIM; ++rr) {
525 auto t_mat = getFTensor2FromArray<3, 3, 3>(OP::locMat, 3 * rr);
526 double row_base = t_w * (t_row_base_fun(J) * t_normal(J));
527 auto t_col_base_fun = data.getFTensor1N<3>(gg, 0);
528 for (int cc = 0; cc != nb_dofs / SPACE_DIM; ++cc) {
529 double col_base = t_col_base_fun(J) * t_normalized_normal(J);
530 t_mat(i, j) += (row_base * col_base) * t_dgap(i, j);
531 ++t_mat;
532 ++t_col_base_fun;
533 }
534 ++t_row_base_fun;
535 }
536 for (; rr != nb_base_functions; ++rr)
537 ++t_row_base_fun;
538 };
539
540 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
541 *gcPtr));
542
543 next();
544 }
545
547 };
548
551
552 if (!this->timeScalingFun.empty())
553 this->locMat *= this->timeScalingFun(this->getFEMethod()->ts_t);
554 if (!this->feScalingFun.empty())
555 this->locMat *= this->feScalingFun(this->getFEMethod());
556 // assemble local matrix
557 CHKERR this->matSetValuesHook(this, data, data, this->locMat);
558
560 }
561
562private:
563};
564
566
569
572
573 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
574
575 FTENSOR_INDEX(3, i);
576 FTENSOR_INDEX(3, J);
577 FTENSOR_INDEX(3, K);
578
579 int nb_integration_pts = OP::getGaussPts().size2();
580
581 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
582 auto t_lambda = getFTensor2FromMat<3, 3>(lambdaPtr->at(get_sense_index()));
583 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
584 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
585 auto t_face_normal = getFTensor1NormalsAtGaussPts();
586 auto t_w = OP::getFTensor0IntegrationWeight();
587
588 auto next = [&]() {
589 ++t_P;
590 ++t_face_normal;
591 ++t_kappa;
592 ++t_delta_kappa;
593 ++t_lambda;
594 ++t_w;
595 };
596
597 double face_dissipation = 0.0;
598 double face_grad_dissipation = 0.0;
599 auto face_handle = OP::getFEEntityHandle();
600 CHKERR getMoab().tag_get_data(dissipationTag, &face_handle, 1,
601 &face_dissipation);
602 CHKERR getMoab().tag_get_data(gradDissipationTag, &face_handle, 1,
603 &face_grad_dissipation);
604
605 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
607 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
608 FTensor::Tensor1<double, 3> t_normalized_normal;
609 t_normalized_normal(J) = t_normal(J);
610 t_normalized_normal.normalize();
612 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
613
614 auto dJ_dkappa = [](auto &t_delta_kappa, auto &t_traction,
615 auto &t_normalized_normal, auto &t_kappa, auto gc) {
616 double kappa = static_cast<double>(t_kappa);
617 double delta_kappa = static_cast<double>(t_delta_kappa);
619 t_traction, t_normalized_normal, OpGetParameters::beta);
621 delta_kappa, teff, kappa, gc);
622 double m_grad =
624 delta_kappa, teff, kappa, gc);
625 return boost::make_tuple(m, m_grad);
626 };
627
628 auto dr_kappa = [](auto &t_delta_kappa, auto &t_traction,
629 auto &t_normalized_normal, auto &t_kappa, auto gc) {
630 double kappa = static_cast<double>(t_kappa);
631 double delta_kappa = static_cast<double>(t_delta_kappa);
632 double kappa_plus_delta = kappa + delta_kappa;
633 double diff_alpha =
634 GriffithCohesiveLaw::getDiffAlpha(kappa_plus_delta, gc);
636 t_traction, t_normalized_normal, diff_alpha,
638 FTENSOR_INDEX(3, i);
639 FTensor::Tensor1<double, 3> t_gap_double;
640 t_gap_double(i) = -t_gap(i);
641 return t_gap_double;
642 };
643
644 auto [J, dJ] = dJ_dkappa(t_delta_kappa, t_traction, t_normalized_normal,
645 t_kappa, *gcPtr);
646 face_dissipation += t_w * J * t_normal.l2();
647 face_grad_dissipation += t_w * dJ * t_normal.l2();
648
649 auto t_dr = dr_kappa(t_delta_kappa, t_traction, t_normalized_normal,
650 t_kappa, *gcPtr);
652 t_l(i) = t_lambda(i, K) * t_normal(K);
653 face_grad_dissipation -= t_w * t_l(i) * t_dr(i);
654
655 next();
656 }
657
658 CHKERR getMoab().tag_set_data(dissipationTag, &face_handle, 1,
659 &face_dissipation);
660 CHKERR getMoab().tag_set_data(gradDissipationTag, &face_handle, 1,
661 &face_grad_dissipation);
662
664 };
665
670
671protected:
672};
673
675
678
681
682 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
683
684 FTENSOR_INDEX(3, i);
685 FTENSOR_INDEX(3, J);
686
687 int nb_integration_pts = OP::getGaussPts().size2();
688 auto nb_dofs = data.getIndices().size();
689
690 int nb_base_functions = data.getN().size2() / 3;
691 auto t_row_base_fun = data.getFTensor1N<3>();
692 auto t_w = OP::getFTensor0IntegrationWeight();
693
694 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
695 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
696 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
697 auto t_face_normal = getFTensor1NormalsAtGaussPts();
698
699 auto next = [&]() {
700 ++t_P;
701 ++t_face_normal;
702 ++t_kappa;
703 ++t_delta_kappa;
704 ++t_w;
705 };
706
707 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
709 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
710 FTensor::Tensor1<double, 3> t_normalized_normal;
711 t_normalized_normal(J) = t_normal(J);
712 t_normalized_normal.normalize();
714 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
715
716 auto dJ_dtraction = [](auto &t_traction, auto &t_normalized_normal,
717 auto &t_delta_kappa, auto &t_kappa, auto gc) {
718 double kappa = static_cast<double>(t_kappa);
719 double delta_kappa = static_cast<double>(t_delta_kappa);
720 FTENSOR_INDEX(3, i);
721 FTensor::Tensor1<double, 3> t_traction_double;
722 t_traction_double(i) = t_traction(i);
723 auto t_dM =
725 t_traction_double, delta_kappa, kappa, t_normalized_normal, gc,
727 FTensor::Tensor1<double, 3> t_dJ_double;
728 t_dJ_double(i) = t_dM(i);
729 return t_dJ_double;
730 };
731
732 auto assemble = [&](const FTensor::Tensor1<double, 3> &&t_dJ) {
733
734 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
735 int bb = 0;
736 for (; bb != nb_dofs / SPACE_DIM; ++bb) {
737 double row_base = t_w * (t_row_base_fun(J) * t_normal(J));
738 t_nf(i) += row_base * t_dJ(i);
739 ++t_nf;
740 ++t_row_base_fun;
741 }
742 for (; bb != nb_base_functions; ++bb)
743 ++t_row_base_fun;
744 };
745
746 assemble(dJ_dtraction(t_traction, t_normalized_normal, t_delta_kappa,
747 t_kappa, *gcPtr));
748
749 next();
750 }
751
753 };
754
757 return VecSetValues<AssemblyTypeSelector<A>>(vec_dJdu, data, OP::locF,
758 ADD_VALUES);
760 }
761
762protected:
763};
764
766
768
770
772 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
773 Tag tag, TagGetType tag_get_type,
774 boost::shared_ptr<VectorDouble> tag_data_ptr,
775 boost::shared_ptr<Range> ents_ptr = nullptr)
776 : OpBrokenBaseCohesive(broken_base_side_data, ents_ptr), tagHandle(tag),
777 tagGetType(tag_get_type), tagDataPtr(tag_data_ptr) {}
778
781
782 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
783 auto &v = *tagDataPtr;
784
785 switch (tagGetType) {
786 case TagSET:
787 break;
788 case TagGET:
789 v.resize(1);
790 v.clear();
791 break;
792 }
793
794 auto get_data = [&](auto &v) {
796
797 auto &moab = getMoab();
798 auto fe_ent = getFEEntityHandle();
799
800 int size;
801 double *data;
802 rval = moab.tag_get_by_ptr(tagHandle, &fe_ent, 1, (const void **)&data,
803 &size);
804 if (
805
806 rval == MB_SUCCESS && size > 0 && size != v.size()
807
808 ) {
809 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
810 "Inconsistent size of tag data");
811 } else {
812 if (rval != MB_SUCCESS || tagGetType == TagSET) {
813 int tag_size[1];
814 tag_size[0] = v.size();
815 void const *tag_data[] = {&v[0]};
816 CHKERR moab.tag_set_by_ptr(tagHandle, &fe_ent, 1, tag_data, tag_size);
817 } else {
818 CHKERR moab.tag_get_by_ptr(tagHandle, &fe_ent, 1,
819 (const void **)&data, &size);
820 std::copy(data, data + size, v.begin());
821 }
822 }
823
825 };
826
827 CHKERR get_data(v);
828
830 }
831
836
837private:
839 boost::shared_ptr<VectorDouble> tagDataPtr;
841};
842
844 EshelbianCore &ep,
845 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
846 SmartPetscObj<Vec> lambda_vec = SmartPetscObj<Vec>()) {
847
848 using EleOnSide =
850 using SideEleOp = EleOnSide::UserDataOperator;
851
852 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
853 pip, {HDIV, H1, L2}, ep.materialH1Positions, ep.frontAdjEdges);
854
855 auto domain_side_flux = [&](auto &pip) {
856 // flux
857 auto broken_data_ptr =
858 boost::make_shared<std::vector<BrokenBaseSideData>>();
860 broken_data_ptr));
861 auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
863 ep.piolaStress, flux_mat_ptr));
864 pip.push_back(new OpSetFlux<SideEleOp>(broken_data_ptr, flux_mat_ptr));
865
866 return broken_data_ptr;
867 };
868
869 auto get_lambda = [&](auto &pip) {
870 boost::shared_ptr<std::array<MatrixDouble, 2>> array_lambda_ptr;
871 if (lambda_vec) {
872 array_lambda_ptr = boost::make_shared<std::array<MatrixDouble, 2>>();
873 auto lambda_mat_ptr = boost::make_shared<MatrixDouble>();
875 ep.piolaStress, lambda_mat_ptr, boost::make_shared<double>(1.0),
876 lambda_vec));
878 auto op = new OP(NOSPACE, OP::OPSPACE);
879 op->doWorkRhsHook =
880 [lambda_mat_ptr, array_lambda_ptr](
881 DataOperator *base_op_ptr, int side, EntityType type,
884 auto op_ptr = static_cast<OP *>(base_op_ptr);
885 auto get_sense_index = [op_ptr]() {
886 return (op_ptr->getSkeletonSense() == 1) ? 0 : 1;
887 };
888 array_lambda_ptr->at(get_sense_index()) = *lambda_mat_ptr;
890 };
891 pip.push_back(op);
892 }
893 return array_lambda_ptr;
894 };
895
896 return std::make_pair(domain_side_flux(pip), get_lambda(pip));
897};
898
900 EshelbianCore &ep,
901 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
902 boost::shared_ptr<Range> interface_range_ptr,
903 SmartPetscObj<Vec> lambda_vec = SmartPetscObj<Vec>()) {
904
905 auto &m_field = ep.mField;
906
907 using BoundaryEle =
909 using EleOnSide =
911 using SideEleOp = EleOnSide::UserDataOperator;
912 using BdyEleOp = BoundaryEle::UserDataOperator;
913
914 auto face_side = [&]() {
915 // First: Iterate over skeleton FEs adjacent to Domain FEs
916 // Note: BoundaryEle, i.e. uses skeleton interation rule
917 auto op_loop_skeleton_side =
919 interface_range_ptr, Sev::noisy);
920 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = [](int, int, int) {
921 return -1;
922 };
923 op_loop_skeleton_side->getSideFEPtr()->setRuleHook =
924 set_integration_at_front_face;
925
926 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
927 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
928 ep.frontAdjEdges);
929
930 return op_loop_skeleton_side;
931 };
932
933 auto op_loop_skeleton_side = face_side();
934
935 // Second: Iterate over domain FEs adjacent to skelton, particularly one
936 // domain element.
937 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
938 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
939 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
940 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
941 boost::make_shared<CGGUserPolynomialBase>();
942 auto [broken_data_flux_ptr, lambda_ptr] = pushCohesiveOpsDomainImpl(
943 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
944 auto op_reset_broken_data = new ForcesAndSourcesCore::UserDataOperator(
946 op_reset_broken_data->doWorkRhsHook =
947 [broken_data_flux_ptr](DataOperator *base_op_ptr, int side,
948 EntityType type,
950 broken_data_flux_ptr->resize(0);
951 return 0;
952 };
953 op_loop_skeleton_side->getOpPtrVector().push_back(op_reset_broken_data);
954 op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
955
956 auto kappa_ptr = boost::make_shared<VectorDouble>();
957 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
958 auto gc_ptr = boost::make_shared<double>();
959
960 op_loop_skeleton_side->getOpPtrVector().push_back(
961 new OpGetParameters(gc_ptr));
962 op_loop_skeleton_side->getOpPtrVector().push_back(
963 new OpGetTag(broken_data_flux_ptr, get_kappa_tag(m_field.get_moab()),
964 TagGetType::TagGET, kappa_ptr));
965 op_loop_skeleton_side->getOpPtrVector().push_back(new OpGetTag(
966 broken_data_flux_ptr, get_delta_kappa_tag(m_field.get_moab()),
967 TagGetType::TagGET, kappa_delta_ptr));
968
969 return std::make_tuple(op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr,
970 kappa_ptr, kappa_delta_ptr);
971};
972
974 EshelbianCore &ep,
975 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
976 boost::shared_ptr<Range> interface_range_ptr,
977 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
979
980 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
981 kappa_delta_ptr] =
982 pushCohesiveOpsImpl(ep, set_integration_at_front_face,
983 interface_range_ptr);
984
985 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
986 op_loop_skeleton_side->getOpPtrVector().push_back(
988 u_gamma_ptr));
989 op_loop_skeleton_side->getOpPtrVector().push_back(new OpCohesiveRhs(
990 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
991
992 pip.push_back(op_loop_skeleton_side);
993
995};
996
998 EshelbianCore &ep,
999 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1000 boost::shared_ptr<Range> interface_range_ptr,
1001 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1003
1004 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
1005 kappa_delta_ptr] =
1006 pushCohesiveOpsImpl(ep, set_integration_at_front_face,
1007 interface_range_ptr);
1008 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1009 op_loop_skeleton_side->getOpPtrVector().push_back(
1011 u_gamma_ptr));
1012 op_loop_skeleton_side->getOpPtrVector().push_back(new OpCohesiveLhs_dPdP(
1013 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1014
1015 pip.push_back(op_loop_skeleton_side);
1016
1018};
1019
1021 EshelbianCore &ep,
1022 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1023 SmartPetscObj<Vec> lambda_vec) {
1024 auto &m_field = ep.mField;
1025
1026 using EleOnSide =
1028 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1029 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
1030 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1031 boost::make_shared<CGGUserPolynomialBase>();
1032 auto [broken_data_flux_ptr, lambda_ptr] = pushCohesiveOpsDomainImpl(
1033 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
1034
1035 auto op_reset_broken_data = new ForcesAndSourcesCore::UserDataOperator(
1037 op_reset_broken_data->doWorkRhsHook =
1038 [broken_data_flux_ptr](DataOperator *base_op_ptr, int side,
1039 EntityType type,
1041 broken_data_flux_ptr->resize(0);
1042 return 0;
1043 };
1044 pip.push_back(op_reset_broken_data);
1045 pip.push_back(op_loop_domain_side);
1046
1047 auto tag_dissipation = get_tag(m_field.get_moab(), "COHESIVE_DISSIPATION", 1);
1048 auto tag_grad_dissipation =
1049 get_tag(m_field.get_moab(), "COHESIVE_DISSIPATION_GRAD", 1);
1050
1051 auto kappa_ptr = boost::make_shared<VectorDouble>();
1052 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
1053 auto gc_ptr = boost::make_shared<double>();
1054
1055 pip.push_back(new OpGetParameters(gc_ptr, Sev::noisy));
1056 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1057 get_kappa_tag(m_field.get_moab()),
1058 TagGetType::TagGET, kappa_ptr));
1059 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1060 get_delta_kappa_tag(m_field.get_moab()),
1061 TagGetType::TagGET, kappa_delta_ptr));
1062
1063 pip.push_back(new OpCohesive_dJ_dkappa(
1064 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr, lambda_ptr,
1065 tag_dissipation, tag_grad_dissipation));
1066
1067 return std::make_pair(tag_dissipation, tag_grad_dissipation);
1068}
1069
1071 EshelbianCore &ep,
1072 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1073 auto &m_field = ep.mField;
1074
1075 using EleOnSide =
1077 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1078 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
1079 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1080 boost::make_shared<CGGUserPolynomialBase>();
1081 auto [broken_data_flux_ptr, lambda_ptr] =
1082 pushCohesiveOpsDomainImpl(ep, op_loop_domain_side->getOpPtrVector());
1083
1084 auto op_reset_broken_data = new ForcesAndSourcesCore::UserDataOperator(
1086 op_reset_broken_data->doWorkRhsHook =
1087 [broken_data_flux_ptr](DataOperator *base_op_ptr, int side,
1088 EntityType type,
1090 broken_data_flux_ptr->resize(0);
1091 return 0;
1092 };
1093 pip.push_back(op_reset_broken_data);
1094 pip.push_back(op_loop_domain_side);
1095
1096 auto kappa_ptr = boost::make_shared<VectorDouble>();
1097 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
1098 auto gc_ptr = boost::make_shared<double>();
1099
1100 pip.push_back(new OpGetParameters(gc_ptr, Sev::noisy));
1101 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1102 get_kappa_tag(m_field.get_moab()),
1103 TagGetType::TagGET, kappa_ptr));
1104 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1105 get_delta_kappa_tag(m_field.get_moab()),
1106 TagGetType::TagGET, kappa_delta_ptr));
1107 auto dJ_dx_vec = createDMVector(ep.dmElastic);
1108 pip.push_back(new OpCohesive_dJ_dP(broken_data_flux_ptr, gc_ptr, kappa_ptr,
1109 kappa_delta_ptr, lambda_ptr, Tag(), Tag(),
1110 dJ_dx_vec));
1111
1112 return dJ_dx_vec;
1113}
1114
1116 EshelbianCore &ep,
1117 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1118 SmartPetscObj<Vec> lambda_vec,
1119 CommInterface::EntitiesPetscVector &dissipation_vec,
1120 CommInterface::EntitiesPetscVector &grad_dissipation_vec) {
1122
1123 using BoundaryEle =
1125 using BdyEleOp = BoundaryEle::UserDataOperator;
1126
1127 auto get_face_ele = [&]() {
1128 auto fe_ptr = boost::make_shared<BoundaryEle>(ep.mField);
1129 fe_ptr->getRuleHook = [](int, int, int) { return -1; };
1130 fe_ptr->setRuleHook = set_integration_at_front_face;
1131 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
1132 fe_ptr->getOpPtrVector(), {L2}, ep.materialH1Positions,
1133 ep.frontAdjEdges);
1134
1135 auto interface_face = [&](FEMethod *fe_method_ptr) {
1136 auto ent = fe_method_ptr->getFEEntityHandle();
1137 if (
1138
1139 ep.interfaceFaces->find(ent) != ep.interfaceFaces->end()
1140
1141 ) {
1142 return true;
1143 };
1144 return false;
1145 };
1146
1147 fe_ptr->exeTestHook = interface_face;
1148
1149 return fe_ptr;
1150 };
1151
1152 auto face_fe = get_face_ele();
1153 auto [tag_dissipation, tag_grad_dissipation] =
1154 pushCohesive_dJ_dkappa_Impl(ep, face_fe->getOpPtrVector(), lambda_vec);
1155
1156 constexpr double zero = 0.0;
1157 CHKERR ep.mField.get_moab().tag_clear_data(tag_dissipation,
1158 *(ep.interfaceFaces), &zero);
1159 CHKERR ep.mField.get_moab().tag_clear_data(tag_grad_dissipation,
1160 *(ep.interfaceFaces), &zero);
1162 ep.dmElastic, ep.skeletonElement, face_fe, 0, ep.mField.get_comm_size());
1164 ep.mField.get_moab(), dissipation_vec, tag_dissipation);
1166 ep.mField.get_moab(), grad_dissipation_vec, tag_grad_dissipation);
1167
1169}
1170
1172 EshelbianCore &ep,
1173 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1174 SmartPetscObj<KSP> ksp, SmartPetscObj<Vec> lambda_vec) {
1176
1177 using BoundaryEle =
1179 using BdyEleOp = BoundaryEle::UserDataOperator;
1180
1181 auto get_face_ele = [&]() {
1182 auto fe_ptr = boost::make_shared<BoundaryEle>(ep.mField);
1183 fe_ptr->getRuleHook = [](int, int, int) { return -1; };
1184 fe_ptr->setRuleHook = set_integration_at_front_face;
1185 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
1186 fe_ptr->getOpPtrVector(), {L2}, ep.materialH1Positions,
1187 ep.frontAdjEdges);
1188
1189 auto interface_face = [&](FEMethod *fe_method_ptr) {
1190 auto ent = fe_method_ptr->getFEEntityHandle();
1191 if (
1192
1193 ep.interfaceFaces->find(ent) != ep.interfaceFaces->end()
1194
1195 ) {
1196 return true;
1197 };
1198 return false;
1199 };
1200
1201 fe_ptr->exeTestHook = interface_face;
1202
1203 return fe_ptr;
1204 };
1205
1206 auto face_fe = get_face_ele();
1207 auto dJ_dx = pushCohesive_dJ_dx_Impl(ep, face_fe->getOpPtrVector());
1209 ep.dmElastic, ep.skeletonElement, face_fe, 0, ep.mField.get_comm_size());
1210 CHKERR VecAssemblyBegin(dJ_dx);
1211 CHKERR VecAssemblyEnd(dJ_dx);
1212 CHKERR VecGhostUpdateBegin(dJ_dx, ADD_VALUES, SCATTER_REVERSE);
1213 CHKERR VecGhostUpdateEnd(dJ_dx, ADD_VALUES, SCATTER_REVERSE);
1214 CHKERR VecGhostUpdateBegin(dJ_dx, INSERT_VALUES, SCATTER_FORWARD);
1215 CHKERR VecGhostUpdateEnd(dJ_dx, INSERT_VALUES, SCATTER_FORWARD);
1216 double dJ_dx_norm2;
1217 CHKERR VecNorm(dJ_dx, NORM_2, &dJ_dx_norm2);
1218 MOFEM_LOG("EP", Sev::inform)
1219 << "evaluateCohesiveLambdaImpl: Norm of dJ/dx vector: " << dJ_dx_norm2;
1220 constexpr double tol = 1e-16;
1221 if (dJ_dx_norm2 < tol) {
1222 CHKERR VecZeroEntries(lambda_vec);
1223 } else {
1224 CHKERR KSPSolveTranspose(ksp, dJ_dx, lambda_vec);
1225 }
1226 CHKERR VecGhostUpdateBegin(lambda_vec, INSERT_VALUES, SCATTER_FORWARD);
1227 CHKERR VecGhostUpdateEnd(lambda_vec, INSERT_VALUES, SCATTER_FORWARD);
1228 double lambda_norm2;
1229 CHKERR VecNorm(lambda_vec, NORM_2, &lambda_norm2);
1230 MOFEM_LOG("EP", Sev::inform)
1231 << "evaluateCohesiveLambdaImpl: Norm of lambda vector: " << lambda_norm2;
1232
1234}
1235
1240
1241 CHKERR TSSetFromOptions(ts);
1242 CHKERR TSSetStepNumber(ts, 0);
1243 CHKERR TSSetTime(ts, 0);
1244 CHKERR TSSetSolution(ts, x);
1245 double dt;
1246 CHKERR TSGetTimeStep(ts, &dt);
1247 MOFEM_LOG("EP", Sev::inform)
1248 << "evaluatePrimalProblemCohesiveImpl: Time step dt: " << dt;
1249 CHKERR TSSolve(ts, PETSC_NULLPTR);
1250 CHKERR TSSetTimeStep(ts, dt);
1251
1252 CHKERR DMoFEMMeshToLocalVector(ep.dmElastic, x, INSERT_VALUES,
1253 SCATTER_FORWARD);
1254 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
1255 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
1256
1257 double x_norm2;
1258 CHKERR VecNorm(x, NORM_2, &x_norm2);
1259 MOFEM_LOG("EP", Sev::inform)
1260 << "evaluatePrimalProblemCohesiveImpl: Norm of displacement vector: "
1261 << x_norm2;
1262
1264}
1265
1268 EshelbianCore *ep,
1269 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1270 SmartPetscObj<TS> time_solver)
1271 : ep_ptr(ep), setIntegrationAtFrontFace(set_integration_at_front_face),
1272 timeSolver(time_solver),
1273
1274 kspSolVec(createDMVector(ep->dmElastic)),
1275 lambdaVec(createDMVector(ep->dmElastic)),
1276 kappaVec(CommInterface::createEntitiesPetscVector(
1277 ep->mField.get_comm(), ep->mField.get_moab(),
1278 [&ep](Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1279 Sev::noisy)),
1281 ep->mField.get_comm(), ep->mField.get_moab(),
1282 [&ep](Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1283 Sev::noisy)),
1285 ep->mField.get_comm(), ep->mField.get_moab(),
1286 [&ep](Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1287 Sev::noisy)) {}
1288
1290 return vectorDuplicate(kappaVec.second);
1291 }
1292
1296
1300
1304
1305private:
1309
1316 PetscReal *f,
1317 Vec g, void *ctx);
1318};
1319
1320boost::shared_ptr<CohesiveTAOCtx> createCohesiveTAOCtx(
1321 EshelbianCore *ep_ptr,
1322 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1323 SmartPetscObj<TS> time_solver) {
1324 return boost::make_shared<CohesiveTAOCtxImpl>(
1325 ep_ptr, set_integration_at_front_face, time_solver);
1326}
1327
1329 PetscReal *f, Vec g,
1330 void *ctx) {
1331
1333 auto cohesive_ctx = static_cast<CohesiveTAOCtxImpl *>(ctx);
1334 auto &ep = *(cohesive_ctx->ep_ptr);
1335
1336#ifndef NDEBUG
1337 CHKERR VecView(delta_kappa, PETSC_VIEWER_STDOUT_WORLD);
1338#endif
1339 // Set delta_kappa values to the tag
1340 CHKERR VecCopy(delta_kappa, cohesive_ctx->kappaVec.second);
1341 CHKERR VecGhostUpdateBegin(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1342 SCATTER_FORWARD);
1343 CHKERR VecGhostUpdateEnd(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1344 SCATTER_FORWARD);
1345 CHKERR CommInterface::setTagFromVector(
1346 ep.mField.get_moab(), cohesive_ctx->kappaVec,
1348
1349 // solve primal problem
1350 CHKERR evaluatePrimalProblemCohesiveImpl(ep, cohesive_ctx->timeSolver,
1351 cohesive_ctx->kspSolVec);
1352 // solve adjoint problem to get lambda
1353 auto ksp = snesGetKSP(tsGetSNES(cohesive_ctx->timeSolver));
1354 CHKERR evaluateCohesiveLambdaImpl(ep, cohesive_ctx->setIntegrationAtFrontFace,
1355 ksp, cohesive_ctx->lambdaVec);
1356 // evaluate dissipation and its gradient
1358 ep, cohesive_ctx->setIntegrationAtFrontFace, cohesive_ctx->lambdaVec,
1359 cohesive_ctx->dissipationVec, cohesive_ctx->gradDissipationVec);
1360
1361 CHKERR VecSum(cohesive_ctx->dissipationVec.second, f);
1362 CHKERR VecCopy(cohesive_ctx->gradDissipationVec.second, g);
1363
1364 MOFEM_LOG("EP", Sev::inform)
1365 << "Cohesive objective function (negative total dissipation): " << *f;
1366
1368}
1369
1377
1378}; // namespace EshelbianPlasticity
Eshelbian plasticity interface.
Add cohesive elements to Eshelbian plasticity module.
#define FTENSOR_INDEX(DIM, I)
constexpr int SPACE_DIM
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Kronecker Delta class.
Tensor1< T, Tensor_Dim > normalize()
@ LASTBASE
Definition definitions.h:69
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ L2
field with C-1 continuity
Definition definitions.h:88
@ H1
continuous field
Definition definitions.h:85
@ NOSPACE
Definition definitions.h:83
@ HDIV
field with continuous normal traction
Definition definitions.h:87
#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_IMPOSSIBLE_CASE
Definition definitions.h:35
@ 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.
double kappa
constexpr auto t_kd
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:557
#define MOFEM_LOG(channel, severity)
Log.
SeverityLevel
Severity levels.
FTensor::Index< 'i', SPACE_DIM > i
double dt
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'J', DIM1 > J
Definition level_set.cpp:30
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
double tol
Tag get_delta_kappa_tag(moab::Interface &moab)
static Tag get_tag(moab::Interface &moab, std::string tag_name, int size)
static MoFEMErrorCode evaluateDissipationAndGradImpl(EshelbianCore &ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, SmartPetscObj< Vec > lambda_vec, CommInterface::EntitiesPetscVector &dissipation_vec, CommInterface::EntitiesPetscVector &grad_dissipation_vec)
static auto pushCohesiveOpsImpl(EshelbianCore &ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, boost::shared_ptr< Range > interface_range_ptr, SmartPetscObj< Vec > lambda_vec=SmartPetscObj< Vec >())
boost::shared_ptr< CohesiveTAOCtx > createCohesiveTAOCtx(EshelbianCore *ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, SmartPetscObj< TS > time_solver)
static MoFEMErrorCode evaluatePrimalProblemCohesiveImpl(EshelbianCore &ep, SmartPetscObj< TS > ts, SmartPetscObj< Vec > x)
MoFEMErrorCode pushCohesiveOpsLhs(EshelbianCore &ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, boost::shared_ptr< Range > interface_range_ptr, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
static Tag get_variable_tag(moab::Interface &moab, std::string tag_name)
static auto pushCohesiveOpsDomainImpl(EshelbianCore &ep, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, SmartPetscObj< Vec > lambda_vec=SmartPetscObj< Vec >())
Tag get_kappa_tag(moab::Interface &moab)
MoFEMErrorCode pushCohesiveOpsRhs(EshelbianCore &ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, boost::shared_ptr< Range > interface_range_ptr, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
ForcesAndSourcesCore::UserDataOperator UserDataOperator
static auto pushCohesive_dJ_dkappa_Impl(EshelbianCore &ep, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, SmartPetscObj< Vec > lambda_vec)
MoFEMErrorCode initializeCohesiveKappaField(EshelbianCore &ep)
static MoFEMErrorCode evaluateCohesiveLambdaImpl(EshelbianCore &ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, SmartPetscObj< KSP > ksp, SmartPetscObj< Vec > lambda_vec)
PipelineManager::ElementsAndOpsByDim< SPACE_DIM >::FaceSideEle EleOnSide
static auto pushCohesive_dJ_dx_Impl(EshelbianCore &ep, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip)
MoFEMErrorCode cohesiveEvaluateObjectiveAndGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
static MoFEMErrorCodeGeneric< moab::ErrorCode > rval
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
constexpr double t
plate stiffness
Definition plate.cpp:58
constexpr double g
FTensor::Index< 'm', 3 > m
boost::shared_ptr< Range > frontAdjEdges
const std::string skeletonElement
MoFEM::Interface & mField
const std::string materialH1Positions
const std::string elementVolumeName
const std::string piolaStress
boost::shared_ptr< Range > interfaceFaces
const std::string hybridSpatialDisp
SmartPetscObj< DM > dmElastic
Elastic problem.
SmartPetscObj< Vec > duplicateGradientVec() override
SmartPetscObj< Vec > duplicateKappaVec() override
CommInterface::EntitiesPetscVector & getKappaVec() override
CohesiveTAOCtxImpl(EshelbianCore *ep, ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face, SmartPetscObj< TS > time_solver)
ForcesAndSourcesCore::GaussHookFun setIntegrationAtFrontFace
CommInterface::EntitiesPetscVector gradDissipationVec
friend MoFEMErrorCode cohesiveEvaluateObjectiveAndGradient(Tao tao, Vec x, PetscReal *f, Vec g, void *ctx)
SmartPetscObj< Vec > duplicateDissipationVec() override
CommInterface::EntitiesPetscVector kappaVec
CommInterface::EntitiesPetscVector dissipationVec
static auto calculateDissipationSurplus(const T &delta_kappa, const T t_eff, const T &kappa, double gc)
static T getAlpha(const T &kappa, double gc, double min_stiffness)
static auto invTau(const T &tau, double Gf)
static auto getDiffTau(const T &tau, const T &kappa, double gc)
static auto calculateGap(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double alpha, double beta)
static auto calculateDissipationSurplusDiffTraction(const FTensor::Tensor1< T, 3 > &t_traction, const T &delta_kappa, const T &kappa, FTensor::Tensor1< double, 3 > &t_n_normalize, double gc, double beta)
static auto calculateEffectiveTraction(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double beta)
static auto calculateY(const T t_eff, const T &kappa, double gc)
static auto calculateDissipationSurplusDiffKappa(const T &delta_kappa, const T &t_eff, const T &kappa, double gc)
static T getTau(const T &k, double gc)
static T getDiffAlpha(const T &kappa, double gc)
static auto calculateDissipation(const T &delta_kappa, const T t_eff, const T &kappa, double gc)
static auto calculateDiffGapDTraction(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double alpha, double beta)
boost::shared_ptr< MatrixDouble > fluxMatPtr
OpBrokenBaseCohesive(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_flux_data_ptr, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
MoFEMErrorCode aSsemble(EntData &data)
Assemble local vector into global vector.
MoFEMErrorCode iNtegrate(EntData &data)
Class dedicated to integrate operator.
boost::shared_ptr< double > tatalDissipationGrad
MoFEMErrorCode iNtegrate(EntData &data)
Class dedicated to integrate operator.
boost::shared_ptr< double > gcPtr
OpCohesiveRhs(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, boost::shared_ptr< double > gc_ptr, boost::shared_ptr< VectorDouble > kappa_ptr, boost::shared_ptr< VectorDouble > kappa_delta_ptr, boost::shared_ptr< std::array< MatrixDouble, 2 > > lambda_ptr=nullptr, Tag dissipation_tags=Tag(), Tag grad_dissipation_tags=Tag(), SmartPetscObj< Vec > vec_dJ_dx=SmartPetscObj< Vec >(), boost::shared_ptr< Range > ents_ptr=nullptr)
boost::shared_ptr< VectorDouble > kappaDeltaPtr
boost::shared_ptr< VectorDouble > kappaPtr
boost::shared_ptr< MatrixDouble > uGammaPtr
boost::shared_ptr< double > totalDissipation
boost::shared_ptr< std::array< MatrixDouble, 2 > > lambdaPtr
MoFEMErrorCode aSsemble(EntData &data)
Assemble local vector into global vector.
MoFEMErrorCode iNtegrate(EntData &data)
Class dedicated to integrate operator.
MoFEMErrorCode iNtegrate(EntData &data)
Class dedicated to integrate operator.
MoFEMErrorCode aSsemble(EntData &data)
Assemble local vector into global vector.
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
OpGetParameters(boost::shared_ptr< double > gc_ptr, Sev severity=Sev::inform)
MoFEMErrorCode iNtegrate(EntData &data)
Class dedicated to integrate operator.
OpGetTag(boost::shared_ptr< std::vector< BrokenBaseSideData > > broken_base_side_data, Tag tag, TagGetType tag_get_type, boost::shared_ptr< VectorDouble > tag_data_ptr, boost::shared_ptr< Range > ents_ptr=nullptr)
MoFEMErrorCode aSsemble(EntData &data)
Assemble local vector into global vector.
boost::shared_ptr< VectorDouble > tagDataPtr
Managing BitRefLevels.
static MoFEMErrorCode updateEntitiesPetscVector(moab::Interface &moab, EntitiesPetscVector &vec, Tag tag, UpdateGhosts update_gosts=defaultUpdateGhosts)
Exchange data between vector and data.
static EntitiesPetscVector createEntitiesPetscVector(MPI_Comm comm, moab::Interface &moab, int dim, const int nb_coeffs, Sev sev=Sev::verbose, int root_rank=0)
Create a ghost vector for exchanging data.
std::pair< std::pair< Range, Range >, SmartPetscObj< Vec > > EntitiesPetscVector
virtual int get_comm_size() const =0
virtual moab::Interface & get_moab()=0
virtual MPI_Comm & get_comm() const =0
base operator to do operations at Gauss Pt. level
Data on single entity (This is passed as argument to DataOperator::doWork)
FieldApproximationBase & getBase()
Get approximation base.
MatrixDouble & getN(const FieldApproximationBase base)
get base functions this return matrix (nb. of rows is equal to nb. of Gauss pts, nb....
FTensor::Tensor1< FTensor::PackPtr< double *, Tensor_Dim >, Tensor_Dim > getFTensor1N(FieldApproximationBase base)
Get base functions for Hdiv/Hcurl spaces.
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
Structure for user loop methods on finite elements.
@ OPSPACE
operator do Work is execute on space data
boost::function< MoFEMErrorCode(ForcesAndSourcesCore *fe_raw_ptr, int order_row, int order_col, int order_data)> GaussHookFun
VectorDouble locF
local entity vector
TimeFun timeScalingFun
assumes that time variable is set
boost::shared_ptr< Range > entsPtr
Entities on which element is run.
static MatSetValuesHook matSetValuesHook
int nbRows
number of dofs on rows
int nbIntegrationPts
number of integration points
virtual MoFEMErrorCode aSsemble(EntData &row_data, EntData &col_data, const bool trans)
Assemble local matrix into global matrix.
virtual MoFEMErrorCode iNtegrate(EntData &row_data, EntData &col_data)
Integrate grad-grad operator.
MatrixDouble locMat
local entity block matrix
virtual size_t getNbOfBaseFunctions(EntitiesFieldData::EntData &data)
Get number of base functions.
FEFun feScalingFun
set by fe entity handle
int nbRowBaseFunctions
number or row base functions
boost::shared_ptr< std::vector< BrokenBaseSideData > > brokenBaseSideData
Operator for broken loop side.
Calculate tenor field using vectorial base, i.e. Hdiv/Hcurl.
Specialization for double precision vector field values calculation.
Element used to execute operators on side of the element.
Template struct for dimension-specific finite element types.
intrusive_ptr for managing petsc objects
BoundaryEle::UserDataOperator BdyEleOp