19#include <boost/math/constants/constants.hpp>
28 static inline constexpr double eps = 1e-8;
30 static inline constexpr double pi = boost::math::constants::pi<double>();
32 boost::math::constants::root_pi<double>();
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);
41 return -tau * (1.0 + (1.0 / (2.0 *
kappa)));
53 return (tau -
kappa * diff_tau) / (tau * tau);
56 template <
typename T>
static auto invTau(
const T &tau,
double Gf) {
57 const T z = (2.0 * Gf * Gf) / (
pi * tau * tau);
64 return 0.5 * t_eff * t_eff * diff_alpha;
69 const T &
kappa,
double gc) {
71 return Y * delta_kappa;
76 const T &
kappa,
double gc) {
79 return -M * delta_kappa;
85 const T &
kappa,
double gc) {
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 * 1
i;
91 std::complex<T> cpx_M =
93 return cpx_M.imag() /
eps;
103 std::complex<T> cpx_delta = delta_kappa;
104 std::complex<T> cpx_kappa =
kappa;
106 for (
auto jj = 0; jj != 3; ++jj) {
107 t_cpx_traction(
i) = t_traction(
i);
108 t_cpx_traction(jj) +=
eps * 1
i;
109 std::complex<T> cpx_teff =
111 std::complex<T> cpx_M =
113 t_diff_traction(jj) = cpx_M.imag() /
eps;
115 return t_diff_traction;
118 template <
typename T>
121 double alpha,
double beta) {
126 t_P(
i,
j) = t_n_normalize(
i) * t_n_normalize(
j);
130 t_normal(
i) = t_P(
i,
j) * t_traction(
j);
132 t_tangential(
i) = t_Q(
i,
j) * t_traction(
j);
133 T teff = (t_n_normalize(
i) * t_normal(
i) / 2.) +
134 std::sqrt((t_normal(
i) * t_normal(
i)) / 2. +
135 (1. / beta) * t_tangential(
i) * t_tangential(
i));
139 (teff * (t_n_normalize(
j) * t_P(
j,
i))) +
140 t_normal(
j) * t_P(
j,
i) +
141 (1. / beta) * t_tangential(
j) * t_Q(
j,
i)
144 return std::make_pair(teff, t_gap);
147 template <
typename T>
151 double alpha,
double beta) {
156 for (
auto jj = 0; jj != 3; ++jj) {
157 t_cpx_traction(
i) = t_traction(
i);
158 t_cpx_traction(jj) +=
eps * 1
i;
159 auto [teff_cpx, t_cpx_gap] =
160 calculateGap(t_cpx_traction, t_n_normalize, alpha, beta);
161 for (
auto ii = 0; ii != 3; ++ii) {
162 auto v = t_cpx_gap(ii).imag();
163 t_dgap(ii, jj) =
v /
eps;
169 template <
typename T>
174 auto [teff, t_gap] =
calculateGap(t_traction, t_n_normalize, 1.0, beta);
203 PetscOptionsBegin(PETSC_COMM_WORLD,
"interface_",
"",
"none");
205 CHKERR PetscOptionsScalar(
"-gc",
"Griffith energy release rate",
"",
207 CHKERR PetscOptionsScalar(
"-min_stiffness",
"Minimal interface stiffness",
209 CHKERR PetscOptionsScalar(
"-strength",
"Strength of interface",
"",
211 CHKERR PetscOptionsScalar(
"-min_kappa",
212 "Minimal kappa to avoid singularity",
"",
214 CHKERR PetscOptionsScalar(
"-kappa0",
"Characteristic length kappa0",
"",
216 CHKERR PetscOptionsScalar(
"-beta",
"Cohesive tangential coupling",
"",
222 <<
"Interface Griffith energy release rate Gc -interface_gc = "
225 <<
"Interface min stiffness -interface_min_stiffness = "
228 <<
"Interface strength -interface_strength = " <<
strength;
230 <<
"Interface minimal kappa -interface_min_kappa = " <<
min_kappa;
232 <<
"Interface characteristic length kappa0 -interface_kappa0 = "
236 double kappa_strength =
240 <<
"Adjusted interface min kappa to capture strength = " <<
min_kappa;
254 const int def_size = 1;
256 moab.tag_get_handle(tag_name.c_str(), def_size, MB_TYPE_DOUBLE, tag,
257 MB_TAG_CREAT | MB_TAG_VARLEN | MB_TAG_SPARSE, NULL);
258 if (
rval == MB_ALREADY_ALLOCATED)
264static Tag get_tag(moab::Interface &moab, std::string tag_name,
int size) {
266 rval = moab.tag_get_handle(tag_name.c_str(), size, MB_TYPE_DOUBLE, tag,
267 MB_TAG_CREAT | MB_TAG_SPARSE, NULL);
268 if (
rval == MB_ALREADY_ALLOCATED)
275 return get_tag(moab,
"KAPPA", 1);
279 return get_tag(moab,
"DELTA_KAPPA", 1);
288 Assembly<A>::OpBrokenBase;
290 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_flux_data_ptr,
291 boost::shared_ptr<Range> ents_ptr =
nullptr)
292 :
OP(broken_flux_data_ptr, ents_ptr) {}
300 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
305 if (!brokenBaseSideData) {
310 auto do_work_rhs = [
this](
int row_side, EntityType row_type,
314 auto base = row_data.
getBase();
317 "row base not set properly");
326 OP::nbIntegrationPts = OP::getGaussPts().size2();
328 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
330 OP::locF.resize(OP::nbRows,
false);
332 OP::locMat.resize(OP::nbRows, OP::nbRows,
false);
338 CHKERR this->aSsemble(row_data);
343 switch (OP::opType) {
345 for (
auto &bd : *brokenBaseSideData) {
347 boost::shared_ptr<MatrixDouble>(brokenBaseSideData, &bd.getFlux());
349 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData());
356 (std::string(
"wrong op type ") +
357 OpBaseDerivativesBase::OpTypeNames[OP::opType])
373 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
374 boost::shared_ptr<double> gc_ptr,
375 boost::shared_ptr<VectorDouble> kappa_ptr,
376 boost::shared_ptr<VectorDouble> kappa_delta_ptr,
377 boost::shared_ptr<std::array<MatrixDouble, 2>> lambda_ptr =
nullptr,
378 Tag dissipation_tags =
Tag(),
Tag grad_dissipation_tags =
Tag(),
380 boost::shared_ptr<Range> ents_ptr =
nullptr)
389 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
395 int nb_integration_pts = OP::getGaussPts().size2();
398 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
400 auto t_face_normal = getFTensor1NormalsAtGaussPts();
402 int nb_base_functions = data.
getN().size2() / 3;
404 auto t_w = OP::getFTensor0IntegrationWeight();
414 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
418 t_normalized_normal(
J) = t_normal(
J);
421 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
423 auto fracture = [
this](
auto &t_traction,
auto &t_normalized_normal,
424 auto &t_kappa,
auto &t_delta_kappa,
auto gc) {
425 double kappa =
static_cast<double>(t_kappa + t_delta_kappa);
432 t_gap_double(
i) = -t_gap(
i);
437 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
439 for (; bb != nb_dofs /
SPACE_DIM; ++bb) {
440 double row_base = t_w * (t_row_base_fun(
J) * t_normal(
J));
441 t_nf(
i) += row_base * t_gap(
i);
445 for (; bb != nb_base_functions; ++bb)
449 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
463 boost::shared_ptr<std::array<MatrixDouble, 2>>
lambdaPtr;
479 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
486 int nb_integration_pts = OP::getGaussPts().size2();
489 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
491 auto t_face_normal = getFTensor1NormalsAtGaussPts();
493 int nb_base_functions = data.
getN().size2() / 3;
495 auto t_w = OP::getFTensor0IntegrationWeight();
505 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
509 t_normalized_normal(
J) = t_normal(
J);
512 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
514 auto fracture = [
this](
auto &t_traction,
auto &t_normalized_normal,
515 auto &t_kappa,
auto &t_delta_kappa,
auto gc) {
516 double kappa =
static_cast<double>(t_kappa + t_delta_kappa);
524 t_dgap_double(
i,
j) = -t_dgap(
i,
j);
525 return t_dgap_double;
530 for (; rr != nb_dofs /
SPACE_DIM; ++rr) {
531 auto t_mat = getFTensor2FromArray<3, 3, 3>(OP::locMat, 3 * rr);
532 double row_base = t_w * (t_row_base_fun(
J) * t_normal(
J));
534 for (
int cc = 0; cc != nb_dofs /
SPACE_DIM; ++cc) {
535 double col_base = t_col_base_fun(
J) * t_normalized_normal(
J);
536 t_mat(
i,
j) += (row_base * col_base) * t_dgap(
i,
j);
542 for (; rr != nb_base_functions; ++rr)
546 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
558 if (!this->timeScalingFun.empty())
559 this->locMat *= this->timeScalingFun(this->getFEMethod()->ts_t);
560 if (!this->feScalingFun.empty())
561 this->locMat *= this->feScalingFun(this->getFEMethod());
563 CHKERR this->matSetValuesHook(
this, data, data, this->locMat);
579 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
585 int nb_integration_pts = OP::getGaussPts().size2();
588 auto t_lambda = getFTensor2FromMat<3, 3>(
lambdaPtr->at(get_sense_index()));
589 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
591 auto t_face_normal = getFTensor1NormalsAtGaussPts();
592 auto t_w = OP::getFTensor0IntegrationWeight();
603 double face_dissipation = 0.0;
604 double face_grad_dissipation = 0.0;
605 auto face_handle = OP::getFEEntityHandle();
609 &face_grad_dissipation);
611 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
615 t_normalized_normal(
J) = t_normal(
J);
618 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
620 auto dJ_dkappa = [](
auto &t_delta_kappa,
auto &t_traction,
621 auto &t_normalized_normal,
auto &t_kappa,
auto gc) {
622 double kappa =
static_cast<double>(t_kappa);
623 double delta_kappa =
static_cast<double>(t_delta_kappa);
627 delta_kappa, teff,
kappa, gc);
630 delta_kappa, teff,
kappa, gc);
631 return boost::make_tuple(
m, m_grad);
634 auto dr_kappa = [](
auto &t_delta_kappa,
auto &t_traction,
635 auto &t_normalized_normal,
auto &t_kappa,
auto gc) {
636 double kappa =
static_cast<double>(t_kappa);
637 double delta_kappa =
static_cast<double>(t_delta_kappa);
638 double kappa_plus_delta =
kappa + delta_kappa;
645 t_gap_double(
i) = -t_gap(
i);
649 auto [
J, dJ] = dJ_dkappa(t_delta_kappa, t_traction, t_normalized_normal,
651 face_dissipation += t_w *
J * t_normal.
l2();
652 face_grad_dissipation += t_w * dJ * t_normal.
l2();
654 auto t_dr = dr_kappa(t_delta_kappa, t_traction, t_normalized_normal,
657 t_l(
i) = t_lambda(
i,
K) * t_normal(
K);
658 face_grad_dissipation -= t_w * t_l(
i) * t_dr(
i);
666 &face_grad_dissipation);
691 int nb_integration_pts = OP::getGaussPts().size2();
694 int nb_base_functions = data.
getN().size2() / 3;
696 auto t_w = OP::getFTensor0IntegrationWeight();
699 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
701 auto t_face_normal = getFTensor1NormalsAtGaussPts();
711 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
715 t_normalized_normal(
J) = t_normal(
J);
718 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
720 auto dJ_dtraction = [](
auto &t_traction,
auto &t_normalized_normal,
721 auto &t_delta_kappa,
auto &t_kappa,
auto gc) {
722 double kappa =
static_cast<double>(t_kappa);
723 double delta_kappa =
static_cast<double>(t_delta_kappa);
726 t_traction_double(
i) = t_traction(
i);
729 t_traction_double, delta_kappa,
kappa, t_normalized_normal, gc,
732 t_dJ_double(
i) = t_dM(
i);
738 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
740 for (; bb != nb_dofs /
SPACE_DIM; ++bb) {
741 double row_base = t_w * (t_row_base_fun(
J) * t_normal(
J));
742 t_nf(
i) += row_base * t_dJ(
i);
746 for (; bb != nb_base_functions; ++bb)
750 assemble(dJ_dtraction(t_traction, t_normalized_normal, t_delta_kappa,
762 return VecSetValues<AssemblyTypeSelector<A>>(
vec_dJdu, data, OP::locF,
777 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
779 boost::shared_ptr<VectorDouble> tag_data_ptr,
780 boost::shared_ptr<Range> ents_ptr =
nullptr)
787 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
799 auto get_data = [&](
auto &
v) {
802 auto &moab = getMoab();
803 auto fe_ent = getFEEntityHandle();
807 rval = moab.tag_get_by_ptr(
tagHandle, &fe_ent, 1, (
const void **)&data,
811 rval == MB_SUCCESS && size > 0 && size !=
v.size()
815 "Inconsistent size of tag data");
819 tag_size[0] =
v.size();
820 void const *tag_data[] = {&
v[0]};
824 (
const void **)&data, &size);
825 std::copy(data, data + size,
v.begin());
850 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
855 using SideEleOp = EleOnSide::UserDataOperator;
857 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
860 auto domain_side_flux = [&](
auto &pip) {
862 auto broken_data_ptr =
863 boost::make_shared<std::vector<BrokenBaseSideData>>();
866 auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
871 return broken_data_ptr;
874 auto get_lambda = [&](
auto &pip) {
875 boost::shared_ptr<std::array<MatrixDouble, 2>> array_lambda_ptr;
877 array_lambda_ptr = boost::make_shared<std::array<MatrixDouble, 2>>();
878 auto lambda_mat_ptr = boost::make_shared<MatrixDouble>();
880 ep.
piolaStress, lambda_mat_ptr, boost::make_shared<double>(1.0),
885 [lambda_mat_ptr, array_lambda_ptr](
889 auto op_ptr =
static_cast<OP *
>(base_op_ptr);
890 auto get_sense_index = [op_ptr]() {
891 return (op_ptr->getSkeletonSense() == 1) ? 0 : 1;
893 array_lambda_ptr->at(get_sense_index()) = *lambda_mat_ptr;
898 return array_lambda_ptr;
901 return std::make_pair(domain_side_flux(pip), get_lambda(pip));
907 boost::shared_ptr<Range> interface_range_ptr,
910 auto &m_field = ep.
mField;
916 using SideEleOp = EleOnSide::UserDataOperator;
917 using BdyEleOp = BoundaryEle::UserDataOperator;
919 auto face_side = [&]() {
922 auto op_loop_skeleton_side =
924 interface_range_ptr, Sev::noisy);
925 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = [](int, int, int) {
928 op_loop_skeleton_side->getSideFEPtr()->setRuleHook =
929 set_integration_at_front_face;
931 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
935 return op_loop_skeleton_side;
938 auto op_loop_skeleton_side = face_side();
945 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
946 boost::make_shared<CGGUserPolynomialBase>();
948 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
951 op_reset_broken_data->doWorkRhsHook =
952 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
955 broken_data_flux_ptr->resize(0);
958 op_loop_skeleton_side->getOpPtrVector().push_back(op_reset_broken_data);
959 op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
961 auto kappa_ptr = boost::make_shared<VectorDouble>();
962 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
963 auto gc_ptr = boost::make_shared<double>();
965 op_loop_skeleton_side->getOpPtrVector().push_back(
967 op_loop_skeleton_side->getOpPtrVector().push_back(
970 op_loop_skeleton_side->getOpPtrVector().push_back(
new OpGetTag(
974 return std::make_tuple(op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr,
975 kappa_ptr, kappa_delta_ptr);
981 boost::shared_ptr<Range> interface_range_ptr,
982 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
985 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
988 interface_range_ptr);
990 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
991 op_loop_skeleton_side->getOpPtrVector().push_back(
994 op_loop_skeleton_side->getOpPtrVector().push_back(
new OpCohesiveRhs(
995 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
997 pip.push_back(op_loop_skeleton_side);
1005 boost::shared_ptr<Range> interface_range_ptr,
1006 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1009 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
1012 interface_range_ptr);
1013 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1014 op_loop_skeleton_side->getOpPtrVector().push_back(
1018 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1020 pip.push_back(op_loop_skeleton_side);
1027 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1029 auto &m_field = ep.
mField;
1035 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1036 boost::make_shared<CGGUserPolynomialBase>();
1038 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
1042 op_reset_broken_data->doWorkRhsHook =
1043 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
1046 broken_data_flux_ptr->resize(0);
1049 pip.push_back(op_reset_broken_data);
1050 pip.push_back(op_loop_domain_side);
1052 auto tag_dissipation =
get_tag(m_field.get_moab(),
"COHESIVE_DISSIPATION", 1);
1053 auto tag_grad_dissipation =
1054 get_tag(m_field.get_moab(),
"COHESIVE_DISSIPATION_GRAD", 1);
1056 auto kappa_ptr = boost::make_shared<VectorDouble>();
1057 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
1058 auto gc_ptr = boost::make_shared<double>();
1061 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1064 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1069 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr, lambda_ptr,
1070 tag_dissipation, tag_grad_dissipation));
1072 return std::make_pair(tag_dissipation, tag_grad_dissipation);
1077 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1078 auto &m_field = ep.
mField;
1084 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1085 boost::make_shared<CGGUserPolynomialBase>();
1086 auto [broken_data_flux_ptr, lambda_ptr] =
1091 op_reset_broken_data->doWorkRhsHook =
1092 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
1095 broken_data_flux_ptr->resize(0);
1098 pip.push_back(op_reset_broken_data);
1099 pip.push_back(op_loop_domain_side);
1101 auto kappa_ptr = boost::make_shared<VectorDouble>();
1102 auto kappa_delta_ptr = boost::make_shared<VectorDouble>();
1103 auto gc_ptr = boost::make_shared<double>();
1106 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1109 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1113 pip.push_back(
new OpCohesive_dJ_dP(broken_data_flux_ptr, gc_ptr, kappa_ptr,
1114 kappa_delta_ptr, lambda_ptr,
Tag(),
Tag(),
1130 using BdyEleOp = BoundaryEle::UserDataOperator;
1132 auto get_face_ele = [&]() {
1133 auto fe_ptr = boost::make_shared<BoundaryEle>(ep.
mField);
1134 fe_ptr->getRuleHook = [](int, int, int) {
return -1; };
1135 fe_ptr->setRuleHook = set_integration_at_front_face;
1136 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
1140 auto interface_face = [&](
FEMethod *fe_method_ptr) {
1141 auto ent = fe_method_ptr->getFEEntityHandle();
1152 fe_ptr->exeTestHook = interface_face;
1157 auto face_fe = get_face_ele();
1158 auto [tag_dissipation, tag_grad_dissipation] =
1161 constexpr double zero = 0.0;
1170 ep.
mField.
get_moab(), grad_dissipation_vec, tag_grad_dissipation);
1183 using BdyEleOp = BoundaryEle::UserDataOperator;
1185 auto get_face_ele = [&]() {
1186 auto fe_ptr = boost::make_shared<BoundaryEle>(ep.
mField);
1187 fe_ptr->getRuleHook = [](int, int, int) {
return -1; };
1188 fe_ptr->setRuleHook = set_integration_at_front_face;
1189 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
1193 auto interface_face = [&](
FEMethod *fe_method_ptr) {
1194 auto ent = fe_method_ptr->getFEEntityHandle();
1205 fe_ptr->exeTestHook = interface_face;
1210 auto face_fe = get_face_ele();
1214 CHKERR VecAssemblyBegin(dJ_dx);
1215 CHKERR VecAssemblyEnd(dJ_dx);
1216 CHKERR VecGhostUpdateBegin(dJ_dx, ADD_VALUES, SCATTER_REVERSE);
1217 CHKERR VecGhostUpdateEnd(dJ_dx, ADD_VALUES, SCATTER_REVERSE);
1218 CHKERR VecGhostUpdateBegin(dJ_dx, INSERT_VALUES, SCATTER_FORWARD);
1219 CHKERR VecGhostUpdateEnd(dJ_dx, INSERT_VALUES, SCATTER_FORWARD);
1221 CHKERR VecNorm(dJ_dx, NORM_2, &dJ_dx_norm2);
1223 <<
"evaluateCohesiveLambdaImpl: Norm of dJ/dx vector: " << dJ_dx_norm2;
1224 constexpr double tol = 1e-16;
1225 if (dJ_dx_norm2 <
tol) {
1226 CHKERR VecZeroEntries(lambda_vec);
1228 CHKERR KSPSolveTranspose(ksp, dJ_dx, lambda_vec);
1230 CHKERR VecGhostUpdateBegin(lambda_vec, INSERT_VALUES, SCATTER_FORWARD);
1231 CHKERR VecGhostUpdateEnd(lambda_vec, INSERT_VALUES, SCATTER_FORWARD);
1232 double lambda_norm2;
1233 CHKERR VecNorm(lambda_vec, NORM_2, &lambda_norm2);
1235 <<
"evaluateCohesiveLambdaImpl: Norm of lambda vector: " << lambda_norm2;
1245 CHKERR TSSetFromOptions(ts);
1246 CHKERR TSSetStepNumber(ts, 0);
1248 CHKERR TSSetSolution(ts, x);
1252 <<
"evaluatePrimalProblemCohesiveImpl: Time step dt: " <<
dt;
1253 CHKERR TSSolve(ts, PETSC_NULLPTR);
1258 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
1259 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
1262 CHKERR VecNorm(x, NORM_2, &x_norm2);
1264 <<
"evaluatePrimalProblemCohesiveImpl: Norm of displacement vector: "
1281 ep->mField.get_comm(), ep->mField.get_moab(),
1286 [&ep](
Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1290 [&ep](
Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1328 return boost::make_shared<CohesiveTAOCtxImpl>(
1329 ep_ptr, set_integration_at_front_face, time_solver);
1333 PetscReal *f, Vec
g,
1338 auto &ep = *(cohesive_ctx->ep_ptr);
1341 CHKERR VecView(delta_kappa, PETSC_VIEWER_STDOUT_WORLD);
1344 CHKERR VecCopy(delta_kappa, cohesive_ctx->kappaVec.second);
1345 CHKERR VecGhostUpdateBegin(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1347 CHKERR VecGhostUpdateEnd(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1355 cohesive_ctx->kspSolVec);
1359 ksp, cohesive_ctx->lambdaVec);
1362 ep, cohesive_ctx->setIntegrationAtFrontFace, cohesive_ctx->lambdaVec,
1363 cohesive_ctx->dissipationVec, cohesive_ctx->gradDissipationVec);
1365 CHKERR VecSum(cohesive_ctx->dissipationVec.second, f);
1366 CHKERR VecCopy(cohesive_ctx->gradDissipationVec.second,
g);
1369 <<
"Cohesive objective function (negative total dissipation): " << *f;
Eshelbian plasticity interface.
#define FTENSOR_INDEX(DIM, I)
ElementsAndOps< SPACE_DIM >::BoundaryEle BoundaryEle
Tensor1< T, Tensor_Dim > normalize()
#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
@ HDIV
field with continuous normal traction
#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_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
PetscErrorCode DMoFEMMeshToLocalVector(DM dm, Vec l, InsertMode mode, ScatterMode scatter_mode)
set local (or ghosted) vector values on mesh for partition only
PetscErrorCode DMoFEMLoopFiniteElements(DM dm, const char fe_name[], MoFEM::FEMethod *method, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
auto createDMVector(DM dm)
Get smart vector from DM.
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.
#define MOFEM_LOG(channel, severity)
Log.
SeverityLevel
Severity levels.
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
FTensor::Index< 'J', DIM1 > J
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 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
auto snesGetKSP(SNES snes)
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data)
MoFEM::OpBrokenTopoBase OP
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
tools::promote_args< ArgumentType, int >::type lambert_w(const ArgumentType &z)
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 > lambdaVec
SmartPetscObj< Vec > kspSolVec
SmartPetscObj< Vec > duplicateGradientVec() override
SmartPetscObj< Vec > duplicateKappaVec() override
CommInterface::EntitiesPetscVector & getKappaVec() override
SmartPetscObj< TS > timeSolver
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 constexpr double pi
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 constexpr double eps
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 constexpr double root_pi
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)
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
SmartPetscObj< Vec > vec_dJdu
boost::shared_ptr< double > totalDissipation
boost::shared_ptr< std::array< MatrixDouble, 2 > > lambdaPtr
MoFEMErrorCode aSsemble(EntData &data)
MoFEMErrorCode iNtegrate(EntData &data)
MoFEMErrorCode iNtegrate(EntData &data)
MoFEMErrorCode aSsemble(EntData &data)
static double min_stiffness
MoFEMErrorCode getOptions()
MoFEMErrorCode doWork(int row_side, EntityType row_type, EntitiesFieldData::EntData &row_data)
boost::shared_ptr< double > gcPtr
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
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.
default operator for TRI element
@ 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
default operator for TET element
BoundaryEle::UserDataOperator BdyEleOp