19#include <boost/math/constants/constants.hpp>
20#include <boost/math/special_functions/lambert_w.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);
58 return T(0.5) * boost::math::lambert_w(z, 0);
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);
134 t_gap_np1(
i) = alpha * (t_normal(
i) + (t_tangential(
i) / beta));
138 template <
typename T>
142 double alpha,
double beta) {
147 for (
auto jj = 0; jj != 3; ++jj) {
148 t_cpx_traction(
i) = t_traction(
i);
149 t_cpx_traction(jj) +=
eps * 1
i;
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;
159 template <
typename T>
166 T t_eff = std::sqrt(
t(
i) *
t(
i));
196 PetscOptionsBegin(PETSC_COMM_WORLD,
"interface_",
"",
"none");
198 CHKERR PetscOptionsScalar(
"-gc",
"Griffith energy release rate",
"",
200 CHKERR PetscOptionsScalar(
"-min_stiffness",
"Minimal interface stiffness",
202 CHKERR PetscOptionsScalar(
"-strength",
"Strength of interface",
"",
204 CHKERR PetscOptionsScalar(
"-min_kappa",
205 "Minimal kappa to avoid singularity",
"",
207 CHKERR PetscOptionsScalar(
"-kappa0",
"Characteristic length kappa0",
"",
209 CHKERR PetscOptionsScalar(
"-beta",
"Cohesive tangential coupling",
"",
215 <<
"Interface Griffith energy release rate Gc -interface_gc = "
218 <<
"Interface min stiffness -interface_min_stiffness = "
221 <<
"Interface strength -interface_strength = " <<
strength;
223 <<
"Interface minimal kappa -interface_min_kappa = " <<
min_kappa;
225 <<
"Interface characteristic length kappa0 -interface_kappa0 = "
229 double kappa_strength =
233 <<
"Adjusted interface min kappa to capture strength = " <<
min_kappa;
248 const int def_size = 1;
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)
258static Tag get_tag(moab::Interface &moab, std::string tag_name,
int size) {
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)
269 return get_tag(moab,
"KAPPA", 1);
273 return get_tag(moab,
"DELTA_KAPPA", 1);
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) {}
304 auto do_work_rhs = [
this](
int row_side, EntityType row_type,
308 auto base = row_data.
getBase();
311 "row base not set properly");
337 switch (OP::opType) {
343 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData());
350 (std::string(
"wrong op type ") +
351 OpBaseDerivativesBase::OpTypeNames[OP::opType])
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)
383 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
389 int nb_integration_pts = OP::getGaussPts().size2();
392 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
394 auto t_face_normal = getFTensor1NormalsAtGaussPts();
396 int nb_base_functions = data.
getN().size2() / 3;
398 auto t_w = OP::getFTensor0IntegrationWeight();
408 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
412 t_normalized_normal(
J) = t_normal(
J);
415 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
426 t_gap_double(
i) = -t_gap(
i);
431 auto t_nf = getFTensor1FromPtr<3>(&*
OP::locF.begin());
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);
439 for (; bb != nb_base_functions; ++bb)
443 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
457 boost::shared_ptr<std::array<MatrixDouble, 2>>
lambdaPtr;
473 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
480 int nb_integration_pts = OP::getGaussPts().size2();
483 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
485 auto t_face_normal = getFTensor1NormalsAtGaussPts();
487 int nb_base_functions = data.
getN().size2() / 3;
489 auto t_w = OP::getFTensor0IntegrationWeight();
499 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
503 t_normalized_normal(
J) = t_normal(
J);
506 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
518 t_dgap_double(
i,
j) = -t_dgap(
i,
j);
519 return t_dgap_double;
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));
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);
536 for (; rr != nb_base_functions; ++rr)
540 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
573 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
579 int nb_integration_pts = OP::getGaussPts().size2();
582 auto t_lambda = getFTensor2FromMat<3, 3>(
lambdaPtr->at(get_sense_index()));
583 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
585 auto t_face_normal = getFTensor1NormalsAtGaussPts();
586 auto t_w = OP::getFTensor0IntegrationWeight();
597 double face_dissipation = 0.0;
598 double face_grad_dissipation = 0.0;
599 auto face_handle = OP::getFEEntityHandle();
603 &face_grad_dissipation);
605 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
609 t_normalized_normal(
J) = t_normal(
J);
612 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
621 delta_kappa, teff,
kappa, gc);
624 delta_kappa, teff,
kappa, gc);
625 return boost::make_tuple(
m, m_grad);
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;
636 t_traction, t_normalized_normal, diff_alpha,
640 t_gap_double(
i) = -t_gap(
i);
644 auto [
J, dJ] = dJ_dkappa(t_delta_kappa, t_traction, t_normalized_normal,
646 face_dissipation += t_w *
J * t_normal.
l2();
647 face_grad_dissipation += t_w * dJ * t_normal.
l2();
649 auto t_dr = dr_kappa(t_delta_kappa, t_traction, t_normalized_normal,
652 t_l(
i) = t_lambda(
i,
K) * t_normal(
K);
653 face_grad_dissipation -= t_w * t_l(
i) * t_dr(
i);
661 &face_grad_dissipation);
682 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
687 int nb_integration_pts = OP::getGaussPts().size2();
690 int nb_base_functions = data.
getN().size2() / 3;
692 auto t_w = OP::getFTensor0IntegrationWeight();
695 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
697 auto t_face_normal = getFTensor1NormalsAtGaussPts();
707 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
711 t_normalized_normal(
J) = t_normal(
J);
714 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
722 t_traction_double(
i) = t_traction(
i);
725 t_traction_double, delta_kappa,
kappa, t_normalized_normal, gc,
728 t_dJ_double(
i) = t_dM(
i);
734 auto t_nf = getFTensor1FromPtr<3>(&*
OP::locF.begin());
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);
742 for (; bb != nb_base_functions; ++bb)
746 assemble(dJ_dtraction(t_traction, t_normalized_normal, t_delta_kappa,
772 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
774 boost::shared_ptr<VectorDouble> tag_data_ptr,
775 boost::shared_ptr<Range> ents_ptr =
nullptr)
782 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
794 auto get_data = [&](
auto &
v) {
797 auto &moab = getMoab();
798 auto fe_ent = getFEEntityHandle();
802 rval = moab.tag_get_by_ptr(
tagHandle, &fe_ent, 1, (
const void **)&data,
806 rval == MB_SUCCESS && size > 0 && size !=
v.size()
810 "Inconsistent size of tag data");
814 tag_size[0] =
v.size();
815 void const *tag_data[] = {&
v[0]};
819 (
const void **)&data, &size);
820 std::copy(data, data + size,
v.begin());
845 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
850 using SideEleOp = EleOnSide::UserDataOperator;
852 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
855 auto domain_side_flux = [&](
auto &pip) {
857 auto broken_data_ptr =
858 boost::make_shared<std::vector<BrokenBaseSideData>>();
861 auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
866 return broken_data_ptr;
869 auto get_lambda = [&](
auto &pip) {
870 boost::shared_ptr<std::array<MatrixDouble, 2>> array_lambda_ptr;
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),
880 [lambda_mat_ptr, array_lambda_ptr](
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;
888 array_lambda_ptr->at(get_sense_index()) = *lambda_mat_ptr;
893 return array_lambda_ptr;
896 return std::make_pair(domain_side_flux(pip), get_lambda(pip));
902 boost::shared_ptr<Range> interface_range_ptr,
905 auto &m_field = ep.
mField;
911 using SideEleOp = EleOnSide::UserDataOperator;
912 using BdyEleOp = BoundaryEle::UserDataOperator;
914 auto face_side = [&]() {
917 auto op_loop_skeleton_side =
919 interface_range_ptr, Sev::noisy);
920 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = [](int, int, int) {
923 op_loop_skeleton_side->getSideFEPtr()->setRuleHook =
924 set_integration_at_front_face;
926 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
930 return op_loop_skeleton_side;
933 auto op_loop_skeleton_side = face_side();
940 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
941 boost::make_shared<CGGUserPolynomialBase>();
943 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
946 op_reset_broken_data->doWorkRhsHook =
947 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
950 broken_data_flux_ptr->resize(0);
953 op_loop_skeleton_side->getOpPtrVector().push_back(op_reset_broken_data);
954 op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
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>();
960 op_loop_skeleton_side->getOpPtrVector().push_back(
962 op_loop_skeleton_side->getOpPtrVector().push_back(
965 op_loop_skeleton_side->getOpPtrVector().push_back(
new OpGetTag(
969 return std::make_tuple(op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr,
970 kappa_ptr, kappa_delta_ptr);
976 boost::shared_ptr<Range> interface_range_ptr,
977 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
980 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
983 interface_range_ptr);
985 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
986 op_loop_skeleton_side->getOpPtrVector().push_back(
989 op_loop_skeleton_side->getOpPtrVector().push_back(
new OpCohesiveRhs(
990 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
992 pip.push_back(op_loop_skeleton_side);
1000 boost::shared_ptr<Range> interface_range_ptr,
1001 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1004 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
1007 interface_range_ptr);
1008 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1009 op_loop_skeleton_side->getOpPtrVector().push_back(
1013 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1015 pip.push_back(op_loop_skeleton_side);
1022 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1024 auto &m_field = ep.
mField;
1030 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1031 boost::make_shared<CGGUserPolynomialBase>();
1033 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
1037 op_reset_broken_data->doWorkRhsHook =
1038 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
1041 broken_data_flux_ptr->resize(0);
1044 pip.push_back(op_reset_broken_data);
1045 pip.push_back(op_loop_domain_side);
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);
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>();
1056 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1059 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1064 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr, lambda_ptr,
1065 tag_dissipation, tag_grad_dissipation));
1067 return std::make_pair(tag_dissipation, tag_grad_dissipation);
1072 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1073 auto &m_field = ep.
mField;
1079 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1080 boost::make_shared<CGGUserPolynomialBase>();
1081 auto [broken_data_flux_ptr, lambda_ptr] =
1086 op_reset_broken_data->doWorkRhsHook =
1087 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
1090 broken_data_flux_ptr->resize(0);
1093 pip.push_back(op_reset_broken_data);
1094 pip.push_back(op_loop_domain_side);
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>();
1101 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1104 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1108 pip.push_back(
new OpCohesive_dJ_dP(broken_data_flux_ptr, gc_ptr, kappa_ptr,
1109 kappa_delta_ptr, lambda_ptr,
Tag(),
Tag(),
1125 using BdyEleOp = BoundaryEle::UserDataOperator;
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(
1135 auto interface_face = [&](
FEMethod *fe_method_ptr) {
1136 auto ent = fe_method_ptr->getFEEntityHandle();
1147 fe_ptr->exeTestHook = interface_face;
1152 auto face_fe = get_face_ele();
1153 auto [tag_dissipation, tag_grad_dissipation] =
1156 constexpr double zero = 0.0;
1166 ep.
mField.
get_moab(), grad_dissipation_vec, tag_grad_dissipation);
1179 using BdyEleOp = BoundaryEle::UserDataOperator;
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(
1189 auto interface_face = [&](
FEMethod *fe_method_ptr) {
1190 auto ent = fe_method_ptr->getFEEntityHandle();
1201 fe_ptr->exeTestHook = interface_face;
1206 auto face_fe = get_face_ele();
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);
1217 CHKERR VecNorm(dJ_dx, NORM_2, &dJ_dx_norm2);
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);
1224 CHKERR KSPSolveTranspose(ksp, dJ_dx, lambda_vec);
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);
1231 <<
"evaluateCohesiveLambdaImpl: Norm of lambda vector: " << lambda_norm2;
1241 CHKERR TSSetFromOptions(ts);
1242 CHKERR TSSetStepNumber(ts, 0);
1244 CHKERR TSSetSolution(ts, x);
1248 <<
"evaluatePrimalProblemCohesiveImpl: Time step dt: " <<
dt;
1249 CHKERR TSSolve(ts, PETSC_NULLPTR);
1254 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
1255 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
1258 CHKERR VecNorm(x, NORM_2, &x_norm2);
1260 <<
"evaluatePrimalProblemCohesiveImpl: Norm of displacement vector: "
1277 ep->mField.get_comm(), ep->mField.get_moab(),
1282 [&ep](
Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1286 [&ep](
Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1324 return boost::make_shared<CohesiveTAOCtxImpl>(
1325 ep_ptr, set_integration_at_front_face, time_solver);
1329 PetscReal *f, Vec
g,
1334 auto &ep = *(cohesive_ctx->ep_ptr);
1337 CHKERR VecView(delta_kappa, PETSC_VIEWER_STDOUT_WORLD);
1340 CHKERR VecCopy(delta_kappa, cohesive_ctx->kappaVec.second);
1341 CHKERR VecGhostUpdateBegin(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1343 CHKERR VecGhostUpdateEnd(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1345 CHKERR CommInterface::setTagFromVector(
1351 cohesive_ctx->kspSolVec);
1353 auto ksp = snesGetKSP(tsGetSNES(cohesive_ctx->timeSolver));
1355 ksp, cohesive_ctx->lambdaVec);
1358 ep, cohesive_ctx->setIntegrationAtFrontFace, cohesive_ctx->lambdaVec,
1359 cohesive_ctx->dissipationVec, cohesive_ctx->gradDissipationVec);
1361 CHKERR VecSum(cohesive_ctx->dissipationVec.second, f);
1362 CHKERR VecCopy(cohesive_ctx->gradDissipationVec.second,
g);
1365 <<
"Cohesive objective function (negative total dissipation): " << *f;
Eshelbian plasticity interface.
Add cohesive elements to Eshelbian plasticity module.
#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
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
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
constexpr double t
plate stiffness
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)
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
SmartPetscObj< Vec > vec_dJdu
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.
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)
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
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.
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
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
default operator for TET element
BoundaryEle::UserDataOperator BdyEleOp