#ifndef __HENCKY_OPS_HPP__
#define __HENCKY_OPS_HPP__
constexpr double eps = std::numeric_limits<float>::epsilon();
auto f = [](
double v) {
return 0.5 * std::log(
v); };
auto d_f = [](
double v) {
return 0.5 /
v; };
auto dd_f = [](
double v) {
return -0.5 / (
v *
v); };
struct isEq {
static inline auto check(
const double &
a,
const double &b) {
}
};
inline auto is_eq(
const double &
a,
const double &b) {
};
template <
int DIM>
inline auto get_uniq_nb(
double *ptr) {
std::array<double, DIM> tmp;
std::copy(ptr, &ptr[DIM], tmp.begin());
std::sort(tmp.begin(), tmp.end());
isEq::absMax = std::max(std::abs(tmp[0]), std::abs(tmp[DIM - 1]));
return std::distance(tmp.begin(), std::unique(tmp.begin(), tmp.end(),
is_eq));
};
template <int DIM>
std::max(std::max(std::abs(eig(0)), std::abs(eig(1))), std::abs(eig(2)));
if (
is_eq(eig(0), eig(1))) {
}
else if (
is_eq(eig(0), eig(2))) {
}
else if (
is_eq(eig(1), eig(2))) {
}
eigen_vec(
i, 0), eigen_vec(
i, 1), eigen_vec(
i, 2),
eigen_vec(
j, 0), eigen_vec(
j, 1), eigen_vec(
j, 2),
eigen_vec(
k, 0), eigen_vec(
k, 1), eigen_vec(
k, 2)};
{
eigen_vec(
i,
j) = eigen_vec_c(
i,
j);
}
};
struct CommonData :
public boost::enable_shared_from_this<CommonData> {
boost::shared_ptr<MatrixDouble>
matDPtr;
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
matLogC);
}
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &
matTangent);
}
};
template <int DIM, IntegrationType I, typename DomainEleOp>
struct OpCalculateEigenValsImpl;
template <int DIM, IntegrationType I, typename DomainEleOp>
struct OpCalculateLogCImpl;
template <int DIM, IntegrationType I, typename DomainEleOp>
struct OpCalculateLogC_dCImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpCalculateHenckyStressImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpCalculateHenckyThermalStressImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpCalculateHenckyThermalStressImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpCalculateHenckyPlasticStressImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpCalculatePiolaStressImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpHenckyTangentImpl;
template <int DIM, IntegrationType I, typename DomainEleOp, int S>
struct OpCalculateHenckyThermalStressdTImpl;
template <int DIM, typename DomainEleOp>
OpCalculateEigenValsImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data)
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
commonDataPtr->matEigVal.resize(DIM, nb_gauss_pts, false);
commonDataPtr->matEigVec.resize(DIM * DIM, nb_gauss_pts, false);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
for (int ii = 0; ii != DIM; ii++)
for (int jj = 0; jj != DIM; jj++)
eigen_vec(ii, jj) = C(ii, jj);
CHKERR computeEigenValuesSymmetric(eigen_vec, eig);
for (auto ii = 0; ii != DIM; ++ii)
eig(ii) = std::max(std::numeric_limits<double>::epsilon(), eig(ii));
if constexpr (DIM == 3) {
if (nb_uniq == 2) {
}
}
t_eig_vec(
i,
j) = eigen_vec(
i,
j);
++t_grad;
++t_eig_val;
++t_eig_vec;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename DomainEleOp>
boost::shared_ptr<CommonData> common_data)
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matLogC.resize(
size_symm, nb_gauss_pts,
false);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
++t_eig_val;
++t_eig_vec;
++t_logC;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename DomainEleOp>
OpCalculateLogC_dCImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data)
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
++t_logC_dC;
++t_eig_val;
++t_eig_vec;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename DomainEleOp, int S>
struct OpCalculateHenckyStressImpl<DIM, GAUSS,
DomainEleOp, S>
OpCalculateHenckyStressImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data)
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matHenckyStress.resize(
size_symm, nb_gauss_pts,
false);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
t_T(
i,
j) = t_D(
i,
j,
k,
l) * t_logC(
k,
l);
++t_logC;
++t_T;
++t_D;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename DomainEleOp, int S>
struct OpCalculateHenckyThermalStressImpl<DIM, GAUSS,
DomainEleOp, S>
OpCalculateHenckyThermalStressImpl(
const std::string
field_name, boost::shared_ptr<VectorDouble> temperature,
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<VectorDouble> coeff_expansion_ptr,
boost::shared_ptr<double> ref_temp_ptr)
commonDataPtr(common_data), coeffExpansionPtr(coeff_expansion_ptr),
refTempPtr(ref_temp_ptr) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matHenckyStress.resize(
size_symm, nb_gauss_pts,
false);
commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
commonDataPtr->matSecondPiolaStress.resize(
size_symm, nb_gauss_pts,
false);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
auto t_S =
getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
auto t_temp = getFTensor0FromVec(*tempPtr);
t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
}
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
#ifdef HENCKY_SMALL_STRAIN
(t_grad(
k,
l) - t_coeff_exp(
k,
l) * (t_temp - (*refTempPtr)));
#else
(t_logC(
k,
l) - t_coeff_exp(
k,
l) * (t_temp - (*refTempPtr)));
t_S(
k,
l) = t_T(
i,
j) * t_logC_dC(
i,
j,
k,
l);
t_P(
i,
l) = t_F(
i,
k) * t_S(
k,
l);
#endif
++t_grad;
++t_logC;
++t_logC_dC;
++t_P;
++t_T;
++t_S;
++t_D;
++t_temp;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
boost::shared_ptr<VectorDouble> tempPtr;
boost::shared_ptr<VectorDouble> coeffExpansionPtr;
boost::shared_ptr<double> refTempPtr;
};
template <int DIM, typename DomainEleOp, int S>
struct OpCalculateHenckyPlasticStressImpl<DIM, GAUSS,
DomainEleOp, S>
OpCalculateHenckyPlasticStressImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<MatrixDouble> mat_D_ptr,
scaleStress(
scale), matDPtr(mat_D_ptr) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
matLogCPlastic = commonDataPtr->matLogCPlastic;
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
auto t_logCPlastic = getFTensor2SymmetricFromMat<DIM>(*matLogCPlastic);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matHenckyStress.resize(
size_symm, nb_gauss_pts,
false);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
t_T(
i,
j) = t_D(
i,
j,
k,
l) * (t_logC(
k,
l) - t_logCPlastic(
k,
l));
t_T(
i,
j) /= scaleStress;
++t_logC;
++t_T;
++t_D;
++t_logCPlastic;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
boost::shared_ptr<MatrixDouble> matLogCPlastic;
const double scaleStress;
};
template <int DIM, typename DomainEleOp, int S>
struct OpCalculatePiolaStressImpl<DIM, GAUSS,
DomainEleOp, S>
OpCalculatePiolaStressImpl(
const std::string
field_name,
boost::shared_ptr<CommonData> common_data)
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
#ifdef HENCKY_SMALL_STRAIN
auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*commonDataPtr->matDPtr);
#endif
auto t_logC = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matLogC);
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
commonDataPtr->matFirstPiolaStress.resize(DIM * DIM, nb_gauss_pts, false);
commonDataPtr->matSecondPiolaStress.resize(
size_symm, nb_gauss_pts,
false);
auto t_P = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matFirstPiolaStress);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
auto t_S =
getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
#ifdef HENCKY_SMALL_STRAIN
t_P(
i,
j) = t_D(
i,
j,
k,
l) * t_grad(
k,
l);
#else
t_S(
k,
l) = t_T(
i,
j) * t_logC_dC(
i,
j,
k,
l);
t_P(
i,
l) = t_F(
i,
k) * t_S(
k,
l);
#endif
++t_grad;
++t_logC;
++t_logC_dC;
++t_P;
++t_T;
++t_S;
#ifdef HENCKY_SMALL_STRAIN
++t_D;
#endif
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename DomainEleOp, int S>
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<MatrixDouble> mat_D_ptr = nullptr)
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
if (mat_D_ptr)
matDPtr = mat_D_ptr;
else
matDPtr = commonDataPtr->matDPtr;
}
MoFEMErrorCode doWork(
int side, EntityType type,
EntData &data) {
const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
commonDataPtr->matTangent.resize(DIM * DIM * DIM * DIM, nb_gauss_pts);
auto dP_dF =
getFTensor4FromMat<DIM, DIM, DIM, DIM, 1>(commonDataPtr->matTangent);
auto t_D = getFTensor4DdgFromMat<DIM, DIM, S>(*matDPtr);
auto t_eig_val = getFTensor1FromMat<DIM>(commonDataPtr->matEigVal);
auto t_eig_vec = getFTensor2FromMat<DIM, DIM>(commonDataPtr->matEigVec);
auto t_T = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matHenckyStress);
auto t_S =
getFTensor2SymmetricFromMat<DIM>(commonDataPtr->matSecondPiolaStress);
auto t_grad = getFTensor2FromMat<DIM, DIM>(*(commonDataPtr->matGradPtr));
auto t_logC_dC = getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->matLogCdC);
for (size_t gg = 0; gg != nb_gauss_pts; ++gg) {
#ifdef HENCKY_SMALL_STRAIN
#else
t_T, nb_uniq);
P_D_P_plus_TL(
i,
j,
k,
l) =
(t_logC_dC(
i,
j, o, p) * t_D(o, p,
m,
n)) * t_logC_dC(
m,
n,
k,
l);
P_D_P_plus_TL(
i,
j,
k,
l) *= 0.5;
t_F(
i,
k) * (P_D_P_plus_TL(
k,
j, o, p) * dC_dF(o, p,
m,
n));
#endif
++dP_dF;
++t_grad;
++t_eig_val;
++t_eig_vec;
++t_logC_dC;
++t_S;
++t_T;
++t_D;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
};
template <int DIM, typename AssemblyDomainEleOp, int S>
OpCalculateHenckyThermalStressdTImpl(
const std::string row_field_name, const std::string col_field_name,
boost::shared_ptr<CommonData> elastic_common_data_ptr,
boost::shared_ptr<VectorDouble> coeff_expansion_ptr);
MoFEMErrorCode iNtegrate(EntitiesFieldData::EntData &row_data,
EntitiesFieldData::EntData &col_data);
private:
boost::shared_ptr<CommonData> elasticCommonDataPtr;
boost::shared_ptr<VectorDouble> coeffExpansionPtr;
};
template <int DIM, typename AssemblyDomainEleOp, int S>
OpCalculateHenckyThermalStressdTImpl<DIM, GAUSS, AssemblyDomainEleOp, S>::
OpCalculateHenckyThermalStressdTImpl(
const std::string row_field_name, const std::string col_field_name,
boost::shared_ptr<CommonData> elastic_common_data_ptr,
boost::shared_ptr<VectorDouble> coeff_expansion_ptr)
AssemblyDomainEleOp::OPROWCOL),
elasticCommonDataPtr(elastic_common_data_ptr),
coeffExpansionPtr(coeff_expansion_ptr) {
this->sYmm = false;
}
template <int DIM, typename AssemblyDomainEleOp, int S>
MoFEMErrorCode
OpCalculateHenckyThermalStressdTImpl<DIM, GAUSS, AssemblyDomainEleOp, S>::
iNtegrate(EntitiesFieldData::EntData &row_data,
EntitiesFieldData::EntData &col_data) {
auto &locMat = AssemblyDomainEleOp::locMat;
const auto nb_integration_pts = row_data.getN().size1();
const auto nb_row_base_functions = row_data.getN().size2();
auto t_w = this->getFTensor0IntegrationWeight();
auto t_row_diff_base = row_data.getFTensor1DiffN<DIM>();
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*elasticCommonDataPtr->matDPtr);
auto t_grad =
getFTensor2FromMat<DIM, DIM>(*(elasticCommonDataPtr->matGradPtr));
auto t_logC_dC =
getFTensor4DdgFromMat<DIM, DIM>(elasticCommonDataPtr->matLogCdC);
t_coeff_exp(d, d) = (*coeffExpansionPtr)[d];
}
t_eigen_strain(
i,
j) = (t_D(
i,
j,
k,
l) * t_coeff_exp(
k,
l));
for (auto gg = 0; gg != nb_integration_pts; ++gg) {
double alpha = this->getMeasure() * t_w;
auto rr = 0;
for (; rr != AssemblyDomainEleOp::nbRows / DIM; ++rr) {
auto t_mat = getFTensor1FromMat<DIM, 1>(locMat, rr * DIM);
auto t_col_base = col_data.getFTensor0N(gg, 0);
for (auto cc = 0; cc != AssemblyDomainEleOp::nbCols; cc++) {
#ifdef HENCKY_SMALL_STRAIN
(t_row_diff_base(
j) * t_eigen_strain(
i,
j)) * (t_col_base * alpha);
#else
t_mat(
i) -= (t_row_diff_base(
j) *
(t_F(
i, o) * ((t_D(
m,
n,
k,
l) * t_coeff_exp(
k,
l)) *
t_logC_dC(
m,
n, o,
j)))) *
(t_col_base * alpha);
#endif
++t_mat;
++t_col_base;
}
++t_row_diff_base;
}
for (; rr != nb_row_base_functions; ++rr)
++t_row_diff_base;
++t_w;
++t_grad;
++t_logC_dC;
++t_D;
}
}
template <typename DomainEleOp> struct HenckyIntegrators {
template <int DIM, IntegrationType I>
template <int DIM, IntegrationType I>
template <int DIM, IntegrationType I>
template <int DIM, IntegrationType I, int S>
OpCalculateHenckyStressImpl<DIM, I, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S>
OpCalculateHenckyThermalStressImpl<DIM, I, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S>
OpCalculateHenckyPlasticStressImpl<DIM, I, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S>
OpCalculatePiolaStressImpl<DIM, GAUSS, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S>
template <int DIM, IntegrationType I, typename AssemblyDomainEleOp, int S>
OpCalculateHenckyThermalStressdTImpl<DIM, GAUSS, AssemblyDomainEleOp, S>;
};
template <int DIM>
MoFEMErrorCode
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string block_name,
boost::shared_ptr<MatrixDouble> mat_D_Ptr, Sev sev,
PetscBool plane_strain_flag = PETSC_FALSE;
CHKERR PetscOptionsGetBool(PETSC_NULLPTR,
"",
"-plane_strain",
&plane_strain_flag, PETSC_NULLPTR);
std::vector<const CubitMeshSets *> meshset_vec_ptr,
double scale, PetscBool plane_strain_flag)
planeStrainFlag(plane_strain_flag) {
"Can not get data from block");
}
MoFEMErrorCode doWork(int side, EntityType type,
EntitiesFieldData::EntData &data) {
for (auto &b : blockData) {
if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
CHKERR getMatDPtr(matDPtr, b.bulkModulusK * scaleYoungModulus,
b.shearModulusG * scaleYoungModulus,
planeStrainFlag);
}
}
CHKERR getMatDPtr(matDPtr, bulkModulusKDefault * scaleYoungModulus,
shearModulusGDefault * scaleYoungModulus,
planeStrainFlag);
}
private:
boost::shared_ptr<MatrixDouble> matDPtr;
const double scaleYoungModulus;
const PetscBool planeStrainFlag;
double bulkModulusK;
double shearModulusG;
};
double bulkModulusKDefault;
double shearModulusGDefault;
std::vector<BlockData> blockData;
MoFEMErrorCode
std::vector<const CubitMeshSets *> meshset_vec_ptr,
Sev sev) {
for (
auto m : meshset_vec_ptr) {
std::vector<double> block_data;
CHKERR m->getAttributes(block_data);
if (block_data.size() != 2) {
"Expected that block has two attribute");
}
auto get_block_ents = [&]() {
m_field.
get_moab().get_entities_by_handle(
m->meshset, ents,
true);
return ents;
};
blockData.push_back(
}
}
MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
auto set_material_stiffness = [&]() {
: 1;
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*mat_D_ptr);
};
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
set_material_stiffness();
}
};
double nu = 0.3;
PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "none");
CHKERR PetscOptionsScalar(
"-young_modulus",
"Young modulus",
"",
E, &
E,
PETSC_NULLPTR);
CHKERR PetscOptionsScalar(
"-poisson_ratio",
"poisson ratio",
"", nu, &nu,
PETSC_NULLPTR);
PetscOptionsEnd();
pip.push_back(new OpMatBlocks(
m_field.
getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
(boost::format("%s(.*)") % block_name).str()
)),
));
}
template <int DIM, IntegrationType I, typename DomainEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string
field_name, std::string block_name, Sev sev,
double scale = 1) {
auto common_ptr = boost::make_shared<HenckyOps::CommonData>();
common_ptr->matDPtr = boost::make_shared<MatrixDouble>();
common_ptr->matGradPtr = boost::make_shared<MatrixDouble>();
common_ptr->matDPtr, sev,
scale),
"addMatBlockOps");
using H = HenckyIntegrators<DomainEleOp>;
pip.push_back(new OpCalculateVectorFieldGradient<DIM, DIM>(
pip.push_back(new typename H::template OpCalculateEigenVals<DIM, I>(
pip.push_back(
new typename H::template OpCalculateLogC<DIM, I>(
field_name, common_ptr));
pip.push_back(new typename H::template OpCalculateLogC_dC<DIM, I>(
pip.push_back(new typename H::template OpCalculateHenckyStress<DIM, I, 0>(
pip.push_back(new typename H::template OpCalculatePiolaStress<DIM, I, 0>(
return common_ptr;
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string
field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
Sev sev) {
using B =
typename FormsIntegrators<DomainEleOp>::template Assembly<
pip.push_back(
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string
field_name, std::string block_name, Sev sev,
double scale = 1) {
common_ptr, sev);
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string
field_name, boost::shared_ptr<HenckyOps::CommonData> common_ptr,
Sev sev) {
using B =
typename FormsIntegrators<DomainEleOp>::template Assembly<
using OpKPiola =
typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
using H = HenckyIntegrators<DomainEleOp>;
pip.push_back(new typename H::template OpHenckyTangent<DIM, I, 0>(
pip.push_back(
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string
field_name, std::string block_name, Sev sev,
double scale = 1) {
common_ptr, sev);
}
}
#endif
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
Kronecker Delta class symmetric.
#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()
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
auto getMat(A &&t_val, B &&t_vec, Fun< double > f)
Get the Mat object.
auto getDiffMat(A &&t_val, B &&t_vec, Fun< double > f, Fun< double > d_f, const int nb)
Get the Diff Mat object.
auto getDiffDiffMat(A &&t_val, B &&t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, C &&t_S, const int nb)
Get the Diff Diff Mat object.
auto get_uniq_nb(double *ptr)
auto sort_eigen_vals(FTensor::Tensor1< double, DIM > &eig, FTensor::Tensor2< double, DIM, DIM > &eigen_vec)
MoFEMErrorCode opFactoryDomainLhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< HenckyOps::CommonData > common_ptr, Sev sev)
MoFEMErrorCode opFactoryDomainRhs(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, boost::shared_ptr< HenckyOps::CommonData > common_ptr, Sev sev)
MoFEMErrorCode addMatBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, boost::shared_ptr< MatrixDouble > mat_D_Ptr, Sev sev, double scale=1)
auto commonDataFactory(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string field_name, std::string block_name, Sev sev, double scale=1)
auto is_eq(const double &a, const double &b)
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
double young_modulus
Young modulus.
double poisson_ratio
Poisson ratio.
constexpr auto field_name
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
FTensor::Index< 'm', 3 > m
auto getMatHenckyStress()
auto getMatFirstPiolaStress()
boost::shared_ptr< MatrixDouble > matLogCPlastic
MatrixDouble matFirstPiolaStress
boost::shared_ptr< MatrixDouble > matDPtr
MatrixDouble matSecondPiolaStress
MatrixDouble matHenckyStress
boost::shared_ptr< MatrixDouble > matGradPtr
OpCalculateHenckyThermalStressdTImpl< DIM, GAUSS, AssemblyDomainEleOp, S > OpCalculateHenckyThermalStressdT
OpCalculatePiolaStressImpl< DIM, GAUSS, DomainEleOp, S > OpCalculatePiolaStress
OpCalculateLogCImpl< DIM, I, DomainEleOp > OpCalculateLogC
OpCalculateHenckyThermalStressImpl< DIM, I, DomainEleOp, S > OpCalculateHenckyThermalStress
OpCalculateHenckyPlasticStressImpl< DIM, I, DomainEleOp, S > OpCalculateHenckyPlasticStress
OpCalculateLogC_dCImpl< DIM, I, DomainEleOp > OpCalculateLogC_dC
OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S > OpHenckyTangent
OpCalculateEigenValsImpl< DIM, I, DomainEleOp > OpCalculateEigenVals
OpCalculateHenckyStressImpl< DIM, I, DomainEleOp, S > OpCalculateHenckyStress
static auto check(const double &a, const double &b)
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
Data on single entity (This is passed as argument to DataOperator::doWork)
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
FormsIntegrators< DomainEleOp >::Assembly< A >::OpBase AssemblyDomainEleOp