v0.14.0
HenckyOps.hpp
/**
* \file HenckyOps.hpp
* \example HenckyOps.hpp
*
* @copyright Copyright (c) 2023
*/
#ifndef __HENCKY_OPS_HPP__
#define __HENCKY_OPS_HPP__
namespace HenckyOps {
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) {
return std::abs(a - b) / absMax < eps;
}
static double absMax;
};
double isEq::absMax = 1;
inline auto is_eq(const double &a, const double &b) {
return isEq::check(a, 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)));
int i = 0, j = 1, k = 2;
if (is_eq(eig(0), eig(1))) {
i = 0;
j = 2;
k = 1;
} else if (is_eq(eig(0), eig(2))) {
i = 0;
j = 1;
k = 2;
} else if (is_eq(eig(1), eig(2))) {
i = 1;
j = 0;
k = 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)};
FTensor::Tensor1<double, 3> eig_c{eig(i), eig(j), eig(k)};
{
eig(i) = eig_c(i);
eigen_vec(i, j) = eigen_vec_c(i, j);
}
};
struct CommonData : public boost::enable_shared_from_this<CommonData> {
boost::shared_ptr<MatrixDouble> matGradPtr;
boost::shared_ptr<MatrixDouble> matDPtr;
boost::shared_ptr<MatrixDouble> matLogCPlastic;
inline auto getMatFirstPiolaStress() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
inline auto getMatHenckyStress() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(),
}
inline auto getMatLogC() {
return boost::shared_ptr<MatrixDouble>(shared_from_this(), &matLogC);
}
inline auto getMatTangent() {
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 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, typename DomainEleOp>
struct OpCalculateEigenValsImpl<DIM, GAUSS, DomainEleOp> : public DomainEleOp {
OpCalculateEigenValsImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
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) {
F(i, j) = t_grad(i, j) + t_kd(i, j);
C(i, j) = F(k, i) ^ F(k, j);
for (int ii = 0; ii != DIM; ii++)
for (int jj = 0; jj != DIM; jj++)
eigen_vec(ii, jj) = C(ii, jj);
for (auto ii = 0; ii != DIM; ++ii)
eig(ii) = std::max(std::numeric_limits<double>::epsilon(), eig(ii));
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
if constexpr (DIM == 3) {
if (nb_uniq == 2) {
sort_eigen_vals<DIM>(eig, eigen_vec);
}
}
t_eig_val(i) = eig(i);
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>
struct OpCalculateLogCImpl<DIM, GAUSS, DomainEleOp> : public DomainEleOp {
OpCalculateLogCImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
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) {
eig(i) = t_eig_val(i);
eigen_vec(i, j) = t_eig_vec(i, j);
auto logC = EigenMatrix::getMat(eig, eigen_vec, f);
t_logC(i, j) = logC(i, j);
++t_eig_val;
++t_eig_vec;
++t_logC;
}
}
private:
boost::shared_ptr<CommonData> commonDataPtr;
};
template <int DIM, typename DomainEleOp>
struct OpCalculateLogC_dCImpl<DIM, GAUSS, DomainEleOp> : public DomainEleOp {
OpCalculateLogC_dCImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
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->matLogCdC.resize(size_symm * size_symm, nb_gauss_pts, false);
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) {
eig(i) = t_eig_val(i);
eigen_vec(i, j) = t_eig_vec(i, j);
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
auto dlogC_dC = EigenMatrix::getDiffMat(eig, eigen_vec, f, d_f, nb_uniq);
dlogC_dC(i, j, k, l) *= 2;
t_logC_dC(i, j, k, l) = dlogC_dC(i, j, k, l);
++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>
: public DomainEleOp {
OpCalculateHenckyStressImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
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 OpCalculateHenckyPlasticStressImpl<DIM, GAUSS, DomainEleOp, S>
: public DomainEleOp {
OpCalculateHenckyPlasticStressImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<MatrixDouble> mat_D_ptr,
const double scale = 1)
: DomainEleOp(field_name, DomainEleOp::OPROW), commonDataPtr(common_data),
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) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
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>
: public DomainEleOp {
OpCalculatePiolaStressImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data)
: DomainEleOp(field_name, DomainEleOp::OPROW),
commonDataPtr(common_data) {
std::fill(&DomainEleOp::doEntities[MBEDGE],
&DomainEleOp::doEntities[MBMAXTYPE], false);
}
MoFEMErrorCode doWork(int side, EntityType type, EntData &data) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
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_F(i, j) = t_grad(i, j) + t_kd(i, j);
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>
struct OpHenckyTangentImpl<DIM, GAUSS, DomainEleOp, S> : public DomainEleOp {
OpHenckyTangentImpl(const std::string field_name,
boost::shared_ptr<CommonData> common_data,
boost::shared_ptr<MatrixDouble> mat_D_ptr = nullptr)
: DomainEleOp(field_name, DomainEleOp::OPROW),
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) {
constexpr auto t_kd = FTensor::Kronecker_Delta<int>();
// const size_t nb_gauss_pts = matGradPtr->size2();
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
dP_dF(i, j, k, l) = t_D(i, j, k, l);
#else
t_F(i, j) = t_grad(i, j) + t_kd(i, j);
eig(i) = t_eig_val(i);
eigen_vec(i, j) = t_eig_vec(i, j);
T(i, j) = t_T(i, j);
// rare case when two eigen values are equal
auto nb_uniq = get_uniq_nb<DIM>(&eig(0));
dC_dF(i, j, k, l) = (t_kd(i, l) * t_F(k, j)) + (t_kd(j, l) * t_F(k, i));
auto TL =
EigenMatrix::getDiffDiffMat(eig, eigen_vec, f, d_f, dd_f, T, nb_uniq);
TL(i, j, k, l) *= 4;
P_D_P_plus_TL(i, j, k, l) =
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;
dP_dF(i, j, m, n) = t_kd(i, m) * (t_kd(k, n) * t_S(k, j));
dP_dF(i, j, m, n) +=
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 <typename DomainEleOp> struct HenckyIntegrators {
template <int DIM, IntegrationType I>
using OpCalculateEigenVals = OpCalculateEigenValsImpl<DIM, I, DomainEleOp>;
template <int DIM, IntegrationType I>
using OpCalculateLogC = OpCalculateLogCImpl<DIM, I, DomainEleOp>;
template <int DIM, IntegrationType I>
using OpCalculateLogC_dC = OpCalculateLogC_dCImpl<DIM, I, DomainEleOp>;
template <int DIM, IntegrationType I, int S>
OpCalculateHenckyStressImpl<DIM, I, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S = 0>
OpCalculateHenckyPlasticStressImpl<DIM, I, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S>
OpCalculatePiolaStressImpl<DIM, GAUSS, DomainEleOp, S>;
template <int DIM, IntegrationType I, int S>
using OpHenckyTangent = OpHenckyTangentImpl<DIM, GAUSS, DomainEleOp, S>;
};
template <int DIM>
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string block_name,
boost::shared_ptr<MatrixDouble> mat_D_Ptr, Sev sev,
double scale = 1) {
struct OpMatBlocks : public DomainEleOp {
OpMatBlocks(boost::shared_ptr<MatrixDouble> m, double bulk_modulus_K,
double shear_modulus_G, MoFEM::Interface &m_field, Sev sev,
std::vector<const CubitMeshSets *> meshset_vec_ptr,
double scale)
: DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), matDPtr(m),
bulkModulusKDefault(bulk_modulus_K),
shearModulusGDefault(shear_modulus_G), scaleYoungModulus(scale) {
CHK_THROW_MESSAGE(extractBlockData(m_field, meshset_vec_ptr, sev),
"Can not get data from block");
}
MoFEMErrorCode doWork(int side, EntityType type,
for (auto &b : blockData) {
if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
CHKERR getMatDPtr(matDPtr, b.bulkModulusK * scaleYoungModulus,
b.shearModulusG * scaleYoungModulus);
}
}
CHKERR getMatDPtr(matDPtr, bulkModulusKDefault * scaleYoungModulus,
shearModulusGDefault * scaleYoungModulus);
}
private:
boost::shared_ptr<MatrixDouble> matDPtr;
const double scaleYoungModulus;
struct BlockData {
double bulkModulusK;
double shearModulusG;
Range blockEnts;
};
double bulkModulusKDefault;
double shearModulusGDefault;
std::vector<BlockData> blockData;
extractBlockData(MoFEM::Interface &m_field,
std::vector<const CubitMeshSets *> meshset_vec_ptr,
Sev sev) {
for (auto m : meshset_vec_ptr) {
MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock") << *m;
std::vector<double> block_data;
CHKERR m->getAttributes(block_data);
if (block_data.size() != 2) {
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
"Expected that block has two attribute");
}
auto get_block_ents = [&]() {
Range ents;
m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
return ents;
};
double young_modulus = block_data[0];
double poisson_ratio = block_data[1];
double bulk_modulus_K = young_modulus / (3 * (1 - 2 * poisson_ratio));
MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
<< "E = " << young_modulus << " nu = " << poisson_ratio;
blockData.push_back(
{bulk_modulus_K, shear_modulus_G, get_block_ents()});
}
}
MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
double bulk_modulus_K, double shear_modulus_G) {
//! [Calculate elasticity tensor]
auto set_material_stiffness = [&]() {
double A = (DIM == 2)
: 1;
auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*mat_D_ptr);
t_D(i, j, k, l) =
2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) * t_kd(i, j) *
t_kd(k, l);
};
//! [Calculate elasticity tensor]
constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
mat_D_ptr->resize(size_symm * size_symm, 1);
set_material_stiffness();
}
};
double E = 1.0;
double nu = 0.3;
CHKERR PetscOptionsBegin(PETSC_COMM_WORLD, "", "", "none");
CHKERR PetscOptionsScalar("-young_modulus", "Young modulus", "", E, &E,
PETSC_NULL);
CHKERR PetscOptionsScalar("-poisson_ratio", "poisson ratio", "", nu, &nu,
PETSC_NULL);
ierr = PetscOptionsEnd();
double bulk_modulus_K = E / (3 * (1 - 2 * nu));
double shear_modulus_G = E / (2 * (1 + nu));
pip.push_back(new OpMatBlocks(
mat_D_Ptr, bulk_modulus_K, shear_modulus_G, m_field, sev,
// Get blockset using regular expression
m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
(boost::format("%s(.*)") % block_name).str()
)),
));
}
template <int DIM, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
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>();
CHK_THROW_MESSAGE(addMatBlockOps<DIM>(m_field, pip, block_name,
common_ptr->matDPtr, sev, scale),
"addMatBlockOps");
using H = HenckyIntegrators<DomainEleOp>;
field_name, common_ptr->matGradPtr));
pip.push_back(new typename H::template OpCalculateEigenVals<DIM, I>(
field_name, common_ptr));
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>(
field_name, common_ptr));
// Assumes constant D matrix per entity
pip.push_back(new typename H::template OpCalculateHenckyStress<DIM, I, 0>(
field_name, common_ptr));
pip.push_back(new typename H::template OpCalculatePiolaStress<DIM, I, 0>(
field_name, common_ptr));
return common_ptr;
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
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<
A>::template LinearForm<I>;
typename B::template OpGradTimesTensor<1, DIM, DIM>;
pip.push_back(
new OpInternalForcePiola("U", common_ptr->getMatFirstPiolaStress()));
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string field_name, std::string block_name, Sev sev, double scale = 1) {
auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
m_field, pip, field_name, block_name, sev, scale);
CHKERR opFactoryDomainRhs<DIM, A, I, DomainEleOp>(m_field, pip, field_name,
common_ptr, sev);
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
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<
A>::template BiLinearForm<I>;
using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
using H = HenckyIntegrators<DomainEleOp>;
// Assumes constant D matrix per entity
pip.push_back(
new typename H::template OpHenckyTangent<DIM, I, 0>(field_name, common_ptr));
pip.push_back(
new OpKPiola(field_name, field_name, common_ptr->getMatTangent()));
}
template <int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp>
MoFEM::Interface &m_field,
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string field_name, std::string block_name, Sev sev, double scale = 1) {
auto common_ptr = commonDataFactory<DIM, I, DomainEleOp>(
m_field, pip, field_name, block_name, sev, scale);
CHKERR opFactoryDomainLhs<DIM, A, I, DomainEleOp>(m_field, pip, field_name,
common_ptr, sev);
}
} // namespace HenckyOps
#endif // __HENCKY_OPS_HPP__
NOSPACE
@ NOSPACE
Definition: definitions.h:83
EigenMatrix::getMat
FTensor::Tensor2_symmetric< double, 3 > getMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f)
Get the Mat object.
Definition: MatrixFunction.cpp:53
MoFEMFunctionReturnHot
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:460
MoFEM::UnknownInterface::getInterface
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.
Definition: UnknownInterface.hpp:93
MoFEM::EntitiesFieldData::EntData
Data on single entity (This is passed as argument to DataOperator::doWork)
Definition: EntitiesFieldData.hpp:128
HenckyOps::CommonData::matLogCPlastic
boost::shared_ptr< MatrixDouble > matLogCPlastic
Definition: HenckyOps.hpp:83
HenckyOps::HenckyIntegrators::OpCalculateLogC
OpCalculateLogCImpl< DIM, I, DomainEleOp > OpCalculateLogC
Definition: HenckyOps.hpp:574
HenckyOps::CommonData::getMatFirstPiolaStress
auto getMatFirstPiolaStress()
Definition: HenckyOps.hpp:94
HenckyOps::HenckyIntegrators::OpCalculateHenckyPlasticStress
OpCalculateHenckyPlasticStressImpl< DIM, I, DomainEleOp, S > OpCalculateHenckyPlasticStress
Definition: HenckyOps.hpp:585
FTensor::Tensor1
Definition: Tensor1_value.hpp:8
MOFEM_LOG_CHANNEL
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
Definition: LogManager.hpp:284
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:121
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
E
HenckyOps::isEq::absMax
static double absMax
Definition: HenckyOps.hpp:23
HenckyOps::CommonData::matEigVec
MatrixDouble matEigVec
Definition: HenckyOps.hpp:86
MoFEM::Types::MatrixDouble
UBlasMatrix< double > MatrixDouble
Definition: Types.hpp:77
HenckyOps::d_f
auto d_f
Definition: HenckyOps.hpp:16
HenckyOps::CommonData::getMatLogC
auto getMatLogC()
Definition: HenckyOps.hpp:104
HenckyOps::opFactoryDomainLhs
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)
Definition: HenckyOps.hpp:807
FTensor::Kronecker_Delta
Kronecker Delta class.
Definition: Kronecker_Delta.hpp:15
OpKPiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradTensorGrad< 1, SPACE_DIM, SPACE_DIM, 1 > OpKPiola
[Only used for dynamics]
Definition: seepage.cpp:64
A
constexpr AssemblyType A
Definition: operators_tests.cpp:30
bulk_modulus_K
double bulk_modulus_K
Definition: dynamic_first_order_con_law.cpp:96
HenckyOps::CommonData::matEigVal
MatrixDouble matEigVal
Definition: HenckyOps.hpp:85
HenckyOps
Definition: HenckyOps.hpp:11
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
scale
double scale
Definition: plastic.cpp:119
HenckyOps::CommonData::matLogC
MatrixDouble matLogC
Definition: HenckyOps.hpp:87
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
HenckyOps::isEq::check
static auto check(const double &a, const double &b)
Definition: HenckyOps.hpp:20
HenckyOps::is_eq
auto is_eq(const double &a, const double &b)
Definition: HenckyOps.hpp:28
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
DomainEleOp
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
Definition: child_and_parent.cpp:36
OpInternalForcePiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition: seepage.cpp:66
BlockData
Definition: ElasticityMixedFormulation.hpp:12
a
constexpr double a
Definition: approx_sphere.cpp:30
HenckyOps::HenckyIntegrators::OpCalculateEigenVals
OpCalculateEigenValsImpl< DIM, I, DomainEleOp > OpCalculateEigenVals
Definition: HenckyOps.hpp:571
convert.type
type
Definition: convert.py:64
HenckyOps::eps
constexpr double eps
Definition: HenckyOps.hpp:13
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:122
HenckyOps::CommonData::getMatHenckyStress
auto getMatHenckyStress()
Definition: HenckyOps.hpp:99
MoFEM::GAUSS
@ GAUSS
Definition: FormsIntegrators.hpp:136
HenckyOps::CommonData::matFirstPiolaStress
MatrixDouble matFirstPiolaStress
Definition: HenckyOps.hpp:89
HenckyOps::CommonData::matSecondPiolaStress
MatrixDouble matSecondPiolaStress
Definition: HenckyOps.hpp:90
HenckyOps::dd_f
auto dd_f
Definition: HenckyOps.hpp:17
FTensor::Tensor4
Definition: Tensor4_value.hpp:18
size_symm
constexpr auto size_symm
Definition: plastic.cpp:42
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
BiLinearForm
OpGradTimesTensor
FormsIntegrators< DomainEleOp >::Assembly< A >::LinearForm< I >::OpGradTimesTensor< 1, FIELD_DIM, SPACE_DIM > OpGradTimesTensor
Definition: operators_tests.cpp:48
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
field_name
constexpr auto field_name
Definition: poisson_2d_homogeneous.cpp:13
FTensor::Index< 'i', DIM >
HenckyOps::HenckyIntegrators::OpHenckyTangent
OpHenckyTangentImpl< DIM, GAUSS, DomainEleOp, S > OpHenckyTangent
Definition: HenckyOps.hpp:592
convert.n
n
Definition: convert.py:82
HenckyOps::HenckyIntegrators::OpCalculateHenckyStress
OpCalculateHenckyStressImpl< DIM, I, DomainEleOp, S > OpCalculateHenckyStress
Definition: HenckyOps.hpp:581
HenckyOps::CommonData::matLogCdC
MatrixDouble matLogCdC
Definition: HenckyOps.hpp:88
v
const double v
phase velocity of light in medium (cm/ns)
Definition: initial_diffusion.cpp:40
Range
DomainEleOp
HenckyOps::sort_eigen_vals
auto sort_eigen_vals(FTensor::Tensor1< double, DIM > &eig, FTensor::Tensor2< double, DIM, DIM > &eigen_vec)
Definition: HenckyOps.hpp:41
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
HenckyOps::CommonData::matGradPtr
boost::shared_ptr< MatrixDouble > matGradPtr
Definition: HenckyOps.hpp:81
HenckyOps::CommonData::matDPtr
boost::shared_ptr< MatrixDouble > matDPtr
Definition: HenckyOps.hpp:82
shear_modulus_G
double shear_modulus_G
Definition: dynamic_first_order_con_law.cpp:97
HenckyOps::CommonData::getMatTangent
auto getMatTangent()
Definition: HenckyOps.hpp:108
HenckyOps::opFactoryDomainRhs
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)
Definition: HenckyOps.hpp:774
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
HenckyOps::HenckyIntegrators::OpCalculateLogC_dC
OpCalculateLogC_dCImpl< DIM, I, DomainEleOp > OpCalculateLogC_dC
Definition: HenckyOps.hpp:577
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
FTensor::Ddg
Definition: Ddg_value.hpp:7
CommonData
Definition: continuity_check_on_skeleton_with_simple_2d_for_h1.cpp:22
MoFEM::Exceptions::ierr
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
Definition: Exceptions.hpp:76
EigenMatrix::getDiffDiffMat
FTensor::Ddg< double, 3, 3 > getDiffDiffMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f, Fun< double > d_f, Fun< double > dd_f, FTensor::Tensor2< double, 3, 3 > &t_S, const int nb)
Definition: MatrixFunction.cpp:78
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
HenckyOps::commonDataFactory
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)
Definition: HenckyOps.hpp:741
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
HenckyOps::addMatBlockOps
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)
Definition: HenckyOps.hpp:597
HenckyOps::CommonData::matTangent
MatrixDouble matTangent
Definition: HenckyOps.hpp:92
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
EigenMatrix::getDiffMat
FTensor::Ddg< double, 3, 3 > getDiffMat(Val< double, 3 > &t_val, Vec< double, 3 > &t_vec, Fun< double > f, Fun< double > d_f, const int nb)
Get the Diff Mat object.
Definition: MatrixFunction.cpp:64
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
HenckyOps::CommonData::matHenckyStress
MatrixDouble matHenckyStress
Definition: HenckyOps.hpp:91
HenckyOps::HenckyIntegrators::OpCalculatePiolaStress
OpCalculatePiolaStressImpl< DIM, GAUSS, DomainEleOp, S > OpCalculatePiolaStress
Definition: HenckyOps.hpp:589
HenckyOps::get_uniq_nb
auto get_uniq_nb(double *ptr)
Definition: HenckyOps.hpp:32
H
double H
Hardening.
Definition: plastic.cpp:124
OpCalculateVectorFieldGradient
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
l
FTensor::Index< 'l', 3 > l
Definition: matrix_function.cpp:21
MoFEM::computeEigenValuesSymmetric
MoFEMErrorCode computeEigenValuesSymmetric(const MatrixDouble &mat, VectorDouble &eig, MatrixDouble &eigen_vec)
compute eigenvalues of a symmetric matrix using lapack dsyev
Definition: Templates.hpp:1452
F
@ F
Definition: free_surface.cpp:394