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)));
54 T dnom =
kappa * min_stiffness + tau;
55 return (tau -
kappa * diff_tau) / (dnom * dnom);
58 template <
typename T>
static auto invTau(
const T &tau,
double Gf) {
59 const T z = (2.0 * Gf * Gf) / (
pi * tau * tau);
65 double min_stiffness) {
67 return 0.5 * t_eff * t_eff * diff_alpha;
72 const T &
kappa,
double gc,
73 double min_stiffness) {
75 return Y * delta_kappa;
80 const T &
kappa,
double gc,
81 double min_stiffness) {
84 return -M * delta_kappa;
90 const T &
kappa,
double gc,
91 double min_stiffness) {
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 * 1
i;
97 std::complex<T> cpx_M =
99 return cpx_M.imag() /
eps;
102 template <
typename T>
106 double beta,
double min_stiffness) {
109 std::complex<T> cpx_delta = delta_kappa;
110 std::complex<T> cpx_kappa =
kappa;
112 for (
auto jj = 0; jj != 3; ++jj) {
113 t_cpx_traction(
i) = t_traction(
i);
114 t_cpx_traction(jj) +=
eps * 1
i;
115 std::complex<T> cpx_teff =
118 cpx_delta, cpx_teff, cpx_kappa, gc, min_stiffness);
119 t_diff_traction(jj) = cpx_M.imag() /
eps;
121 return t_diff_traction;
124 template <
typename T>
127 double alpha,
double beta,
128 bool sign_sensitive =
true) {
133 t_P(
i,
j) = t_n_normalize(
i) * t_n_normalize(
j);
137 t_normal(
i) = t_P(
i,
j) * t_traction(
j);
139 t_tangential(
i) = t_Q(
i,
j) * t_traction(
j);
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;
147 if (std::real(s) > std::numeric_limits<double>::epsilon()) {
150 t_n_normalize(
i) * teff / 2. +
152 (1.0 / 4.0) * (teff / s) * t_normal(
i) +
154 (1.0 / beta) * (teff / s) * t_tangential(
i)
158 return std::make_pair(teff, t_gap);
160 T teff = std::sqrt((t_normal(
i) * t_normal(
i)) +
161 (1. / beta) * t_tangential(
i) * t_tangential(
i));
165 t_normal(
i) + (1.0 / beta) * t_tangential(
i)
168 return std::make_pair(teff, t_gap);
172 template <
typename T>
176 double alpha,
double beta,
177 bool sign_sensitive =
true) {
182 for (
auto jj = 0; jj != 3; ++jj) {
183 t_cpx_traction(
i) = t_traction(
i);
184 t_cpx_traction(jj) +=
eps * 1
i;
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;
195 template <
typename T>
200 auto [teff, t_gap] =
calculateGap(t_traction, t_n_normalize, 1.0, beta);
229 PetscOptionsBegin(PETSC_COMM_WORLD,
"interface_",
"",
"none");
231 CHKERR PetscOptionsScalar(
"-gc",
"Griffith energy release rate",
"",
233 CHKERR PetscOptionsScalar(
"-min_stiffness",
"Minimal interface stiffness",
235 CHKERR PetscOptionsScalar(
"-strength",
"Strength of interface",
"",
237 CHKERR PetscOptionsScalar(
"-min_kappa",
238 "Minimal kappa to avoid singularity",
"",
240 CHKERR PetscOptionsScalar(
"-kappa0",
"Characteristic length kappa0",
"",
242 CHKERR PetscOptionsScalar(
"-beta",
"Cohesive tangential coupling",
"",
248 <<
"Interface Griffith energy release rate Gc -interface_gc = "
251 <<
"Interface min stiffness -interface_min_stiffness = "
254 <<
"Interface strength -interface_strength = " <<
strength;
256 <<
"Interface minimal kappa -interface_min_kappa = " <<
min_kappa;
258 <<
"Interface characteristic length kappa0 -interface_kappa0 = "
261 <<
"Interface tangential coupling -interface_beta = " <<
beta;
264 double kappa_strength =
268 <<
"Adjusted interface min kappa to capture strength = " <<
min_kappa;
280static Tag get_tag(moab::Interface &moab, std::string tag_name,
int size) {
281 std::vector<double> dummy(size, 0.);
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)
292 return get_tag(moab,
"KAPPA", 1);
296 return get_tag(moab,
"DELTA_KAPPA", 1);
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) {}
317 if (OP::entsPtr->find(this->getFEEntityHandle()) == OP::entsPtr->end())
322 if (!brokenBaseSideData) {
327 auto do_work_rhs = [
this](
int row_side, EntityType row_type,
331 auto base = row_data.
getBase();
334 "row base not set properly");
343 OP::nbIntegrationPts = OP::getGaussPts().size2();
345 OP::nbRowBaseFunctions = OP::getNbOfBaseFunctions(row_data);
347 OP::locF.resize(OP::nbRows,
false);
349 OP::locMat.resize(OP::nbRows, OP::nbRows,
false);
353 CHKERR this->iNtegrate(row_data);
355 CHKERR this->aSsemble(row_data);
360 switch (OP::opType) {
362 for (
auto &bd : *brokenBaseSideData) {
364 boost::shared_ptr<MatrixDouble>(brokenBaseSideData, &bd.getFlux());
366 CHKERR do_work_rhs(bd.getSide(), bd.getType(), bd.getData());
373 (std::string(
"wrong op type ") +
374 OpBaseDerivativesBase::OpTypeNames[OP::opType])
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)
406 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
412 int nb_integration_pts = OP::getGaussPts().size2();
415 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
417 auto t_face_normal = getFTensor1NormalsAtGaussPts();
419 int nb_base_functions = data.
getN().size2() / 3;
421 auto t_w = OP::getFTensor0IntegrationWeight();
431 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
435 t_normalized_normal(
J) = t_normal(
J);
438 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
450 t_gap_double(
i) = -t_gap(
i);
455 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
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);
463 for (; bb != nb_base_functions; ++bb)
467 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
481 boost::shared_ptr<std::array<MatrixDouble, 2>>
lambdaPtr;
497 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
504 int nb_integration_pts = OP::getGaussPts().size2();
507 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
509 auto t_face_normal = getFTensor1NormalsAtGaussPts();
511 int nb_base_functions = data.
getN().size2() / 3;
513 auto t_w = OP::getFTensor0IntegrationWeight();
523 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
527 t_normalized_normal(
J) = t_normal(
J);
530 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
543 t_dgap_double(
i,
j) = -t_dgap(
i,
j);
544 return t_dgap_double;
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));
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);
561 for (; rr != nb_base_functions; ++rr)
565 assemble(fracture(t_traction, t_normalized_normal, t_kappa, t_delta_kappa,
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());
582 CHKERR this->matSetValuesHook(
this, data, data, this->locMat);
598 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
604 int nb_integration_pts = OP::getGaussPts().size2();
607 auto t_lambda = getFTensor2FromMat<3, 3>(
lambdaPtr->at(get_sense_index()));
608 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
610 auto t_face_normal = getFTensor1NormalsAtGaussPts();
611 auto t_w = OP::getFTensor0IntegrationWeight();
622 double face_dissipation = 0.0;
623 double face_grad_dissipation = 0.0;
624 auto face_handle = OP::getFEEntityHandle();
628 &face_grad_dissipation);
630 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
634 t_normalized_normal(
J) = t_normal(
J);
637 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
650 return boost::make_tuple(
m, m_grad);
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;
664 t_gap_double(
i) = -t_gap(
i);
668 auto [
J, dJ] = dJ_dkappa(t_delta_kappa, t_traction, t_normalized_normal,
670 face_dissipation += t_w *
J * t_normal.
l2();
671 face_grad_dissipation += t_w * dJ * t_normal.
l2();
673 auto t_dr = dr_kappa(t_delta_kappa, t_traction, t_normalized_normal,
676 t_l(
i) = t_lambda(
i,
K) * t_normal(
K);
677 face_grad_dissipation -= t_w * t_l(
i) * t_dr(
i);
685 &face_grad_dissipation);
710 int nb_integration_pts = OP::getGaussPts().size2();
713 int nb_base_functions = data.
getN().size2() / 3;
715 auto t_w = OP::getFTensor0IntegrationWeight();
718 auto t_kappa = getFTensor0FromVec<0>(*
kappaPtr);
720 auto t_face_normal = getFTensor1NormalsAtGaussPts();
730 for (
auto gg = 0; gg != nb_integration_pts; ++gg) {
734 t_normalized_normal(
J) = t_normal(
J);
737 t_traction(
i) = t_P(
i,
J) * t_normalized_normal(
J);
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);
745 t_traction_double(
i) = t_traction(
i);
748 t_traction_double, delta_kappa,
kappa, t_normalized_normal, gc,
751 t_dJ_double(
i) = t_dM(
i);
757 auto t_nf = getFTensor1FromPtr<3>(&*OP::locF.begin());
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);
765 for (; bb != nb_base_functions; ++bb)
769 assemble(dJ_dtraction(t_traction, t_normalized_normal, t_delta_kappa,
781 return VecSetValues<AssemblyTypeSelector<A>>(
vec_dJdu, data, OP::locF,
796 boost::shared_ptr<std::vector<BrokenBaseSideData>> broken_base_side_data,
798 boost::shared_ptr<VectorDouble> tag_data_ptr,
799 boost::shared_ptr<Range> ents_ptr =
nullptr)
806 auto get_sense_index = [
this]() {
return (
faceSense == 1) ? 0 : 1; };
818 auto get_data = [&](
auto &
v) {
821 auto &moab = getMoab();
822 auto fe_ent = getFEEntityHandle();
826 rval = moab.tag_get_by_ptr(
tagHandle, &fe_ent, 1, (
const void **)&data,
830 rval == MB_SUCCESS && size > 0 && size !=
v.size()
834 "Inconsistent size of tag data");
838 tag_size[0] =
v.size();
839 void const *tag_data[] = {&
v[0]};
843 (
const void **)&data, &size);
844 std::copy(data, data + size,
v.begin());
869 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
874 using SideEleOp = EleOnSide::UserDataOperator;
876 EshelbianPlasticity::AddHOOps<SPACE_DIM, SPACE_DIM, SPACE_DIM>::add(
879 auto domain_side_flux = [&](
auto &pip) {
881 auto broken_data_ptr =
882 boost::make_shared<std::vector<BrokenBaseSideData>>();
885 auto flux_mat_ptr = boost::make_shared<MatrixDouble>();
890 return broken_data_ptr;
893 auto get_lambda = [&](
auto &pip) {
894 boost::shared_ptr<std::array<MatrixDouble, 2>> array_lambda_ptr;
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),
904 [lambda_mat_ptr, array_lambda_ptr](
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;
912 array_lambda_ptr->at(get_sense_index()) = *lambda_mat_ptr;
917 return array_lambda_ptr;
920 return std::make_pair(domain_side_flux(pip), get_lambda(pip));
926 boost::shared_ptr<Range> interface_range_ptr,
929 auto &m_field = ep.
mField;
935 using SideEleOp = EleOnSide::UserDataOperator;
936 using BdyEleOp = BoundaryEle::UserDataOperator;
938 auto face_side = [&]() {
941 auto op_loop_skeleton_side =
943 interface_range_ptr, Sev::noisy);
944 op_loop_skeleton_side->getSideFEPtr()->getRuleHook = [](int, int, int) {
947 op_loop_skeleton_side->getSideFEPtr()->setRuleHook =
948 set_integration_at_front_face;
950 EshelbianPlasticity::AddHOOps<SPACE_DIM - 1, SPACE_DIM, SPACE_DIM>::add(
954 return op_loop_skeleton_side;
957 auto op_loop_skeleton_side = face_side();
964 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
965 boost::make_shared<CGGUserPolynomialBase>(
nullptr,
true);
967 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
970 op_reset_broken_data->doWorkRhsHook =
971 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
974 broken_data_flux_ptr->resize(0);
977 op_loop_skeleton_side->getOpPtrVector().push_back(op_reset_broken_data);
978 op_loop_skeleton_side->getOpPtrVector().push_back(op_loop_domain_side);
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>();
984 op_loop_skeleton_side->getOpPtrVector().push_back(
986 op_loop_skeleton_side->getOpPtrVector().push_back(
989 op_loop_skeleton_side->getOpPtrVector().push_back(
new OpGetTag(
993 return std::make_tuple(op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr,
994 kappa_ptr, kappa_delta_ptr);
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);
1009 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1010 op_loop_skeleton_side->getOpPtrVector().push_back(
1013 op_loop_skeleton_side->getOpPtrVector().push_back(
new OpCohesiveRhs(
1014 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1016 pip.push_back(op_loop_skeleton_side);
1024 boost::shared_ptr<Range> interface_range_ptr,
1025 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1028 auto [op_loop_skeleton_side, broken_data_flux_ptr, gc_ptr, kappa_ptr,
1031 interface_range_ptr);
1032 auto u_gamma_ptr = boost::make_shared<MatrixDouble>();
1033 op_loop_skeleton_side->getOpPtrVector().push_back(
1037 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr));
1039 pip.push_back(op_loop_skeleton_side);
1046 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
1048 auto &m_field = ep.
mField;
1054 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1055 boost::make_shared<CGGUserPolynomialBase>(
nullptr,
true);
1057 ep, op_loop_domain_side->getOpPtrVector(), lambda_vec);
1061 op_reset_broken_data->doWorkRhsHook =
1062 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
1065 broken_data_flux_ptr->resize(0);
1068 pip.push_back(op_reset_broken_data);
1069 pip.push_back(op_loop_domain_side);
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);
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>();
1080 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1083 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1088 broken_data_flux_ptr, gc_ptr, kappa_ptr, kappa_delta_ptr, lambda_ptr,
1089 tag_dissipation, tag_grad_dissipation));
1091 return std::make_pair(tag_dissipation, tag_grad_dissipation);
1096 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip) {
1097 auto &m_field = ep.
mField;
1103 op_loop_domain_side->getSideFEPtr()->getUserPolynomialBase() =
1104 boost::make_shared<CGGUserPolynomialBase>(
nullptr,
true);
1105 auto [broken_data_flux_ptr, lambda_ptr] =
1110 op_reset_broken_data->doWorkRhsHook =
1111 [broken_data_flux_ptr](
DataOperator *base_op_ptr,
int side,
1114 broken_data_flux_ptr->resize(0);
1117 pip.push_back(op_reset_broken_data);
1118 pip.push_back(op_loop_domain_side);
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>();
1125 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1128 pip.push_back(
new OpGetTag(broken_data_flux_ptr,
1132 pip.push_back(
new OpCohesive_dJ_dP(broken_data_flux_ptr, gc_ptr, kappa_ptr,
1133 kappa_delta_ptr, lambda_ptr,
Tag(),
Tag(),
1149 using BdyEleOp = BoundaryEle::UserDataOperator;
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(
1159 auto interface_face = [&](
FEMethod *fe_method_ptr) {
1160 auto ent = fe_method_ptr->getFEEntityHandle();
1171 fe_ptr->exeTestHook = interface_face;
1176 auto face_fe = get_face_ele();
1177 auto [tag_dissipation, tag_grad_dissipation] =
1180 constexpr double zero = 0.0;
1189 ep.
mField.
get_moab(), grad_dissipation_vec, tag_grad_dissipation);
1202 using BdyEleOp = BoundaryEle::UserDataOperator;
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(
1212 auto interface_face = [&](
FEMethod *fe_method_ptr) {
1213 auto ent = fe_method_ptr->getFEEntityHandle();
1224 fe_ptr->exeTestHook = interface_face;
1229 auto face_fe = get_face_ele();
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);
1240 CHKERR VecNorm(dJ_dx, NORM_2, &dJ_dx_norm2);
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);
1247 CHKERR KSPSolveTranspose(ksp, dJ_dx, lambda_vec);
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);
1254 <<
"evaluateCohesiveLambdaImpl: Norm of lambda vector: " << lambda_norm2;
1264 CHKERR TSSetFromOptions(ts);
1265 CHKERR TSSetStepNumber(ts, 0);
1267 CHKERR TSSetSolution(ts, x);
1271 <<
"evaluatePrimalProblemCohesiveImpl: Time step dt: " <<
dt;
1272 CHKERR TSSolve(ts, PETSC_NULLPTR);
1277 CHKERR VecGhostUpdateBegin(x, INSERT_VALUES, SCATTER_FORWARD);
1278 CHKERR VecGhostUpdateEnd(x, INSERT_VALUES, SCATTER_FORWARD);
1281 CHKERR VecNorm(x, NORM_2, &x_norm2);
1283 <<
"evaluatePrimalProblemCohesiveImpl: Norm of displacement vector: "
1300 ep->mField.get_comm(), ep->mField.get_moab(),
1305 [&ep](
Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1309 [&ep](
Range r) { return intersect(*(ep->interfaceFaces), r); }, 1,
1347 return boost::make_shared<CohesiveTAOCtxImpl>(
1348 ep_ptr, set_integration_at_front_face, time_solver);
1352 PetscReal *f, Vec
g,
1357 auto &ep = *(cohesive_ctx->ep_ptr);
1360 CHKERR VecView(delta_kappa, PETSC_VIEWER_STDOUT_WORLD);
1363 CHKERR VecCopy(delta_kappa, cohesive_ctx->kappaVec.second);
1364 CHKERR VecGhostUpdateBegin(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1366 CHKERR VecGhostUpdateEnd(cohesive_ctx->kappaVec.second, INSERT_VALUES,
1374 cohesive_ctx->kspSolVec);
1378 ksp, cohesive_ctx->lambdaVec);
1381 ep, cohesive_ctx->setIntegrationAtFrontFace, cohesive_ctx->lambdaVec,
1382 cohesive_ctx->dissipationVec, cohesive_ctx->gradDissipationVec);
1384 CHKERR VecSum(cohesive_ctx->dissipationVec.second, f);
1385 CHKERR VecCopy(cohesive_ctx->gradDissipationVec.second,
g);
1388 <<
"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, RowColData rc=RowColData::COL)
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, RowColData rc=RowColData::COL)
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 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)
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
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 > 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, double min_stiffness)
static T getAlpha(const T &kappa, double gc, double min_stiffness)
static constexpr double pi
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 constexpr double eps
static T getTau(const T &k, double gc)
static constexpr double root_pi
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 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