8 #ifndef __HENCKY_OPS_HPP__
9 #define __HENCKY_OPS_HPP__
13 constexpr
double eps = std::numeric_limits<float>::epsilon();
15 auto f = [](
double v) {
return 0.5 * std::log(
v); };
16 auto d_f = [](
double v) {
return 0.5 /
v; };
17 auto dd_f = [](
double v) {
return -0.5 / (
v *
v); };
20 static inline auto check(
const double &
a,
const double &b) {
28 inline auto is_eq(
const double &
a,
const double &b) {
33 std::array<double, DIM> tmp;
34 std::copy(ptr, &ptr[DIM], tmp.begin());
35 std::sort(tmp.begin(), tmp.end());
36 isEq::absMax = std::max(std::abs(tmp[0]), std::abs(tmp[DIM - 1]));
37 return std::distance(tmp.begin(), std::unique(tmp.begin(), tmp.end(),
is_eq));
45 std::max(std::max(std::abs(eig(0)), std::abs(eig(1))), std::abs(eig(2)));
47 int i = 0,
j = 1,
k = 2;
49 if (
is_eq(eig(0), eig(1))) {
53 }
else if (
is_eq(eig(0), eig(2))) {
57 }
else if (
is_eq(eig(1), eig(2))) {
64 eigen_vec(
i, 0), eigen_vec(
i, 1), eigen_vec(
i, 2),
66 eigen_vec(
j, 0), eigen_vec(
j, 1), eigen_vec(
j, 2),
68 eigen_vec(
k, 0), eigen_vec(
k, 1), eigen_vec(
k, 2)};
76 eigen_vec(
i,
j) = eigen_vec_c(
i,
j);
80 struct CommonData :
public boost::enable_shared_from_this<CommonData> {
95 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
100 return boost::shared_ptr<MatrixDouble>(shared_from_this(),
105 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
matLogC);
109 return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
matTangent);
113 template <
int DIM, IntegrationType I,
typename DomainEleOp>
116 template <
int DIM, IntegrationType I,
typename DomainEleOp>
119 template <
int DIM, IntegrationType I,
typename DomainEleOp>
122 template <
int DIM, IntegrationType I,
typename DomainEleOp,
int S>
125 template <
int DIM, IntegrationType I,
typename DomainEleOp,
int S>
128 template <
int DIM, IntegrationType I,
typename DomainEleOp,
int S>
131 template <
int DIM, IntegrationType I,
typename DomainEleOp,
int S>
134 template <
int DIM,
typename DomainEleOp>
138 boost::shared_ptr<CommonData> common_data)
140 commonDataPtr(common_data) {
141 std::fill(&DomainEleOp::doEntities[MBEDGE],
142 &DomainEleOp::doEntities[MBMAXTYPE],
false);
154 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
156 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
158 commonDataPtr->matEigVal.resize(DIM, nb_gauss_pts,
false);
159 commonDataPtr->matEigVec.resize(DIM * DIM, nb_gauss_pts,
false);
160 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
161 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
163 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
173 for (
int ii = 0; ii != DIM; ii++)
174 for (
int jj = 0; jj != DIM; jj++)
175 eigen_vec(ii, jj) = C(ii, jj);
178 for (
auto ii = 0; ii != DIM; ++ii)
179 eig(ii) = std::max(std::numeric_limits<double>::epsilon(), eig(ii));
182 auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
183 if constexpr (DIM == 3) {
185 sort_eigen_vals<DIM>(eig, eigen_vec);
189 t_eig_val(
i) = eig(
i);
190 t_eig_vec(
i,
j) = eigen_vec(
i,
j);
204 template <
int DIM,
typename DomainEleOp>
208 boost::shared_ptr<CommonData> common_data)
210 commonDataPtr(common_data) {
211 std::fill(&DomainEleOp::doEntities[MBEDGE],
212 &DomainEleOp::doEntities[MBMAXTYPE],
false);
223 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
224 constexpr
auto size_symm = (DIM * (DIM + 1)) / 2;
225 commonDataPtr->matLogC.resize(
size_symm, nb_gauss_pts,
false);
227 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
228 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
230 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
232 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
236 eig(
i) = t_eig_val(
i);
237 eigen_vec(
i,
j) = t_eig_vec(
i,
j);
239 t_logC(
i,
j) = logC(
i,
j);
253 template <
int DIM,
typename DomainEleOp>
257 boost::shared_ptr<CommonData> common_data)
259 commonDataPtr(common_data) {
260 std::fill(&DomainEleOp::doEntities[MBEDGE],
261 &DomainEleOp::doEntities[MBMAXTYPE],
false);
272 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
273 constexpr
auto size_symm = (DIM * (DIM + 1)) / 2;
275 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
276 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
277 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
279 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
283 eig(
i) = t_eig_val(
i);
284 eigen_vec(
i,
j) = t_eig_vec(
i,
j);
287 auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
289 dlogC_dC(
i,
j,
k,
l) *= 2;
291 t_logC_dC(
i,
j,
k,
l) = dlogC_dC(
i,
j,
k,
l);
305 template <
int DIM,
typename DomainEleOp,
int S>
310 boost::shared_ptr<CommonData> common_data)
312 commonDataPtr(common_data) {
313 std::fill(&DomainEleOp::doEntities[MBEDGE],
314 &DomainEleOp::doEntities[MBMAXTYPE],
false);
328 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
329 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
330 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
331 constexpr
auto size_symm = (DIM * (DIM + 1)) / 2;
332 commonDataPtr->matHenckyStress.resize(
size_symm, nb_gauss_pts,
false);
333 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
335 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
336 t_T(
i,
j) = t_D(
i,
j,
k,
l) * t_logC(
k,
l);
349 template <
int DIM,
typename DomainEleOp,
int S>
354 boost::shared_ptr<CommonData> common_data,
355 boost::shared_ptr<MatrixDouble> mat_D_ptr,
356 const double scale = 1)
358 scaleStress(
scale), matDPtr(mat_D_ptr) {
359 std::fill(&DomainEleOp::doEntities[MBEDGE],
360 &DomainEleOp::doEntities[MBMAXTYPE],
false);
362 matLogCPlastic = commonDataPtr->matLogCPlastic;
376 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
377 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
378 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
379 auto t_logCPlastic = getFTensor2SymmetricFromMat<DIM>(*matLogCPlastic);
380 constexpr
auto size_symm = (DIM * (DIM + 1)) / 2;
381 commonDataPtr->matHenckyStress.resize(
size_symm, nb_gauss_pts,
false);
382 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
384 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
385 t_T(
i,
j) = t_D(
i,
j,
k,
l) * (t_logC(
k,
l) - t_logCPlastic(
k,
l));
386 t_T(
i,
j) /= scaleStress;
403 template <
int DIM,
typename DomainEleOp,
int S>
408 boost::shared_ptr<CommonData> common_data)
410 commonDataPtr(common_data) {
411 std::fill(&DomainEleOp::doEntities[MBEDGE],
412 &DomainEleOp::doEntities[MBMAXTYPE],
false);
426 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
427 #ifdef HENCKY_SMALL_STRAIN
428 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
430 auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
431 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
432 constexpr
auto size_symm = (DIM * (DIM + 1)) / 2;
433 commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts,
false);
434 commonDataPtr->matSecondPiolaStress.resize(
size_symm, nb_gauss_pts,
false);
435 auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
436 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
438 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
439 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
441 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
443 #ifdef HENCKY_SMALL_STRAIN
444 t_P(
i,
j) = t_D(
i,
j,
k,
l) * t_grad(
k,
l);
448 t_S(
k,
l) = t_T(
i,
j) * t_logC_dC(
i,
j,
k,
l);
449 t_P(
i,
l) = t_F(
i,
k) * t_S(
k,
l);
458 #ifdef HENCKY_SMALL_STRAIN
470 template <
int DIM,
typename DomainEleOp,
int S>
473 boost::shared_ptr<CommonData> common_data,
474 boost::shared_ptr<MatrixDouble> mat_D_ptr =
nullptr)
476 commonDataPtr(common_data) {
477 std::fill(&DomainEleOp::doEntities[MBEDGE],
478 &DomainEleOp::doEntities[MBMAXTYPE],
false);
482 matDPtr = commonDataPtr->matDPtr;
499 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
500 commonDataPtr->matTangent.resize(DIM * DIM * DIM * DIM, nb_gauss_pts);
502 getFTensor4FromMat<DIM, DIM, DIM, DIM, 1>(commonDataPtr->matTangent);
504 auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
505 auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
506 auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
507 auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
509 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
510 auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
511 auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
513 for (
size_t gg = 0; gg != nb_gauss_pts; ++gg) {
515 #ifdef HENCKY_SMALL_STRAIN
516 dP_dF(
i,
j,
k,
l) = t_D(
i,
j,
k,
l);
525 eig(
i) = t_eig_val(
i);
526 eigen_vec(
i,
j) = t_eig_vec(
i,
j);
530 auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
540 P_D_P_plus_TL(
i,
j,
k,
l) =
542 (t_logC_dC(
i,
j, o, p) * t_D(o, p,
m,
n)) * t_logC_dC(
m,
n,
k,
l);
543 P_D_P_plus_TL(
i,
j,
k,
l) *= 0.5;
546 t_F(
i,
k) * (P_D_P_plus_TL(
k,
j, o, p) * dC_dF(o, p,
m,
n));
570 template <
int DIM, IntegrationType I>
573 template <
int DIM, IntegrationType I>
576 template <
int DIM, IntegrationType I>
579 template <
int DIM, IntegrationType I,
int S>
583 template <
int DIM, IntegrationType I,
int S = 0>
587 template <
int DIM, IntegrationType I,
int S>
591 template <
int DIM, IntegrationType I,
int S>
598 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
599 std::string block_name,
600 boost::shared_ptr<MatrixDouble> mat_D_Ptr,
Sev sev,
607 std::vector<const CubitMeshSets *> meshset_vec_ptr,
613 "Can not get data from block");
620 for (
auto &b : blockData) {
622 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
623 CHKERR getMatDPtr(matDPtr, b.bulkModulusK * scaleYoungModulus,
624 b.shearModulusG * scaleYoungModulus);
629 CHKERR getMatDPtr(matDPtr, bulkModulusKDefault * scaleYoungModulus,
630 shearModulusGDefault * scaleYoungModulus);
635 boost::shared_ptr<MatrixDouble> matDPtr;
636 const double scaleYoungModulus;
640 double shearModulusG;
644 double bulkModulusKDefault;
645 double shearModulusGDefault;
646 std::vector<BlockData> blockData;
650 std::vector<const CubitMeshSets *> meshset_vec_ptr,
654 for (
auto m : meshset_vec_ptr) {
656 std::vector<double> block_data;
657 CHKERR m->getAttributes(block_data);
658 if (block_data.size() != 2) {
660 "Expected that block has two attribute");
662 auto get_block_ents = [&]() {
665 m_field.
get_moab().get_entities_by_handle(
m->meshset, ents,
true);
684 MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
688 auto set_material_stiffness = [&]() {
694 double A = (DIM == 2)
698 auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*mat_D_ptr);
705 constexpr
auto size_symm = (DIM * (DIM + 1)) / 2;
707 set_material_stiffness();
715 CHKERR PetscOptionsBegin(PETSC_COMM_WORLD,
"",
"",
"none");
716 CHKERR PetscOptionsScalar(
"-young_modulus",
"Young modulus",
"",
E, &
E,
718 CHKERR PetscOptionsScalar(
"-poisson_ratio",
"poisson ratio",
"", nu, &nu,
720 ierr = PetscOptionsEnd();
724 pip.push_back(
new OpMatBlocks(
728 m_field.
getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
730 (boost::format(
"%s(.*)") % block_name).str()
740 template <
int DIM, IntegrationType I,
typename DomainEleOp>
743 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
746 auto common_ptr = boost::make_shared<HenckyOps::CommonData>();
747 common_ptr->matDPtr = boost::make_shared<MatrixDouble>();
748 common_ptr->matGradPtr = boost::make_shared<MatrixDouble>();
751 common_ptr->matDPtr, sev,
scale),
758 pip.push_back(
new typename H::template OpCalculateEigenVals<DIM, I>(
761 new typename H::template OpCalculateLogC<DIM, I>(
field_name, common_ptr));
762 pip.push_back(
new typename H::template OpCalculateLogC_dC<DIM, I>(
765 pip.push_back(
new typename H::template OpCalculateHenckyStress<DIM, I, 0>(
767 pip.push_back(
new typename H::template OpCalculatePiolaStress<DIM, I, 0>(
773 template <
int DIM, AssemblyType A, IntegrationType I,
typename DomainEleOp>
776 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
777 std::string
field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
781 using B =
typename FormsIntegrators<DomainEleOp>::template Assembly<
782 A>::template LinearForm<I>;
791 template <
int DIM, AssemblyType A, IntegrationType I,
typename DomainEleOp>
794 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
798 auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
806 template <
int DIM, AssemblyType A, IntegrationType I,
typename DomainEleOp>
809 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
810 std::string
field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
814 using B =
typename FormsIntegrators<DomainEleOp>::template Assembly<
816 using OpKPiola =
typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
821 new typename H::template OpHenckyTangent<DIM, I, 0>(
field_name, common_ptr));
828 template <
int DIM, AssemblyType A, IntegrationType I,
typename DomainEleOp>
831 boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
835 auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
844 #endif // __HENCKY_OPS_HPP__