v0.15.5
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>
21
22#include <EshelbianCohesive.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>
51 static T getDiffAlpha(const T &kappa, double gc, double min_stiffness) {
52 T tau = getTau(kappa, gc);
53 T diff_tau = getDiffTau(tau, kappa, gc);
54 T dnom = kappa * min_stiffness + tau;
55 return (tau - kappa * diff_tau) / (dnom * dnom);
56 }
57
58 template <typename T> static auto invTau(const T &tau, double Gf) {
59 const T z = (2.0 * Gf * Gf) / (pi * tau * tau);
60 return T(0.5) * boost::math::lambert_w(z, 0); // k=0 == principal branch
61 }
62
63 template <typename T>
64 static auto calculateY(const T t_eff, const T &kappa, double gc,
65 double min_stiffness) {
66 T diff_alpha = getDiffAlpha(kappa, gc, min_stiffness);
67 return 0.5 * t_eff * t_eff * diff_alpha;
68 }
69
70 template <typename T>
71 static auto calculateDissipation(const T &delta_kappa, const T t_eff,
72 const T &kappa, double gc,
73 double min_stiffness) {
74 T Y = calculateY(t_eff, kappa + delta_kappa, gc, min_stiffness);
75 return Y * delta_kappa;
76 }
77
78 template <typename T>
79 static auto calculateDissipationSurplus(const T &delta_kappa, const T t_eff,
80 const T &kappa, double gc,
81 double min_stiffness) {
82 T M = calculateY(t_eff, kappa + delta_kappa, gc, min_stiffness) -
83 getTau(kappa + delta_kappa, gc);
84 return -M * delta_kappa;
85 }
86
87 template <typename T>
88 static auto calculateDissipationSurplusDiffKappa(const T &delta_kappa,
89 const T &t_eff,
90 const T &kappa, double gc,
91 double min_stiffness) {
92
93 std::complex<T> cpx_delta = delta_kappa;
94 std::complex<T> cpx_t_eff = t_eff;
95 std::complex<T> cpx_kappa = kappa;
96 cpx_delta += eps * 1i;
97 std::complex<T> cpx_M =
98 calculateDissipationSurplus(cpx_delta, cpx_t_eff, cpx_kappa, gc, min_stiffness);
99 return cpx_M.imag() / eps;
100 }
101
102 template <typename T>
104 const FTensor::Tensor1<T, 3> &t_traction, const T &delta_kappa,
105 const T &kappa, FTensor::Tensor1<double, 3> &t_n_normalize, double gc,
106 double beta, double min_stiffness) {
107 FTENSOR_INDEX(3, i);
108 FTensor::Tensor1<double, 3> t_diff_traction;
109 std::complex<T> cpx_delta = delta_kappa;
110 std::complex<T> cpx_kappa = kappa;
111 FTensor::Tensor1<std::complex<T>, 3> t_cpx_traction;
112 for (auto jj = 0; jj != 3; ++jj) {
113 t_cpx_traction(i) = t_traction(i);
114 t_cpx_traction(jj) += eps * 1i;
115 std::complex<T> cpx_teff =
116 calculateEffectiveTraction(t_cpx_traction, t_n_normalize, beta);
117 std::complex<T> cpx_M = calculateDissipationSurplus(
118 cpx_delta, cpx_teff, cpx_kappa, gc, min_stiffness);
119 t_diff_traction(jj) = cpx_M.imag() / eps;
120 }
121 return t_diff_traction;
122 }
123
124 template <typename T>
125 static auto calculateGap(const FTensor::Tensor1<T, 3> &t_traction,
126 FTensor::Tensor1<double, 3> &t_n_normalize,
127 double alpha, double beta,
128 bool sign_sensitive = true) {
129 FTENSOR_INDEX(3, i);
130 FTENSOR_INDEX(3, j);
131 constexpr auto t_kd = FTensor::Kronecker_Delta<double>();
133 t_P(i, j) = t_n_normalize(i) * t_n_normalize(j);
135 t_Q(i, j) = t_kd(i, j) - t_P(i, j);
136 FTensor::Tensor1<T, 3> t_normal;
137 t_normal(i) = t_P(i, j) * t_traction(j);
138 FTensor::Tensor1<T, 3> t_tangential;
139 t_tangential(i) = t_Q(i, j) * t_traction(j);
140
141 if (sign_sensitive) {
142 T s = std::sqrt((t_normal(i) * t_normal(i)) / 4. +
143 (1. / beta) * t_tangential(i) * t_tangential(i));
144 T teff = (t_n_normalize(i) * t_normal(i) / 2.) + s;
145
146 FTensor::Tensor1<T, 3> t_gap{0., 0., 0.};
147 if (std::real(s) > std::numeric_limits<double>::epsilon()) {
148 t_gap(i) = alpha * (
149
150 t_n_normalize(i) * teff / 2. +
151
152 (1.0 / 4.0) * (teff / s) * t_normal(i) +
153
154 (1.0 / beta) * (teff / s) * t_tangential(i)
155
156 );
157 }
158 return std::make_pair(teff, t_gap);
159 } else {
160 T teff = std::sqrt((t_normal(i) * t_normal(i)) +
161 (1. / beta) * t_tangential(i) * t_tangential(i));
163 t_gap(i) = alpha * (
164
165 t_normal(i) + (1.0 / beta) * t_tangential(i)
166
167 );
168 return std::make_pair(teff, t_gap);
169 }
170 }
171
172 template <typename T>
173 static auto
175 FTensor::Tensor1<double, 3> &t_n_normalize,
176 double alpha, double beta,
177 bool sign_sensitive = true) {
178 FTENSOR_INDEX(3, i);
180 FTensor::Tensor1<std::complex<T>, 3> t_cpx_delta;
181 FTensor::Tensor1<std::complex<T>, 3> t_cpx_traction;
182 for (auto jj = 0; jj != 3; ++jj) {
183 t_cpx_traction(i) = t_traction(i);
184 t_cpx_traction(jj) += eps * 1i;
185 auto [teff_cpx, t_cpx_gap] = calculateGap(t_cpx_traction, t_n_normalize,
186 alpha, beta, sign_sensitive);
187 for (auto ii = 0; ii != 3; ++ii) {
188 auto v = t_cpx_gap(ii).imag();
189 t_dgap(ii, jj) = v / eps;
190 }
191 }
192 return t_dgap;
193 }
194
195 template <typename T>
196 static auto
198 FTensor::Tensor1<double, 3> &t_n_normalize,
199 double beta) {
200 auto [teff, t_gap] = calculateGap(t_traction, t_n_normalize, 1.0, beta);
201 return teff;
202 }
203};
204
206 OpGetParameters(boost::shared_ptr<double> gc_ptr, Sev severity = Sev::inform)
208 gcPtr(gc_ptr), logSev(severity) {
209 CHK_THROW_MESSAGE(getOptions(), "Failed to get EshelbianCohesive options");
210 }
211
212 MoFEMErrorCode doWork(int row_side, EntityType row_type,
213 EntitiesFieldData::EntData &row_data) {
215 *gcPtr = defaultGc;
217 }
218
219 static inline double strength = -1;
220 static inline double min_kappa = 1e-12;
221 static inline double min_stiffness = 1e-2;
222 static inline double kappa0 = 1;
223 static inline double beta = 1;
224
225private:
228
229 PetscOptionsBegin(PETSC_COMM_WORLD, "interface_", "", "none");
230
231 CHKERR PetscOptionsScalar("-gc", "Griffith energy release rate", "",
232 defaultGc, &defaultGc, PETSC_NULLPTR);
233 CHKERR PetscOptionsScalar("-min_stiffness", "Minimal interface stiffness",
234 "", min_stiffness, &min_stiffness, PETSC_NULLPTR);
235 CHKERR PetscOptionsScalar("-strength", "Strength of interface", "",
236 strength, &strength, PETSC_NULLPTR);
237 CHKERR PetscOptionsScalar("-min_kappa",
238 "Minimal kappa to avoid singularity", "",
239 min_kappa, &min_kappa, PETSC_NULLPTR);
240 CHKERR PetscOptionsScalar("-kappa0", "Characteristic length kappa0", "",
241 kappa0, &kappa0, PETSC_NULLPTR);
242 CHKERR PetscOptionsScalar("-beta", "Cohesive tangential coupling", "",
243 beta, &beta, PETSC_NULLPTR);
244
245 PetscOptionsEnd();
246
247 MOFEM_LOG("EP", logSev)
248 << "Interface Griffith energy release rate Gc -interface_gc = "
249 << defaultGc;
250 MOFEM_LOG("EP", logSev)
251 << "Interface min stiffness -interface_min_stiffness = "
252 << min_stiffness;
253 MOFEM_LOG("EP", logSev)
254 << "Interface strength -interface_strength = " << strength;
255 MOFEM_LOG("EP", logSev)
256 << "Interface minimal kappa -interface_min_kappa = " << min_kappa;
257 MOFEM_LOG("EP", logSev)
258 << "Interface characteristic length kappa0 -interface_kappa0 = "
259 << kappa0;
260 MOFEM_LOG("EP", logSev)
261 << "Interface tangential coupling -interface_beta = " << beta;
262
263 if (strength > 0) {
264 double kappa_strength =
265 GriffithCohesiveLaw::invTau(static_cast<double>(strength), defaultGc);
266 min_kappa = std::max(min_kappa, kappa_strength);
267 MOFEM_LOG("EP", logSev)
268 << "Adjusted interface min kappa to capture strength = " << min_kappa;
269 }
270
272 }
273
274 double defaultGc = 1.0;
275 boost::shared_ptr<double> gcPtr;
277};
278
279
280static Tag get_tag(moab::Interface &moab, std::string tag_name, int size) {
281 std::vector<double> dummy(size, 0.);
282 Tag tag;
283 rval = moab.tag_get_handle(tag_name.c_str(), size, MB_TYPE_DOUBLE, tag,
284 MB_TAG_CREAT | MB_TAG_SPARSE, dummy.data());
285 if (rval == MB_ALREADY_ALLOCATED)
286 rval = MB_SUCCESS;
287 CHK_MOAB_THROW(rval, "Failed to get tag " + tag_name);
288 return tag;
289};
290
291Tag get_kappa_tag(moab::Interface &moab) {
292 return get_tag(moab, "KAPPA", 1);
293}
294
295Tag get_delta_kappa_tag(moab::Interface &moab) {
296 return get_tag(moab, "DELTA_KAPPA", 1);
297}
298
300 : public FormsIntegrators<FaceElementForcesAndSourcesCore::
301 UserDataOperator>::Assembly<A>::OpBrokenBase {
302
303 using OP =
305 Assembly<A>::OpBrokenBase;
307 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_flux_data_ptr,
308 boost::shared_ptr<Range> ents_ptr = nullptr)
309 : OP(broken_flux_data_ptr, ents_ptr) {}
310
311 MoFEMErrorCode doWork(int row_side, EntityType row_type,
312 EntitiesFieldData::EntData &row_data) {
313
315
316 if (OP::entsPtr) {
317 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
319 }
320
321#ifndef NDEBUG
322 if (!brokenBaseSideData) {
323 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE, "space not set");
324 }
325#endif // NDEBUG
326
327 auto do_work_rhs = [this](int row_side, EntityType row_type,
328 EntitiesFieldData::EntData &row_data) {
330#ifndef NDEBUG
331 auto base = row_data.getBase();
332 if (base < 0 && base >= LASTBASE) {
333 SETERRQ(PETSC_COMM_SELF, MOFEM_IMPOSSIBLE_CASE,
334 "row base not set properly");
335 }
336#endif // NDEBUG
337
338 // get number of dofs on row
339 OP::nbRows = row_data.getIndices().size();
340 if (!OP::nbRows)
342 // get number of integration points
343 OP::nbIntegrationPts = OP::getGaussPts().size2();
344 // get row base functions
345 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
346 // resize and clear the right hand side vector
347 OP::locF.resize(OP::nbRows, false);
348 OP::locF.clear();
349 OP::locMat.resize(OP::nbRows, OP::nbRows, false);
350 OP::locMat.clear();
351 if (OP::nbRows) {
352 // integrate local vector
353 CHKERR this->iNtegrate(row_data);
354 // assemble local vector
355 CHKERR this->aSsemble(row_data);
356 }
358 };
359
360 switch (OP::opType) {
361 case OP::OPSPACE:
362 for (auto &bd : *brokenBaseSideData) {
363 fluxMatPtr =
364 boost::shared_ptr<MatrixDouble>(brokenBaseSideData, &bd.getFlux());
365 faceSense = bd.getSense();
366 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData());
367 fluxMatPtr.reset();
368 faceSense = 0;
369 }
370 break;
371 default:
373 (std::string("wrong op type ") +
374 OpBaseDerivativesBase::OpTypeNames[OP::opType])
375 .c_str());
376 }
377
379 }
380
381protected:
382 boost::shared_ptr<MatrixDouble> fluxMatPtr;
383 int faceSense = 0;
384};
385
387
390 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
391 boost::shared_ptr<double> gc_ptr,
392 boost::shared_ptr<VectorDouble> kappa_ptr,
393 boost::shared_ptr<VectorDouble> kappa_delta_ptr,
394 boost::shared_ptr<std::array<MatrixDouble, 2>> lambda_ptr = nullptr,
395 Tag dissipation_tags = Tag(), Tag grad_dissipation_tags = Tag(),
397 boost::shared_ptr<Range> ents_ptr = nullptr)
398 : OpBrokenBaseCohesive(broken_base_side_data, ents_ptr), gcPtr(gc_ptr),
399 kappaPtr(kappa_ptr), kappaDeltaPtr(kappa_delta_ptr),
400 lambdaPtr(lambda_ptr), dissipationTag(dissipation_tags),
401 gradDissipationTag(grad_dissipation_tags), vec_dJdu(vec_dJ_dx) {}
402
405
406 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
407
408 auto nb_dofs = data.getIndices().size();
409 FTENSOR_INDEX(3, i);
410 FTENSOR_INDEX(3, J);
411
412 int nb_integration_pts = OP::getGaussPts().size2();
413
414 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
415 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
416 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
417 auto t_face_normal = getFTensor1NormalsAtGaussPts();
418
419 int nb_base_functions = data.getN().size2() / 3;
420 auto t_row_base_fun = data.getFTensor1N<3>();
421 auto t_w = OP::getFTensor0IntegrationWeight();
422
423 auto next = [&]() {
424 ++t_P;
425 ++t_face_normal;
426 ++t_kappa;
427 ++t_delta_kappa;
428 ++t_w;
429 };
430
431 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
433 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
434 FTensor::Tensor1<double, 3> t_normalized_normal;
435 t_normalized_normal(J) = t_normal(J);
436 t_normalized_normal.normalize();
438 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
439
440 auto fracture = [this](auto &t_traction, auto &t_normalized_normal,
441 auto &t_kappa, auto &t_delta_kappa, auto gc) {
442 double kappa = static_cast<double>(t_kappa + t_delta_kappa);
445 auto [teff, t_gap] = GriffithCohesiveLaw::calculateGap(
446 t_traction, t_normalized_normal, alpha, OpGetParameters::beta,
447 true);
448 FTENSOR_INDEX(3, i);
449 FTensor::Tensor1<double, 3> t_gap_double;
450 t_gap_double(i) = -t_gap(i);
451 return t_gap_double;
452 };
453
454 auto assemble = [&](const FTensor::Tensor1<double, 3> &&t_gap) {
455 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
456 int bb = 0;
457 for (; bb != nb_dofs / SPACE_DIM; ++bb) {
458 double row_base = t_w * (t_row_base_fun(J) * t_normal(J));
459 t_nf(i) += row_base * t_gap(i);
460 ++t_nf;
461 ++t_row_base_fun;
462 }
463 for (; bb != nb_base_functions; ++bb)
464 ++t_row_base_fun;
465 };
466
467 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
468 *gcPtr));
469
470 next();
471 }
472
474 };
475
476protected:
477 boost::shared_ptr<MatrixDouble> uGammaPtr;
478 boost::shared_ptr<double> gcPtr;
479 boost::shared_ptr<VectorDouble> kappaPtr;
480 boost::shared_ptr<VectorDouble> kappaDeltaPtr;
481 boost::shared_ptr<std::array<MatrixDouble, 2>> lambdaPtr;
482 boost::shared_ptr<double> totalDissipation;
483 boost::shared_ptr<double> tatalDissipationGrad;
487};
488
490
493
496
497 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
498
499 auto nb_dofs = data.getIndices().size();
500 FTENSOR_INDEX(3, i);
501 FTENSOR_INDEX(3, j);
502 FTENSOR_INDEX(3, J);
503
504 int nb_integration_pts = OP::getGaussPts().size2();
505
506 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
507 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
508 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
509 auto t_face_normal = getFTensor1NormalsAtGaussPts();
510
511 int nb_base_functions = data.getN().size2() / 3;
512 auto t_row_base_fun = data.getFTensor1N<3>();
513 auto t_w = OP::getFTensor0IntegrationWeight();
514
515 auto next = [&]() {
516 ++t_P;
517 ++t_kappa;
518 ++t_delta_kappa;
519 ++t_face_normal;
520 ++t_w;
521 };
522
523 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
525 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
526 FTensor::Tensor1<double, 3> t_normalized_normal;
527 t_normalized_normal(J) = t_normal(J);
528 t_normalized_normal.normalize();
530 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
531
532 auto fracture = [this](auto &t_traction, auto &t_normalized_normal,
533 auto &t_kappa, auto &t_delta_kappa, auto gc) {
534 double kappa = static_cast<double>(t_kappa + t_delta_kappa);
538 t_traction, t_normalized_normal, alpha, OpGetParameters::beta,
539 true);
540 FTENSOR_INDEX(3, i);
541 FTENSOR_INDEX(3, j);
542 FTensor::Tensor2<double, 3, 3> t_dgap_double;
543 t_dgap_double(i, j) = -t_dgap(i, j);
544 return t_dgap_double;
545 };
546
547 auto assemble = [&](const FTensor::Tensor2<double, 3, 3> &&t_dgap) {
548 int rr = 0;
549 for (; rr != nb_dofs / SPACE_DIM; ++rr) {
550 auto t_mat = getFTensor2FromArray<3, 3, 3>(OP::locMat, 3 * rr);
551 double row_base = t_w * (t_row_base_fun(J) * t_normal(J));
552 auto t_col_base_fun = data.getFTensor1N<3>(gg, 0);
553 for (int cc = 0; cc != nb_dofs / SPACE_DIM; ++cc) {
554 double col_base = t_col_base_fun(J) * t_normalized_normal(J);
555 t_mat(i, j) += (row_base * col_base) * t_dgap(i, j);
556 ++t_mat;
557 ++t_col_base_fun;
558 }
559 ++t_row_base_fun;
560 }
561 for (; rr != nb_base_functions; ++rr)
562 ++t_row_base_fun;
563 };
564
565 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
566 *gcPtr));
567
568 next();
569 }
570
572 };
573
576
577 if (!this->timeScalingFun.empty())
578 this->locMat *= this->timeScalingFun(this->getFEMethod()->ts_t);
579 if (!this->feScalingFun.empty())
580 this->locMat *= this->feScalingFun(this->getFEMethod());
581 // assemble local matrix
582 CHKERR this->matSetValuesHook(this, data, data, this->locMat);
583
585 }
586
587private:
588};
589
591
594
597
598 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
599
600 FTENSOR_INDEX(3, i);
601 FTENSOR_INDEX(3, J);
602 FTENSOR_INDEX(3, K);
603
604 int nb_integration_pts = OP::getGaussPts().size2();
605
606 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
607 auto t_lambda = getFTensor2FromMat<3, 3>(lambdaPtr->at(get_sense_index()));
608 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
609 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
610 auto t_face_normal = getFTensor1NormalsAtGaussPts();
611 auto t_w = OP::getFTensor0IntegrationWeight();
612
613 auto next = [&]() {
614 ++t_P;
615 ++t_face_normal;
616 ++t_kappa;
617 ++t_delta_kappa;
618 ++t_lambda;
619 ++t_w;
620 };
621
622 double face_dissipation = 0.0;
623 double face_grad_dissipation = 0.0;
624 auto face_handle = OP::getFEEntityHandle();
625 CHKERR getMoab().tag_get_data(dissipationTag, &face_handle, 1,
626 &face_dissipation);
627 CHKERR getMoab().tag_get_data(gradDissipationTag, &face_handle, 1,
628 &face_grad_dissipation);
629
630 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
632 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
633 FTensor::Tensor1<double, 3> t_normalized_normal;
634 t_normalized_normal(J) = t_normal(J);
635 t_normalized_normal.normalize();
637 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
638
639 auto dJ_dkappa = [](auto &t_delta_kappa, auto &t_traction,
640 auto &t_normalized_normal, auto &t_kappa, auto gc) {
641 double kappa = static_cast<double>(t_kappa);
642 double delta_kappa = static_cast<double>(t_delta_kappa);
644 t_traction, t_normalized_normal, OpGetParameters::beta);
646 delta_kappa, teff, kappa, gc, OpGetParameters::min_stiffness);
647 double m_grad =
649 delta_kappa, teff, kappa, gc, OpGetParameters::min_stiffness);
650 return boost::make_tuple(m, m_grad);
651 };
652
653 auto dr_kappa = [](auto &t_delta_kappa, auto &t_traction,
654 auto &t_normalized_normal, auto &t_kappa, auto gc) {
655 double kappa = static_cast<double>(t_kappa);
656 double delta_kappa = static_cast<double>(t_delta_kappa);
657 double kappa_plus_delta = kappa + delta_kappa;
658 double diff_alpha = GriffithCohesiveLaw::getDiffAlpha(
659 kappa_plus_delta, gc, OpGetParameters::min_stiffness);
660 auto [teff, t_gap] = GriffithCohesiveLaw::calculateGap(
661 t_traction, t_normalized_normal, diff_alpha, OpGetParameters::beta);
662 FTENSOR_INDEX(3, i);
663 FTensor::Tensor1<double, 3> t_gap_double;
664 t_gap_double(i) = -t_gap(i);
665 return t_gap_double;
666 };
667
668 auto [J, dJ] = dJ_dkappa(t_delta_kappa, t_traction, t_normalized_normal,
669 t_kappa, *gcPtr);
670 face_dissipation += t_w * J * t_normal.l2();
671 face_grad_dissipation += t_w * dJ * t_normal.l2();
672
673 auto t_dr = dr_kappa(t_delta_kappa, t_traction, t_normalized_normal,
674 t_kappa, *gcPtr);
676 t_l(i) = t_lambda(i, K) * t_normal(K);
677 face_grad_dissipation -= t_w * t_l(i) * t_dr(i);
678
679 next();
680 }
681
682 CHKERR getMoab().tag_set_data(dissipationTag, &face_handle, 1,
683 &face_dissipation);
684 CHKERR getMoab().tag_set_data(gradDissipationTag, &face_handle, 1,
685 &face_grad_dissipation);
686
687
689 };
690
695
696protected:
697};
698
700
703
706
707 FTENSOR_INDEX(3, i);
708 FTENSOR_INDEX(3, J);
709
710 int nb_integration_pts = OP::getGaussPts().size2();
711 auto nb_dofs = data.getIndices().size();
712
713 int nb_base_functions = data.getN().size2() / 3;
714 auto t_row_base_fun = data.getFTensor1N<3>();
715 auto t_w = OP::getFTensor0IntegrationWeight();
716
717 auto t_P = getFTensor2FromMat<3, 3>(*(OP::fluxMatPtr));
718 auto t_kappa = getFTensor0FromVec<0>(*kappaPtr);
719 auto t_delta_kappa = getFTensor0FromVec<0>(*kappaDeltaPtr);
720 auto t_face_normal = getFTensor1NormalsAtGaussPts();
721
722 auto next = [&]() {
723 ++t_P;
724 ++t_face_normal;
725 ++t_kappa;
726 ++t_delta_kappa;
727 ++t_w;
728 };
729
730 for (auto gg = 0; gg != nb_integration_pts; ++gg) {
732 t_normal(J) = t_face_normal(J) * (faceSense / 2.);
733 FTensor::Tensor1<double, 3> t_normalized_normal;
734 t_normalized_normal(J) = t_normal(J);
735 t_normalized_normal.normalize();
737 t_traction(i) = t_P(i, J) * t_normalized_normal(J);
738
739 auto dJ_dtraction = [](auto &t_traction, auto &t_normalized_normal,
740 auto &t_delta_kappa, auto &t_kappa, auto gc) {
741 double kappa = static_cast<double>(t_kappa);
742 double delta_kappa = static_cast<double>(t_delta_kappa);
743 FTENSOR_INDEX(3, i);
744 FTensor::Tensor1<double, 3> t_traction_double;
745 t_traction_double(i) = t_traction(i);
746 auto t_dM =
748 t_traction_double, delta_kappa, kappa, t_normalized_normal, gc,
750 FTensor::Tensor1<double, 3> t_dJ_double;
751 t_dJ_double(i) = t_dM(i);
752 return t_dJ_double;
753 };
754
755 auto assemble = [&](const FTensor::Tensor1<double, 3> &&t_dJ) {
756
757 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
758 int bb = 0;
759 for (; bb != nb_dofs / SPACE_DIM; ++bb) {
760 double row_base = t_w * (t_row_base_fun(J) * t_normal(J));
761 t_nf(i) += row_base * t_dJ(i);
762 ++t_nf;
763 ++t_row_base_fun;
764 }
765 for (; bb != nb_base_functions; ++bb)
766 ++t_row_base_fun;
767 };
768
769 assemble(dJ_dtraction(t_traction, t_normalized_normal, t_delta_kappa,
770 t_kappa, *gcPtr));
771
772 // cerr << "AAAA " << OP::locF << endl;
773 next();
774 }
775
777 };
778
781 return VecSetValues<AssemblyTypeSelector<A>>(vec_dJdu, data, OP::locF,
782 ADD_VALUES);
784 }
785
786protected:
787};
788
790
792
794
796 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
797 Tag tag, TagGetType tag_get_type,
798 boost::shared_ptr<VectorDouble> tag_data_ptr,
799 boost::shared_ptr<Range> ents_ptr = nullptr)
800 : OpBrokenBaseCohesive(broken_base_side_data, ents_ptr), tagHandle(tag),
801 tagGetType(tag_get_type), tagDataPtr(tag_data_ptr) {}
802
805
806 auto get_sense_index = [this]() { return (faceSense == 1) ? 0 : 1; };
807 auto &v = *tagDataPtr;
808
809 switch (tagGetType) {
810 case TagSET:
811 break;
812 case TagGET:
813 v.resize(1);
814 v.clear();
815 break;
816 }
817
818 auto get_data = [&](auto &v) {
820
821 auto &moab = getMoab();
822 auto fe_ent = getFEEntityHandle();
823
824 int size;
825 double *data;
826 rval = moab.tag_get_by_ptr(tagHandle, &fe_ent, 1, (const void **)&data,
827 &size);
828 if (
829
830 rval == MB_SUCCESS && size > 0 && size != v.size()
831
832 ) {
833 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
834 "Inconsistent size of tag data");
835 } else {
836 if (rval != MB_SUCCESS || tagGetType == TagSET) {
837 int tag_size[1];
838 tag_size[0] = v.size();
839 void const *tag_data[] = {&v[0]};
840 CHKERR moab.tag_set_by_ptr(tagHandle, &fe_ent, 1, tag_data, tag_size);
841 } else {
842 CHKERR moab.tag_get_by_ptr(tagHandle, &fe_ent, 1,
843 (const void **)&data, &size);
844 std::copy(data, data + size, v.begin());
845 }
846 }
847
849 };
850
851 CHKERR get_data(v);
852
854 }
855
860
861private:
863 boost::shared_ptr<VectorDouble> tagDataPtr;
865};
866
868 EshelbianCore &ep,
869 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
870 SmartPetscObj<Vec> lambda_vec = SmartPetscObj<Vec>()) {
871
872 using EleOnSide =
874 using SideEleOp = EleOnSide::UserDataOperator;
875
876 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
877 pip, {HDIV, H1, L2}, ep.materialH1Positions, ep.frontAdjEdges);
878
879 auto domain_side_flux = [&](auto &pip) {
880 // flux
881 auto broken_data_ptr =
882 boost::make_shared<std::vector<BrokenBaseSideData>>();
884 broken_data_ptr));
885 auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
887 ep.piolaStress, flux_mat_ptr));
888 pip.push_back(new OpSetFlux<SideEleOp>(broken_data_ptr, flux_mat_ptr));
889
890 return broken_data_ptr;
891 };
892
893 auto get_lambda = [&](auto &pip) {
894 boost::shared_ptr<std::array<MatrixDouble, 2>> array_lambda_ptr;
895 if (lambda_vec) {
896 array_lambda_ptr = boost::make_shared<std::array<MatrixDouble, 2>>();
897 auto lambda_mat_ptr = boost::make_shared<MatrixDouble>();
899 ep.piolaStress, lambda_mat_ptr, boost::make_shared<double>(1.0),
900 lambda_vec));
902 auto op = new OP(NOSPACE, OP::OPSPACE);
903 op->doWorkRhsHook =
904 [lambda_mat_ptr, array_lambda_ptr](
905 DataOperator *base_op_ptr, int side, EntityType type,
908 auto op_ptr = static_cast<OP *>(base_op_ptr);
909 auto get_sense_index = [op_ptr]() {
910 return (op_ptr->getSkeletonSense() == 1) ? 0 : 1;
911 };
912 array_lambda_ptr->at(get_sense_index()) = *lambda_mat_ptr;
914 };
915 pip.push_back(op);
916 }
917 return array_lambda_ptr;
918 };
919
920 return std::make_pair(domain_side_flux(pip), get_lambda(pip));
921};
922
924 EshelbianCore &ep,
925 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
926 boost::shared_ptr<Range> interface_range_ptr,
927 SmartPetscObj<Vec> lambda_vec = SmartPetscObj<Vec>()) {
928
929 auto &m_field = ep.mField;
930
931 using BoundaryEle =
933 using EleOnSide =
935 using SideEleOp = EleOnSide::UserDataOperator;
936 using BdyEleOp = BoundaryEle::UserDataOperator;
937
938 auto face_side = [&]() {
939 // First: Iterate over skeleton FEs adjacent to Domain FEs
940 // Note: BoundaryEle, i.e. uses skeleton interation rule
941 auto op_loop_skeleton_side =
943 interface_range_ptr, Sev::noisy);
944 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = [](int, int, int) {
945 return -1;
946 };
947 op_loop_skeleton_side->getSideFEPtr()->setRuleHook =
948 set_integration_at_front_face;
949
950 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
951 op_loop_skeleton_side->getOpPtrVector(), {L2}, ep.materialH1Positions,
952 ep.frontAdjEdges);
953
954 return op_loop_skeleton_side;
955 };
956
957 auto op_loop_skeleton_side = face_side();
958
959 // Second: Iterate over domain FEs adjacent to skelton, particularly one
960 // domain element.
961 // Note: EleOnSide, i.e. uses on domain projected skeleton rule
962 auto op_loop_domain_side = new OpBrokenLoopSide<EleOnSide>(
963 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
964 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
965 boost::make_shared<CGGUserPolynomialBase>(nullptr, true);
966 auto [broken_data_flux_ptr, lambda_ptr] = pushCohesiveOpsDomainImpl(
967 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
968 auto op_reset_broken_data = new ForcesAndSourcesCore::UserDataOperator(
970 op_reset_broken_data->doWorkRhsHook =
971 [broken_data_flux_ptr](DataOperator *base_op_ptr, int side,
972 EntityType type,
974 broken_data_flux_ptr->resize(0);
975 return 0;
976 };
977 op_loop_skeleton_side->getOpPtrVector().push_back(op_reset_broken_data);
978 op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
979
980 auto kappa_ptr = boost::make_shared<VectorDouble>();
981 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
982 auto gc_ptr = boost::make_shared<double>();
983
984 op_loop_skeleton_side->getOpPtrVector().push_back(
985 new OpGetParameters(gc_ptr));
986 op_loop_skeleton_side->getOpPtrVector().push_back(
987 new OpGetTag(broken_data_flux_ptr, get_kappa_tag(m_field.get_moab()),
988 TagGetType::TagGET, kappa_ptr));
989 op_loop_skeleton_side->getOpPtrVector().push_back(new OpGetTag(
990 broken_data_flux_ptr, get_delta_kappa_tag(m_field.get_moab()),
991 TagGetType::TagGET, kappa_delta_ptr));
992
993 return std::make_tuple(op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr,
994 kappa_ptr, kappa_delta_ptr);
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
1009 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1010 op_loop_skeleton_side->getOpPtrVector().push_back(
1012 u_gamma_ptr));
1013 op_loop_skeleton_side->getOpPtrVector().push_back(new OpCohesiveRhs(
1014 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1015
1016 pip.push_back(op_loop_skeleton_side);
1017
1019};
1020
1022 EshelbianCore &ep,
1023 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1024 boost::shared_ptr<Range> interface_range_ptr,
1025 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1027
1028 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
1029 kappa_delta_ptr] =
1030 pushCohesiveOpsImpl(ep, set_integration_at_front_face,
1031 interface_range_ptr);
1032 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1033 op_loop_skeleton_side->getOpPtrVector().push_back(
1035 u_gamma_ptr));
1036 op_loop_skeleton_side->getOpPtrVector().push_back(new OpCohesiveLhs_dPdP(
1037 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1038
1039 pip.push_back(op_loop_skeleton_side);
1040
1042};
1043
1045 EshelbianCore &ep,
1046 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1047 SmartPetscObj<Vec> lambda_vec) {
1048 auto &m_field = ep.mField;
1049
1050 using EleOnSide =
1052 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1053 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
1054 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1055 boost::make_shared<CGGUserPolynomialBase>(nullptr, true);
1056 auto [broken_data_flux_ptr, lambda_ptr] = pushCohesiveOpsDomainImpl(
1057 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
1058
1059 auto op_reset_broken_data = new ForcesAndSourcesCore::UserDataOperator(
1061 op_reset_broken_data->doWorkRhsHook =
1062 [broken_data_flux_ptr](DataOperator *base_op_ptr, int side,
1063 EntityType type,
1065 broken_data_flux_ptr->resize(0);
1066 return 0;
1067 };
1068 pip.push_back(op_reset_broken_data);
1069 pip.push_back(op_loop_domain_side);
1070
1071 auto tag_dissipation = get_tag(m_field.get_moab(), "COHESIVE_DISSIPATION", 1);
1072 auto tag_grad_dissipation =
1073 get_tag(m_field.get_moab(), "COHESIVE_DISSIPATION_GRAD", 1);
1074
1075 auto kappa_ptr = boost::make_shared<VectorDouble>();
1076 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
1077 auto gc_ptr = boost::make_shared<double>();
1078
1079 pip.push_back(new OpGetParameters(gc_ptr, Sev::noisy));
1080 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1081 get_kappa_tag(m_field.get_moab()),
1082 TagGetType::TagGET, kappa_ptr));
1083 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1084 get_delta_kappa_tag(m_field.get_moab()),
1085 TagGetType::TagGET, kappa_delta_ptr));
1086
1087 pip.push_back(new OpCohesive_dJ_dkappa(
1088 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr, lambda_ptr,
1089 tag_dissipation, tag_grad_dissipation));
1090
1091 return std::make_pair(tag_dissipation, tag_grad_dissipation);
1092}
1093
1095 EshelbianCore &ep,
1096 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1097 auto &m_field = ep.mField;
1098
1099 using EleOnSide =
1101 auto op_loop_domain_side = new OpLoopSide<EleOnSide>(
1102 m_field, ep.elementVolumeName, SPACE_DIM, Sev::noisy);
1103 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1104 boost::make_shared<CGGUserPolynomialBase>(nullptr, true);
1105 auto [broken_data_flux_ptr, lambda_ptr] =
1106 pushCohesiveOpsDomainImpl(ep, op_loop_domain_side->getOpPtrVector());
1107
1108 auto op_reset_broken_data = new ForcesAndSourcesCore::UserDataOperator(
1110 op_reset_broken_data->doWorkRhsHook =
1111 [broken_data_flux_ptr](DataOperator *base_op_ptr, int side,
1112 EntityType type,
1114 broken_data_flux_ptr->resize(0);
1115 return 0;
1116 };
1117 pip.push_back(op_reset_broken_data);
1118 pip.push_back(op_loop_domain_side);
1119
1120 auto kappa_ptr = boost::make_shared<VectorDouble>();
1121 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
1122 auto gc_ptr = boost::make_shared<double>();
1123
1124 pip.push_back(new OpGetParameters(gc_ptr, Sev::noisy));
1125 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1126 get_kappa_tag(m_field.get_moab()),
1127 TagGetType::TagGET, kappa_ptr));
1128 pip.push_back(new OpGetTag(broken_data_flux_ptr,
1129 get_delta_kappa_tag(m_field.get_moab()),
1130 TagGetType::TagGET, kappa_delta_ptr));
1131 auto dJ_dx_vec = createDMVector(ep.dmElastic);
1132 pip.push_back(new OpCohesive_dJ_dP(broken_data_flux_ptr, gc_ptr, kappa_ptr,
1133 kappa_delta_ptr, lambda_ptr, Tag(), Tag(),
1134 dJ_dx_vec));
1135
1136 return dJ_dx_vec;
1137}
1138
1140 EshelbianCore &ep,
1141 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1142 SmartPetscObj<Vec> lambda_vec,
1143 CommInterface::EntitiesPetscVector &dissipation_vec,
1144 CommInterface::EntitiesPetscVector &grad_dissipation_vec) {
1146
1147 using BoundaryEle =
1149 using BdyEleOp = BoundaryEle::UserDataOperator;
1150
1151 auto get_face_ele = [&]() {
1152 auto fe_ptr = boost::make_shared<BoundaryEle>(ep.mField);
1153 fe_ptr->getRuleHook = [](int, int, int) { return -1; };
1154 fe_ptr->setRuleHook = set_integration_at_front_face;
1155 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
1156 fe_ptr->getOpPtrVector(), {L2}, ep.materialH1Positions,
1157 ep.frontAdjEdges);
1158
1159 auto interface_face = [&](FEMethod *fe_method_ptr) {
1160 auto ent = fe_method_ptr->getFEEntityHandle();
1161 if (
1162
1163 ep.interfaceFaces->find(ent) != ep.interfaceFaces->end()
1164
1165 ) {
1166 return true;
1167 };
1168 return false;
1169 };
1170
1171 fe_ptr->exeTestHook = interface_face;
1172
1173 return fe_ptr;
1174 };
1175
1176 auto face_fe = get_face_ele();
1177 auto [tag_dissipation, tag_grad_dissipation] =
1178 pushCohesive_dJ_dkappa_Impl(ep, face_fe->getOpPtrVector(), lambda_vec);
1179
1180 constexpr double zero = 0.0;
1181 CHKERR ep.mField.get_moab().tag_clear_data(tag_dissipation,
1182 *(ep.interfaceFaces), &zero);
1183 CHKERR ep.mField.get_moab().tag_clear_data(tag_grad_dissipation,
1184 *(ep.interfaceFaces), &zero);
1187 ep.mField.get_moab(), dissipation_vec, tag_dissipation);
1189 ep.mField.get_moab(), grad_dissipation_vec, tag_grad_dissipation);
1190
1192}
1193
1195 EshelbianCore &ep,
1196 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1197 SmartPetscObj<KSP> ksp, SmartPetscObj<Vec> lambda_vec) {
1199
1200 using BoundaryEle =
1202 using BdyEleOp = BoundaryEle::UserDataOperator;
1203
1204 auto get_face_ele = [&]() {
1205 auto fe_ptr = boost::make_shared<BoundaryEle>(ep.mField);
1206 fe_ptr->getRuleHook = [](int, int, int) { return -1; };
1207 fe_ptr->setRuleHook = set_integration_at_front_face;
1208 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
1209 fe_ptr->getOpPtrVector(), {L2}, ep.materialH1Positions,
1210 ep.frontAdjEdges);
1211
1212 auto interface_face = [&](FEMethod *fe_method_ptr) {
1213 auto ent = fe_method_ptr->getFEEntityHandle();
1214 if (
1215
1216 ep.interfaceFaces->find(ent) != ep.interfaceFaces->end()
1217
1218 ) {
1219 return true;
1220 };
1221 return false;
1222 };
1223
1224 fe_ptr->exeTestHook = interface_face;
1225
1226 return fe_ptr;
1227 };
1228
1229 auto face_fe = get_face_ele();
1230 auto dJ_dx = pushCohesive_dJ_dx_Impl(ep, face_fe->getOpPtrVector());
1232 ep.dmElastic, ep.skeletonElement, face_fe, 0, ep.mField.get_comm_size());
1233 CHKERR VecAssemblyBegin(dJ_dx);
1234 CHKERR VecAssemblyEnd(dJ_dx);
1235 CHKERR VecGhostUpdateBegin(dJ_dx, ADD_VALUES, SCATTER_REVERSE);
1236 CHKERR VecGhostUpdateEnd(dJ_dx, ADD_VALUES, SCATTER_REVERSE);
1237 CHKERR VecGhostUpdateBegin(dJ_dx, INSERT_VALUES, SCATTER_FORWARD);
1238 CHKERR VecGhostUpdateEnd(dJ_dx, INSERT_VALUES, SCATTER_FORWARD);
1239 double dJ_dx_norm2;
1240 CHKERR VecNorm(dJ_dx, NORM_2, &dJ_dx_norm2);
1241 MOFEM_LOG("EP", Sev::inform)
1242 << "evaluateCohesiveLambdaImpl: Norm of dJ/dx vector: " << dJ_dx_norm2;
1243 constexpr double tol = 1e-16;
1244 if (dJ_dx_norm2 < tol) {
1245 CHKERR VecZeroEntries(lambda_vec);
1246 } else {
1247 CHKERR KSPSolveTranspose(ksp, dJ_dx, lambda_vec);
1248 }
1249 CHKERR VecGhostUpdateBegin(lambda_vec, INSERT_VALUES, SCATTER_FORWARD);
1250 CHKERR VecGhostUpdateEnd(lambda_vec, INSERT_VALUES, SCATTER_FORWARD);
1251 double lambda_norm2;
1252 CHKERR VecNorm(lambda_vec, NORM_2, &lambda_norm2);
1253 MOFEM_LOG("EP", Sev::inform)
1254 << "evaluateCohesiveLambdaImpl: Norm of lambda vector: " << lambda_norm2;
1255
1257}
1258
1263
1264 CHKERR TSSetFromOptions(ts);
1265 CHKERR TSSetStepNumber(ts, 0);
1266 CHKERR TSSetTime(ts, 0);
1267 CHKERR TSSetSolution(ts, x);
1268 double dt;
1269 CHKERR TSGetTimeStep(ts, &dt);
1270 MOFEM_LOG("EP", Sev::inform)
1271 << "evaluatePrimalProblemCohesiveImpl: Time step dt: " << dt;
1272 CHKERR TSSolve(ts, PETSC_NULLPTR);
1273 CHKERR TSSetTimeStep(ts, dt);
1274
1275 CHKERR DMoFEMMeshToLocalVector(ep.dmElastic, x, INSERT_VALUES,
1276 SCATTER_FORWARD);
1277 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
1278 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
1279
1280 double x_norm2;
1281 CHKERR VecNorm(x, NORM_2, &x_norm2);
1282 MOFEM_LOG("EP", Sev::inform)
1283 << "evaluatePrimalProblemCohesiveImpl: Norm of displacement vector: "
1284 << x_norm2;
1285
1287}
1288
1291 EshelbianCore *ep,
1292 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1293 SmartPetscObj<TS> time_solver)
1294 : ep_ptr(ep), setIntegrationAtFrontFace(set_integration_at_front_face),
1295 timeSolver(time_solver),
1296
1297 kspSolVec(createDMVector(ep->dmElastic)),
1298 lambdaVec(createDMVector(ep->dmElastic)),
1299 kappaVec(CommInterface::createEntitiesPetscVector(
1300 ep->mField.get_comm(), ep->mField.get_moab(),
1301 [&ep](Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1302 Sev::noisy)),
1304 ep->mField.get_comm(), ep->mField.get_moab(),
1305 [&ep](Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1306 Sev::noisy)),
1308 ep->mField.get_comm(), ep->mField.get_moab(),
1309 [&ep](Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1310 Sev::noisy)) {}
1311
1313 return vectorDuplicate(kappaVec.second);
1314 }
1315
1319
1323
1327
1328private:
1332
1339 PetscReal *f,
1340 Vec g, void *ctx);
1341};
1342
1343boost::shared_ptr<CohesiveTAOCtx> createCohesiveTAOCtx(
1344 EshelbianCore *ep_ptr,
1345 ForcesAndSourcesCore::GaussHookFun set_integration_at_front_face,
1346 SmartPetscObj<TS> time_solver) {
1347 return boost::make_shared<CohesiveTAOCtxImpl>(
1348 ep_ptr, set_integration_at_front_face, time_solver);
1349}
1350
1352 PetscReal *f, Vec g,
1353 void *ctx) {
1354
1356 auto cohesive_ctx = static_cast<CohesiveTAOCtxImpl *>(ctx);
1357 auto &ep = *(cohesive_ctx->ep_ptr);
1358
1359#ifndef NDEBUG
1360 CHKERR VecView(delta_kappa, PETSC_VIEWER_STDOUT_WORLD);
1361#endif
1362 // Set delta_kappa values to the tag
1363 CHKERR VecCopy(delta_kappa, cohesive_ctx->kappaVec.second);
1364 CHKERR VecGhostUpdateBegin(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1365 SCATTER_FORWARD);
1366 CHKERR VecGhostUpdateEnd(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1367 SCATTER_FORWARD);
1369 ep.mField.get_moab(), cohesive_ctx->kappaVec,
1371
1372 // solve primal problem
1373 CHKERR evaluatePrimalProblemCohesiveImpl(ep, cohesive_ctx->timeSolver,
1374 cohesive_ctx->kspSolVec);
1375 // solve adjoint problem to get lambda
1376 auto ksp = snesGetKSP(tsGetSNES(cohesive_ctx->timeSolver));
1377 CHKERR evaluateCohesiveLambdaImpl(ep, cohesive_ctx->setIntegrationAtFrontFace,
1378 ksp, cohesive_ctx->lambdaVec);
1379 // evaluate dissipation and its gradient
1381 ep, cohesive_ctx->setIntegrationAtFrontFace, cohesive_ctx->lambdaVec,
1382 cohesive_ctx->dissipationVec, cohesive_ctx->gradDissipationVec);
1383
1384 CHKERR VecSum(cohesive_ctx->dissipationVec.second, f);
1385 CHKERR VecCopy(cohesive_ctx->gradDissipationVec.second, g);
1386
1387 MOFEM_LOG("EP", Sev::inform)
1388 << "Cohesive objective function (negative total dissipation): " << *f;
1389
1391}
1392
1401
1402}; // namespace EshelbianPlasticity
Eshelbian plasticity interface.
#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, RowColData rc=RowColData::COL)
set local (or ghosted) vector values on mesh for partition only
Definition DMMoFEM.cpp:514
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:576
auto createDMVector(DM dm, RowColData rc=RowColData::COL)
Get smart vector from DM.
Definition DMMoFEM.hpp:1237
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
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 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
auto snesGetKSP(SNES snes)
auto tsGetSNES(TS ts)
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
tools::promote_args< ArgumentType, int >::type lambert_w(const ArgumentType &z)
Definition lambert_w.hpp:40
constexpr double g
FTensor::Index< 'm', 3 > m
boost::shared_ptr< Range > frontAdjEdges
const std::string skeletonElement
MoFEM::Interface & mField
const std::string materialH1Positions
std::vector< Tag > listTagsToTransfer
list of tags to transfer to postprocessor
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, double min_stiffness)
static T getAlpha(const T &kappa, double gc, double min_stiffness)
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, double min_stiffness)
static auto calculateY(const T t_eff, const T &kappa, double gc, double min_stiffness)
static auto calculateGap(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double alpha, double beta, bool sign_sensitive=true)
static auto invTau(const T &tau, double Gf)
static auto getDiffTau(const T &tau, const T &kappa, double gc)
static auto calculateEffectiveTraction(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double beta)
static auto calculateDiffGapDTraction(const FTensor::Tensor1< T, 3 > &t_traction, FTensor::Tensor1< double, 3 > &t_n_normalize, double alpha, double beta, bool sign_sensitive=true)
static T getDiffAlpha(const T &kappa, double gc, double min_stiffness)
static T getTau(const T &k, double gc)
static auto calculateDissipation(const T &delta_kappa, const T t_eff, const T &kappa, double gc, double min_stiffness)
static auto calculateDissipationSurplusDiffKappa(const T &delta_kappa, const T &t_eff, const T &kappa, double gc, double min_stiffness)
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 iNtegrate(EntData &data)
boost::shared_ptr< double > tatalDissipationGrad
MoFEMErrorCode iNtegrate(EntData &data)
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)
MoFEMErrorCode iNtegrate(EntData &data)
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)
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)
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 MoFEMErrorCode setTagFromVector(moab::Interface &moab, EntitiesPetscVector &vec, Tag tag)
std::pair< std::pair< Range, Range >, SmartPetscObj< Vec > > EntitiesPetscVector
static EntitiesPetscVector createEntitiesPetscVector(MPI_Comm comm, moab::Interface &moab, std::function< Range(Range)> get_entities_fun, const int nb_coeffs, Sev sev=Sev::verbose, int root_rank=0, bool get_vertices=true)
Create a ghost vector for exchanging data.
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
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
double tol