v0.14.0
ADOLCPlasticityMaterialModels.hpp
/** \file ADOLCPlasticityMaterialModels.hpp
* \ingroup adoc_plasticity
* \brief Matetial models for plasticity
*
* \example ADOLCPlasticityMaterialModels.hpp
*
*/
#ifndef __ADOLC_MATERIAL_MODELS_HPP
#define __ADOLC_MATERIAL_MODELS_HPP
#ifndef WITH_ADOL_C
#error "MoFEM need to be compiled with ADOL-C"
#endif
namespace ADOLCPlasticity {
/** \brief J2 plasticity (Kinematic Isotropic (Linear) Hardening)
* \ingroup adoc_plasticity
*/
/**
* @brief J2 plasticity (Kinematic Isotropic (Linear) Hardening)
*
* @tparam DIM model dimension (2 - is plane stress)
*/
template <int DIM> struct J2Plasticity;
/**
* @brief J2 (Von Misses) plasticity
* @ingroup adoc_plasticity
*
* Free energy:
*
* \f[
* \psi(\varepsilon^e, \alpha, \beta) =
* \frac{1}{2} \lambda \textrm{tr}[\varepsilon^e]^2
* + \mu \varepsilon^e : \varepsilon^e
* +
* \sigma_y\alpha
* +
* \frac{1}{2} \phi H \alpha^2
* +
* \frac{1}{2} (1-\phi)K \beta^2
* \f]
*
* Isotropic hardening variable \f$\alpha\f$ is in first index of
* internal variables vector. Kinematic hardening variable is in
* the remaining indices of internal variables vector.
*
* Yield function:
*
* Deviator of actual stress:
* \f[
* \eta=\textrm{dev}[\sigma]-\overline{\beta}
* \f]
* where \f$\overline{\beta}\f$ is back stress.
*
* Square of the deviator norm
* \f[
* f = \frac{3}{2} \eta:\eta
* \f]
*
* Yield function:
* \f[
* y = \sqrt{f} - \overline{\alpha}
* \f]
* where \f$f\f$ is defined in \ref eva
*
* This is associated potential model, so flow potential is equal to yield
*
*/
template <> struct J2Plasticity<3> : public ClosestPointProjection {
J2Plasticity() : ClosestPointProjection() {
nbInternalVariables = 7;
activeVariablesW.resize(12 + nbInternalVariables, false);
}
adouble tR;
double mu; //< Lamé parameter
double lambda; //< Lamé parameters
double H; //< Isotropic hardening
double K; //< Kinematic hardening
double phi; //< Combined isotropic/kinematic hardening
double sigmaY; //< Yield stress
enum ModelParams { LAMBDA, MU, ISOH, KINK, PHI, SIGMA_Y, LAST_PARAM };
using ParamsArray = std::array<double, LAST_PARAM>;
boost::shared_ptr<ParamsArray> defaultParamsArrayPtr = nullptr;
boost::shared_ptr<ParamsArray> paramsArrayPtr = nullptr;
boost::shared_ptr<ParamsArray> oldParamsArrayPtr = nullptr;
MoFEMErrorCode getDefaultMaterialParameters() {
double young_modulus = 1;
double poisson_ratio = 0.25;
sigmaY = 1;
H = 0.1;
K = 0;
phi = 1;
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-young_modulus", &young_modulus,
PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-poisson_ratio", &poisson_ratio,
0);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-sigma_y", &sigmaY, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-H", &H, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-K", &K, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-phi", &phi, PETSC_NULL);
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "Young modulus: " << young_modulus;
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "Poisson ratio: " << poisson_ratio;
defaultParamsArrayPtr = boost::make_shared<ParamsArray>();
(*defaultParamsArrayPtr)[LAMBDA] = lambda;
(*defaultParamsArrayPtr)[MU] = mu;
(*defaultParamsArrayPtr)[H] = H;
(*defaultParamsArrayPtr)[K] = K;
(*defaultParamsArrayPtr)[PHI] = phi;
(*defaultParamsArrayPtr)[SIGMA_Y] = sigmaY;
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "LAMBDA: " << (*defaultParamsArrayPtr)[LAMBDA];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "MU: " << (*defaultParamsArrayPtr)[MU];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "H: " << (*defaultParamsArrayPtr)[ISOH];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "K: " << (*defaultParamsArrayPtr)[KINK];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "PHI: " << (*defaultParamsArrayPtr)[PHI];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "SIGMA_Y: " << (*defaultParamsArrayPtr)[SIGMA_Y];
paramsArrayPtr = boost::make_shared<ParamsArray>();
std::copy(defaultParamsArrayPtr->begin(), defaultParamsArrayPtr->end(),
paramsArrayPtr->begin());
oldParamsArrayPtr = boost::make_shared<ParamsArray>();
std::fill(oldParamsArrayPtr->begin(), oldParamsArrayPtr->end(), 0);
};
/**
* @brief Get material parameters from block
*
* @param m_field
* @param pip
* @param block_name
* @param sev
* @return MoFEMErrorCode
*/
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string block_name, Sev sev) {
struct OpMatBlocks : public ForcesAndSourcesCore::UserDataOperator {
OpMatBlocks(boost::shared_ptr<ParamsArray> params_array_ptr,
boost::shared_ptr<ParamsArray> def_params_array_otr,
MoFEM::Interface &m_field, Sev sev,
std::vector<const CubitMeshSets *> meshset_vec_ptr)
: OP(NOSPACE, OP::OPSPACE), paramsArrayPtr(params_array_ptr),
defParamsArray(def_params_array_otr) {
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()) {
std::copy(b.paramsArray.begin(), b.paramsArray.end(),
paramsArrayPtr->begin());
}
}
std::copy(defParamsArray->begin(), defParamsArray->end(),
paramsArrayPtr->begin());
}
private:
boost::shared_ptr<ParamsArray> paramsArrayPtr;
boost::shared_ptr<ParamsArray> defParamsArray;
struct BlockData {
ParamsArray paramsArray;
Range blockEnts;
};
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("ADOLCPlasticityWold", sev, "JP2 parameters") << *m;
std::vector<double> block_data;
CHKERR m->getAttributes(block_data);
if (block_data.size() != 6) {
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;
};
blockData.push_back({*defParamsArray, get_block_ents()});
std::copy(block_data.begin(), block_data.end(),
blockData.back().paramsArray.begin());
const auto young_modulus = blockData.back().paramsArray[LAMBDA];
const auto poisson_ratio = blockData.back().paramsArray[MU];
// Is assumed that uset provide young_modulus and poisson_ratio in
// fist two argiumens of the block
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "Young modulus: " << (blockData.back().paramsArray)[LAMBDA];
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "Poisson ratio: " << (blockData.back().paramsArray)[MU];
blockData.back().paramsArray[LAMBDA] =
blockData.back().paramsArray[MU] = MU(young_modulus, poisson_ratio);
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "LAMBDA: " << (blockData.back().paramsArray)[LAMBDA];
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "MU: " << (blockData.back().paramsArray)[MU];
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "H: " << (blockData.back().paramsArray)[ISOH];
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "K: " << (blockData.back().paramsArray)[KINK];
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "PHI: " << (blockData.back().paramsArray)[PHI];
MOFEM_LOG("ADOLCPlasticityWold", sev)
<< "SIGMA_Y: " << (blockData.back().paramsArray)[SIGMA_Y];
}
}
};
auto cubit_meshsets_vec_ptr =
m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
std::regex((boost::format("%s(.*)") % block_name).str()));
pip.push_back(new OpMatBlocks(paramsArrayPtr, defaultParamsArrayPtr,
m_field, sev, cubit_meshsets_vec_ptr));
}
/**
* @brief Set material parameters at integration point
*
* @param tag
* @param recalculate_elastic_tangent
* @return MoFEMErrorCode
*/
MoFEMErrorCode setParams(short tag, bool &recalculate_elastic_tangent) {
switch (tag) {
case TypesTags::W: {
set_param_vec(tapesTags[tag], paramsArrayPtr->size(),
paramsArrayPtr->data());
for (auto p = 0; p != LAST_PARAM; p++) {
if ((*paramsArrayPtr)[p] != (*oldParamsArrayPtr)[p]) {
recalculate_elastic_tangent = true;
std::copy(paramsArrayPtr->begin(), paramsArrayPtr->end(),
oldParamsArrayPtr->begin());
break;
}
}
} break;
case TypesTags::Y:
break;
default:
SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
"Unknown tag for setParams");
}
}
//! [Free energy J2]
MoFEMErrorCode freeHemholtzEnergy() {
// ADOL-C register parameters (those can varied for integration points)
auto p_lambda = mkparam(lambda); // 0
auto p_mu = mkparam(mu); // 1
auto p_H = mkparam(H); // 2
auto p_K = mkparam(K); // 3
auto p_phi = mkparam(phi); // 4
auto p_sigma_y = mkparam(sigmaY); // 5
// Operator converting Voight and tensor notation
auto t_voight_strain_op = voight_to_strain_op();
// Variable in voight notation
auto t_voight_total_strain = getFTensor1FromPtr<6>(&a_sTrain[0]);
auto t_voight_plastic_strain = getFTensor1FromPtr<6>(&a_plasticStrain[0]);
auto t_voight_internal_tensor =
getFTensor1FromPtr<6>(&a_internalVariables[1]);
// Convert variables to tensor notation
tTotalStrain(i, j) = t_voight_strain_op(i, j, Z) * t_voight_total_strain(Z);
tPlasticStrain(i, j) =
t_voight_strain_op(i, j, Z) * t_voight_plastic_strain(Z);
tInternalTensor(i, j) =
t_voight_strain_op(i, j, Z) * t_voight_internal_tensor(Z);
// Evaluate elastic strain
tElasticStrain(i, j) = tTotalStrain(i, j) - tPlasticStrain(i, j);
tR = tElasticStrain(i, i);
// Strain energy
a_w = 0.5 * p_lambda * tR * tR +
p_mu * (tElasticStrain(i, j) * tElasticStrain(i, j));
// Isotropic strain energy
a_w +=
a_internalVariables[0] * p_sigma_y +
0.5 * p_phi * p_H * (a_internalVariables[0] * a_internalVariables[0]);
// Kinematic strain energy
a_w += (0.5 * (1 - p_phi) * p_K) *
(tInternalTensor(i, j) * tInternalTensor(i, j));
}
//! [Free energy J2]
//! [Yield function J2]
MoFEMErrorCode evalF() {
auto t_voight_stress_op = voight_to_stress_op();
auto t_vioght_stress = getFTensor1FromPtr<6>(&a_sTress[0]);
auto t_back_stress = getFTensor1FromPtr<6>(&a_internalFluxes[1]);
tStress(i, j) = t_voight_stress_op(i, j, Z) * t_vioght_stress(Z);
tBackStress(i, j) = t_voight_stress_op(i, j, Z) * t_back_stress(Z);
constexpr double third = boost::math::constants::third<double>();
tR = tStress(i, i) / 3.;
tDeviator(i, j) = tStress(i, j) - tR * t_kd(i, j) - tBackStress(i, j);
// Fix the constant to mach uniaxial test
constexpr double c = 3. / 2.;
f = c * tDeviator(i, j) * tDeviator(i, j);
}
MoFEMErrorCode yieldFunction() {
CHKERR evalF();
a_y = sqrt(f) - a_internalFluxes[0];
}
//! [Yield function J2]
MoFEMErrorCode flowPotential() {
CHKERR evalF();
a_h = sqrt(f) - a_internalFluxes[0];
}
};
//! [J2 2D]
template <> struct J2Plasticity<2> : public J2Plasticity<3> {
MoFEMErrorCode freeHemholtzEnergy() {
auto p_lambda = mkparam(lambda); // 0
auto p_mu = mkparam(mu); // 1
auto p_H = mkparam(H); // 2
auto p_K = mkparam(K); // 3
auto p_phi = mkparam(phi); // 4
auto p_sigma_y = mkparam(sigmaY); // 5
auto t_voight_strain_op = voight_to_strain_op();
auto t_voight_total_strain = getFTensor1FromPtr<6>(&a_sTrain[0]);
auto t_voight_plastic_strain = getFTensor1FromPtr<6>(&a_plasticStrain[0]);
auto t_voight_internal_tensor =
getFTensor1FromPtr<6>(&a_internalVariables[1]);
tTotalStrain(i, j) = t_voight_strain_op(i, j, Z) * t_voight_total_strain(Z);
tPlasticStrain(i, j) =
t_voight_strain_op(i, j, Z) * t_voight_plastic_strain(Z);
tInternalTensor(i, j) =
t_voight_strain_op(i, j, Z) * t_voight_internal_tensor(Z);
tElasticStrain(i, j) = tTotalStrain(i, j) - tPlasticStrain(i, j);
//! [Plane stress]
// Calculate st rain at ezz enforcing plane strain case
tElasticStrain(2, 2) =
(-tElasticStrain(0, 0) * p_lambda - tElasticStrain(1, 1) * p_lambda) /
(p_lambda + 2 * p_mu);
//! [Plane stress]
tR = tElasticStrain(i, i);
a_w = 0.5 * p_lambda * tR * tR +
p_mu * (tElasticStrain(i, j) * tElasticStrain(i, j));
// plastic part
// isotropic
a_w +=
a_internalVariables[0] * p_sigma_y +
0.5 * p_phi * p_H * (a_internalVariables[0] * a_internalVariables[0]);
// kinematic
a_w += (0.5 * (1 - p_phi) * p_K) *
(tInternalTensor(i, j) * tInternalTensor(i, j));
}
};
//! [J2 2D]
/** \brief Paraboloidal yield criterion
*
* See paper for reference \cite ullah2017three and \cite stier2015comparing
*
* Linear hardening
* \f[
* \psi =
* \frac{1}{2} \lambda \textrm{tr}[\varepsilon]^2 + \mu \varepsilon :
* \varepsilon
* +
* \sigma_{yt}\alpha_0
* +
* \frac{1}{2} H_t \alpha_0^2
* +
* \sigma_{yc}\alpha_1
* +
* \frac{1}{2} H_c \alpha_1^2
* \f]
*
* Exponential hardening
* \f[
* \psi =
* \frac{1}{2} \lambda \textrm{tr}[\varepsilon]^2 + \mu \varepsilon :
* \varepsilon\\
* +
* (\sigma_{yt}+H_t)\alpha_0
* +
* \frac{H_t}{n_t} \exp{(-n_t \alpha_0)}\\
* +
* (\sigma_{yc}+H_c)\alpha_0
* +
* \frac{H_c}{n_c} \exp{(-n_c \alpha_1)}
* \f]
*
* \f[
* I_1 = \textrm{tr} (\boldsymbol{\sigma})
* \f]
*
* \f[
* \eta=\textrm{dev}[\boldsymbol{\sigma}]
* \f]
*
* \f[
* J_2 = \frac{1}{2} \eta:\eta
* \f]
*
* \f[
* y = 6J_2 + 2I_1\left(\overline{\alpha_1}-\overline{\alpha_0}\right) -
* 2\overline{\alpha_0} \,\overline{\alpha_1} \f]
* where
* \f[
* \overline{\alpha_0}=\frac{\partial \psi}{\partial \alpha_0}=\sigma_{yt} + H_t
* \alpha_0
* \f]
*
* \f[
* \overline{\alpha_1}=\frac{\partial \psi}{\partial \alpha_1}=\sigma_{yc} + H_c
* \alpha_1
* \f]
*
* \f[
* \Psi = 6J_2 + 2\alpha I_1 \left(\overline{\alpha_1}-\overline{\alpha_0}\right)
* - 2\overline{\alpha_0} \,\overline{\alpha_1} \f] \f[ \alpha=
* \frac{1-2\nu_p}{1+\nu_p}
* \f]
*
* \ingroup adoc_plasticity
*/
struct ParaboloidalPlasticity : public ClosestPointProjection {
}
double lambda; //< Lamé parameters
double mu; //< Lamé parameter
double nup; //< Plastic Poisson’s ratio
double Ht, Hc; //< Isotropic hardening for tension and compression
double sigmaYt, sigmaYc; //< Yield stress in tension and compression
double nt, nc; //< Hardening parameters for tension and compression
enum ModelParams {
MU,
NUP,
HT,
HC,
NT,
NC,
};
using ParamsArray = std::array<double, LAST_PARAM>;
boost::shared_ptr<ParamsArray> defaultParamsArrayPtr = nullptr;
double young_modulus = 1;
double poisson_ratio = 0.25;
sigmaYt = 1;
sigmaYc = 1;
Ht = 0.1;
Hc = 0.1;
nc = 1.;
nt = 1.;
nup = 0.;
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-young_modulus", &young_modulus,
PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-poisson_ratio", &poisson_ratio,
0);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-sigma_yt", &sigmaYt, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-sigma_yc", &sigmaYc, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-Ht", &Ht, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-Hc", &Hc, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-nt", &nt, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-nc", &nc, PETSC_NULL);
CHKERR PetscOptionsGetScalar(PETSC_NULL, "-nup", &nup, PETSC_NULL);
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "Young modulus: " << young_modulus;
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "Poisson ratio: " << poisson_ratio;
defaultParamsArrayPtr = boost::make_shared<ParamsArray>();
(*defaultParamsArrayPtr)[LAMBDA] = lambda;
(*defaultParamsArrayPtr)[MU] = mu;
(*defaultParamsArrayPtr)[NUP] = nup;
(*defaultParamsArrayPtr)[SIGMA_YT] = sigmaYt;
(*defaultParamsArrayPtr)[SIGMA_YC] = sigmaYc;
(*defaultParamsArrayPtr)[HT] = Ht;
(*defaultParamsArrayPtr)[HC] = Hc;
(*defaultParamsArrayPtr)[NT] = nt;
(*defaultParamsArrayPtr)[NC] = nc;
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "LAMBDA: " << (*defaultParamsArrayPtr)[LAMBDA];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "MU: " << (*defaultParamsArrayPtr)[MU];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "NUP: " << (*defaultParamsArrayPtr)[NUP];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "SIGMA_YT: " << (*defaultParamsArrayPtr)[SIGMA_YT];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "SIGMA_YC: " << (*defaultParamsArrayPtr)[SIGMA_YC];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "HT: " << (*defaultParamsArrayPtr)[HT];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "HC: " << (*defaultParamsArrayPtr)[HC];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "NT: " << (*defaultParamsArrayPtr)[NT];
MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
<< "NC: " << (*defaultParamsArrayPtr)[NC];
};
boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> &pip,
std::string block_name, Sev sev) {
return 0;
}
MoFEMErrorCode setParams(short tag, bool &recalculate_elastic_tangent) {
recalculate_elastic_tangent = true;
return 0;
}
//! [Paraboloidal free energy]
auto t_voight_strain_op = voight_to_strain_op();
auto t_voight_total_strain = getFTensor1FromPtr<6>(&a_sTrain[0]);
auto t_voight_plastic_strain = getFTensor1FromPtr<6>(&a_plasticStrain[0]);
tTotalStrain(i, j) = t_voight_strain_op(i, j, Z) * t_voight_total_strain(Z);
t_voight_strain_op(i, j, Z) * t_voight_plastic_strain(Z);
a_w = 0.5 * lambda * tR * tR +
}
//! [Paraboloidal free energy]
//! [Paraboloidal yield function]
auto t_voight_stress_op = voight_to_stress_op();
auto t_vioght_stress = getFTensor1FromPtr<6>(&a_sTress[0]);
tStress(i, j) = t_voight_stress_op(i, j, Z) * t_vioght_stress(Z);
tR = tStress(i, i);
I1 = tR;
tR /= 3.;
tDeviator(i, j) = tStress(i, j) - tR * t_kd(i, j);
J2 = 0.5 * (tDeviator(i, j) * tDeviator(i, j));
}
a_y = 6 * J2 + 2 * I1 * (a_internalFluxes[1] - a_internalFluxes[0]) -
}
//! [Paraboloidal yield function]
//! [Paraboloidal flow potential]
double alpha =
(1 - 2 * nup) /
(1 + nup); // relation between alpha and plastic poisson's ratio
a_h = 6 * J2 +
2 * alpha * I1 * (a_internalFluxes[1] - a_internalFluxes[0]) -
}
//! [Paraboloidal flow potential]
};
template <typename MODEL>
inline auto createMaterial(
std::array<int, ClosestPointProjection::LAST_TAPE> tapes_tags = {0, 1, 2}) {
auto cp_ptr = boost::make_shared<MODEL>();
cp_ptr->getDefaultMaterialParameters();
cp_ptr->createMatAVecR();
cp_ptr->snesCreate();
cp_ptr->recordTapes();
cp_ptr->tapesTags = tapes_tags;
return cp_ptr;
}
} // namespace ADOLCPlasticity
#endif //__ADOLC_MATERIAL_MODELS_HPP
NOSPACE
@ NOSPACE
Definition: definitions.h:83
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::K
VectorDouble K
Definition: Projection10NodeCoordsOnField.cpp:125
ADOLCPlasticity::ClosestPointProjection::a_h
adouble a_h
Definition: ADOLCPlasticity.hpp:330
ADOLCPlasticity::ParaboloidalPlasticity::defaultParamsArrayPtr
boost::shared_ptr< ParamsArray > defaultParamsArrayPtr
Definition: ADOLCPlasticityMaterialModels.hpp:571
young_modulus
double young_modulus
Young modulus.
Definition: plastic.cpp:121
ADOLCPlasticity::ParaboloidalPlasticity::tElasticStrain
FTensor::Tensor2_symmetric< adouble, 3 > tElasticStrain
Definition: ADOLCPlasticityMaterialModels.hpp:544
ADOLCPlasticity::ParaboloidalPlasticity::SIGMA_YC
@ SIGMA_YC
Definition: ADOLCPlasticityMaterialModels.hpp:562
ADOLCPlasticity::ClosestPointProjection::ClosestPointProjection
ClosestPointProjection()
Definition: ClosetPointProjection.cpp:15
sigmaY
double sigmaY
Yield stress.
Definition: plastic.cpp:123
CHK_THROW_MESSAGE
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
Definition: definitions.h:609
ADOLCPlasticity::ClosestPointProjection::a_plasticStrain
ublas::vector< adouble > a_plasticStrain
Definition: ADOLCPlasticity.hpp:325
ADOLCPlasticity::ParaboloidalPlasticity::HT
@ HT
Definition: ADOLCPlasticityMaterialModels.hpp:563
ADOLCPlasticity::ParaboloidalPlasticity::nt
double nt
Definition: ADOLCPlasticityMaterialModels.hpp:555
ADOLCPlasticity::ParaboloidalPlasticity::tR
adouble tR
Definition: ADOLCPlasticityMaterialModels.hpp:548
ADOLCPlasticity::ParaboloidalPlasticity::j
FTensor::Index< 'j', 3 > j
Definition: ADOLCPlasticityMaterialModels.hpp:539
phi
static double phi
Definition: poisson_2d_dis_galerkin.cpp:30
MoFEM::Exceptions::MoFEMErrorCode
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
Definition: Exceptions.hpp:56
ADOLCPlasticity::ParaboloidalPlasticity::Ht
double Ht
Definition: ADOLCPlasticityMaterialModels.hpp:553
ADOLCPlasticity::ParaboloidalPlasticity::sigmaYt
double sigmaYt
Definition: ADOLCPlasticityMaterialModels.hpp:554
ADOLCPlasticity::ParaboloidalPlasticity::Hc
double Hc
Definition: ADOLCPlasticityMaterialModels.hpp:553
LAMBDA
#define LAMBDA(E, NU)
Definition: fem_tools.h:22
ADOLCPlasticity::ClosestPointProjection::a_sTrain
ublas::vector< adouble > a_sTrain
Definition: ADOLCPlasticity.hpp:324
ADOLCPlasticity::ParaboloidalPlasticity::NUP
@ NUP
Definition: ADOLCPlasticityMaterialModels.hpp:560
W
double W
Definition: free_surface.cpp:166
FTensor::Tensor2_symmetric< adouble, 3 >
ADOLCPlasticity::ParaboloidalPlasticity::flowPotential
MoFEMErrorCode flowPotential()
[Paraboloidal yield function]
Definition: ADOLCPlasticityMaterialModels.hpp:705
ADOLCPlasticity::ParaboloidalPlasticity::mu
double mu
Definition: ADOLCPlasticityMaterialModels.hpp:551
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
ADOLCPlasticity::ParaboloidalPlasticity::NT
@ NT
Definition: ADOLCPlasticityMaterialModels.hpp:565
ADOLCPlasticity::ParaboloidalPlasticity::ModelParams
ModelParams
Definition: ADOLCPlasticityMaterialModels.hpp:557
c
const double c
speed of light (cm/ns)
Definition: initial_diffusion.cpp:39
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
ADOLCPlasticity::ParaboloidalPlasticity::nc
double nc
Definition: ADOLCPlasticityMaterialModels.hpp:555
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
ADOLCPlasticity::J2Plasticity< 3 >::J2Plasticity
J2Plasticity()
Definition: ADOLCPlasticityMaterialModels.hpp:75
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
ADOLCPlasticity::ParaboloidalPlasticity::getDefaultMaterialParameters
MoFEMErrorCode getDefaultMaterialParameters()
Definition: ADOLCPlasticityMaterialModels.hpp:573
BlockData
Definition: ElasticityMixedFormulation.hpp:12
ADOLCPlasticity::ParaboloidalPlasticity::evalF
MoFEMErrorCode evalF()
[Paraboloidal yield function]
Definition: ADOLCPlasticityMaterialModels.hpp:681
ADOLCPlasticity::ParaboloidalPlasticity::I1
adouble I1
[Paraboloidal free energy]
Definition: ADOLCPlasticityMaterialModels.hpp:678
convert.type
type
Definition: convert.py:64
ADOLCPlasticity::ParaboloidalPlasticity::i
FTensor::Index< 'i', 3 > i
Definition: ADOLCPlasticityMaterialModels.hpp:538
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:122
ADOLCPlasticity::createMaterial
auto createMaterial(std::array< int, ClosestPointProjection::LAST_TAPE > tapes_tags={0, 1, 2})
Definition: ADOLCPlasticityMaterialModels.hpp:721
ADOLCPlasticity::voight_to_strain_op
FTensor::Dg< double, 3, 6 > voight_to_strain_op()
Op convert Vight strain vector to strain tensor.
Definition: ADOLCPlasticity.hpp:29
ADOLCPlasticity::ParaboloidalPlasticity::Z
FTensor::Index< 'Z', 6 > Z
Definition: ADOLCPlasticityMaterialModels.hpp:540
ADOLCPlasticity::ParaboloidalPlasticity::ParaboloidalPlasticity
ParaboloidalPlasticity()
Definition: ADOLCPlasticityMaterialModels.hpp:533
ADOLCPlasticity::ParaboloidalPlasticity::tPlasticStrain
FTensor::Tensor2_symmetric< adouble, 3 > tPlasticStrain
Definition: ADOLCPlasticityMaterialModels.hpp:543
ADOLCPlasticity::ParaboloidalPlasticity::addMatBlockOps
MoFEMErrorCode addMatBlockOps(MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, Sev sev)
Get model parameters from blocks.
Definition: ADOLCPlasticityMaterialModels.hpp:641
ADOLCPlasticity::ClosestPointProjection::a_w
adouble a_w
Definition: ADOLCPlasticity.hpp:328
ADOLCPlasticity
Definition: ADOLCPlasticity.hpp:24
ADOLCPlasticity::ParaboloidalPlasticity::nup
double nup
Definition: ADOLCPlasticityMaterialModels.hpp:552
i
FTensor::Index< 'i', SPACE_DIM > i
Definition: hcurl_divergence_operator_2d.cpp:27
ADOLCPlasticity::ClosestPointProjection::a_internalFluxes
ublas::vector< adouble > a_internalFluxes
Definition: ADOLCPlasticity.hpp:327
t_kd
constexpr auto t_kd
Definition: free_surface.cpp:137
ADOLCPlasticity::ClosestPointProjection::nbInternalVariables
int nbInternalVariables
Definition: ADOLCPlasticity.hpp:141
ADOLCPlasticity::ClosestPointProjection::activeVariablesW
VectorAdaptor activeVariablesW
Definition: ADOLCPlasticity.hpp:173
EntData
EntitiesFieldData::EntData EntData
Definition: child_and_parent.cpp:37
FTensor::Index< 'i', 3 >
ADOLCPlasticity::ParaboloidalPlasticity::tDeviator
FTensor::Tensor2_symmetric< adouble, 3 > tDeviator
Definition: ADOLCPlasticityMaterialModels.hpp:546
ADOLCPlasticity::ParaboloidalPlasticity::LAST_PARAM
@ LAST_PARAM
Definition: ADOLCPlasticityMaterialModels.hpp:567
Range
adouble
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
MOFEM_LOG
#define MOFEM_LOG(channel, severity)
Log.
Definition: LogManager.hpp:308
ADOLCPlasticity::ParaboloidalPlasticity::tTotalStrain
FTensor::Tensor2_symmetric< adouble, 3 > tTotalStrain
Definition: ADOLCPlasticityMaterialModels.hpp:542
ADOLCPlasticity::ClosestPointProjection::a_internalVariables
ublas::vector< adouble > a_internalVariables
Definition: ADOLCPlasticity.hpp:326
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
UserDataOperator
ForcesAndSourcesCore::UserDataOperator UserDataOperator
Definition: HookeElement.hpp:75
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
ADOLCPlasticity::ClosestPointProjection::a_y
adouble a_y
Definition: ADOLCPlasticity.hpp:329
mu
double mu
Definition: dynamic_first_order_con_law.cpp:98
lambda
static double lambda
Definition: incompressible_elasticity.cpp:199
ADOLCPlasticity::ParaboloidalPlasticity::LAMBDA
@ LAMBDA
Definition: ADOLCPlasticityMaterialModels.hpp:558
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
ADOLCPlasticity::voight_to_stress_op
FTensor::Dg< double, 3, 6 > voight_to_stress_op()
Op convert Vight stress vector to stress tensor.
Definition: ADOLCPlasticity.hpp:49
ADOLCPlasticity::ClosestPointProjection::a_sTress
ublas::vector< adouble > a_sTress
Definition: ADOLCPlasticity.hpp:327
ADOLCPlasticity::ParaboloidalPlasticity::ParamsArray
std::array< double, LAST_PARAM > ParamsArray
Definition: ADOLCPlasticityMaterialModels.hpp:570
ADOLCPlasticity::ParaboloidalPlasticity::MU
@ MU
Definition: ADOLCPlasticityMaterialModels.hpp:559
ADOLCPlasticity::ParaboloidalPlasticity::sigmaYc
double sigmaYc
Definition: ADOLCPlasticityMaterialModels.hpp:554
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
ADOLCPlasticity::ParaboloidalPlasticity::setParams
MoFEMErrorCode setParams(short tag, bool &recalculate_elastic_tangent)
Set parameters for ADOL-C of constitutive relations.
Definition: ADOLCPlasticityMaterialModels.hpp:647
ADOLCPlasticity::ParaboloidalPlasticity::yieldFunction
MoFEMErrorCode yieldFunction()
Set yield function.
Definition: ADOLCPlasticityMaterialModels.hpp:695
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
ADOLCPlasticity::ParaboloidalPlasticity::SIGMA_YT
@ SIGMA_YT
Definition: ADOLCPlasticityMaterialModels.hpp:561
MoFEM::PetscOptionsGetScalar
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
Definition: DeprecatedPetsc.hpp:162
third
constexpr double third
Definition: EshelbianADOL-C.cpp:17
ADOLCPlasticity::ParaboloidalPlasticity::NC
@ NC
Definition: ADOLCPlasticityMaterialModels.hpp:566
ADOLCPlasticity::ParaboloidalPlasticity::HC
@ HC
Definition: ADOLCPlasticityMaterialModels.hpp:564
ADOLCPlasticity::ParaboloidalPlasticity::tStress
FTensor::Tensor2_symmetric< adouble, 3 > tStress
Definition: ADOLCPlasticityMaterialModels.hpp:545
ADOLCPlasticity::ParaboloidalPlasticity::freeHemholtzEnergy
MoFEMErrorCode freeHemholtzEnergy()
[Paraboloidal free energy]
Definition: ADOLCPlasticityMaterialModels.hpp:653
OP
ADOLCPlasticity::ParaboloidalPlasticity::lambda
double lambda
Definition: ADOLCPlasticityMaterialModels.hpp:550
MoFEMFunctionReturn
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:429
H
double H
Hardening.
Definition: plastic.cpp:124
MoFEMFunctionBegin
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:359
MU
#define MU(E, NU)
Definition: fem_tools.h:23
ADOLCPlasticity::ParaboloidalPlasticity::J2
adouble J2
Definition: ADOLCPlasticityMaterialModels.hpp:678
PlasticOps::addMatBlockOps
MoFEMErrorCode addMatBlockOps(MoFEM::Interface &m_field, std::string block_name, Pip &pip, boost::shared_ptr< MatrixDouble > mat_D_Ptr, boost::shared_ptr< CommonData::BlockParams > mat_params_ptr, double scale_value, Sev sev)
Definition: PlasticOps.hpp:248