v0.15.0
Loading...
Searching...
No Matches
ADOLCPlasticity::HeterogeneousParaboloidalPlasticity Struct Reference

#include "users_modules/adolc-plasticity/src/ADOLCPlasticityMaterialModels.hpp"

Inheritance diagram for ADOLCPlasticity::HeterogeneousParaboloidalPlasticity:
[legend]
Collaboration diagram for ADOLCPlasticity::HeterogeneousParaboloidalPlasticity:
[legend]

Public Member Functions

MoFEMErrorCode addMatBlockOps (MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, Sev sev) override
 Get model parameters from blocks.
 
 ParaboloidalPlasticity ()
 
- Public Member Functions inherited from ADOLCPlasticity::ParaboloidalPlasticity
 ParaboloidalPlasticity ()
 
MoFEMErrorCode getDefaultMaterialParameters ()
 
MoFEMErrorCode addMatBlockOps (MoFEM::Interface &m_field, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string block_name, Sev sev)
 Get model parameters from blocks.
 
MoFEMErrorCode setParams (short tag, bool &recalculate_elastic_tangent)
 Set parameters for ADOL-C of constitutive relations.
 
MoFEMErrorCode freeHemholtzEnergy ()
 [Paraboloidal free energy]
 
MoFEMErrorCode codedHessianW (vector< double * > hessian)
 [Paraboloidal free energy]
 
MoFEMErrorCode evalF ()
 [Paraboloidal yield function]
 
MoFEMErrorCode yieldFunction ()
 Set yield function.
 
MoFEMErrorCode flowPotential ()
 [Paraboloidal yield function]
 
- Public Member Functions inherited from ADOLCPlasticity::ClosestPointProjection
VectorAdaptor getPlasticStrain ()
 
VectorAdaptor getTotalStrain ()
 
VectorAdaptor getInternalVariables ()
 
VectorAdaptor getActiveVariablesYH ()
 
VectorAdaptor getStress ()
 
VectorAdaptor getInternalFluxes ()
 
 ClosestPointProjection ()
 
MoFEMErrorCode recordW ()
 Record strain energy.
 
MoFEMErrorCode recordY ()
 Record yield function.
 
MoFEMErrorCode recordH ()
 Record flow potential.
 
MoFEMErrorCode playW ()
 
MoFEMErrorCode playW_NoHessian ()
 
MoFEMErrorCode playY ()
 
MoFEMErrorCode playY_NoGradient ()
 
MoFEMErrorCode playH ()
 
MoFEMErrorCode playH_NoHessian ()
 
MoFEMErrorCode createMatAVecR ()
 
MoFEMErrorCode evaluatePotentials ()
 
MoFEMErrorCode playPotentials ()
 
MoFEMErrorCode playPotentials_NoHessian ()
 
MoFEMErrorCode calculateR (Vec R)
 
MoFEMErrorCode calculateA ()
 
MoFEMErrorCode snesCreate ()
 Create nested snes.
 
MoFEMErrorCode solveClosestProjection ()
 Solve nonlinear system of equations to find stress, internal fluxes, and Lagrange plastic multiplier.
 
MoFEMErrorCode consistentTangent ()
 Calculate consistent tangent matrix.
 
MoFEMErrorCode recordTapes ()
 Record tapes.
 

Additional Inherited Members

- Public Types inherited from ADOLCPlasticity::ParaboloidalPlasticity
enum  ModelParams {
  LAMBDA , MU , NUP , SIGMA_YT ,
  SIGMA_YC , HT , HC , NT ,
  NC , LAST_PARAM
}
 
using ParamsArray = std::array<double, LAST_PARAM>
 
- Public Types inherited from ADOLCPlasticity::ClosestPointProjection
enum  TypesTags { W = 0 , Y , H , LAST_TAPE }
 
- Public Attributes inherited from ADOLCPlasticity::ParaboloidalPlasticity
FTensor::Index< 'i', 3 > i
 
FTensor::Index< 'j', 3 > j
 
FTensor::Index< 'Z', 6 > Z
 
FTensor::Tensor2_symmetric< adouble, 3 > tTotalStrain
 
FTensor::Tensor2_symmetric< adouble, 3 > tPlasticStrain
 
FTensor::Tensor2_symmetric< adouble, 3 > tElasticStrain
 
FTensor::Tensor2_symmetric< adouble, 3 > tStress
 
FTensor::Tensor2_symmetric< adouble, 3 > tDeviator
 
adouble tR
 
double lambda
 
double mu
 
double nup
 
double Ht
 
double Hc
 
double sigmaYt
 
double sigmaYc
 
double nt
 
double nc
 
boost::shared_ptr< ParamsArraydefaultParamsArrayPtr = nullptr
 
boost::shared_ptr< ParamsArrayparamsArrayPtr = nullptr
 
boost::shared_ptr< ParamsArrayoldParamsArrayPtr = nullptr
 
adouble I1
 
adouble J2
 
- Public Attributes inherited from ADOLCPlasticity::ClosestPointProjection
int nbInternalVariables
 
boost::function< int(int, int, int)> integrationRule
 
VectorDouble internalVariables0
 
VectorDouble plasticStrain0
 
std::array< int, LAST_TAPEtapesTags
 
VectorAdaptor activeVariablesW
 
VectorAdaptor gradientW
 
double w
 
double y
 
double h
 
double deltaGamma
 
MatrixDouble Ep
 
MatrixDouble Cp
 
MatrixDouble Cep
 
ublas::symmetric_matrix< double, ublas::lower > C
 
ublas::symmetric_matrix< double, ublas::lower > D
 
PetscBool implHessianW
 
MatrixDouble hessianW
 
VectorDouble gradientY
 
VectorDouble gradientH
 
MatrixDouble hessianH
 
MatrixDouble partialWStrainPlasticStrain
 
VectorAdaptor partialYSigma
 
VectorAdaptor partialYFlux
 
VectorAdaptor partialHSigma
 
VectorAdaptor partialHFlux
 
ublas::symmetric_matrix< double, ublas::lower > partial2HSigma
 
ublas::symmetric_matrix< double, ublas::lower > partial2HFlux
 
MatrixDouble partial2HSigmaFlux
 
SmartPetscObj< Mat > A
 
SmartPetscObj< Vec > R
 
SmartPetscObj< Vec > Chi
 
SmartPetscObj< Vec > dChi
 
ublas::matrix< double, ublas::column_major > dataA
 
SmartPetscObj< SNES > sNes
 
ublas::vector< adoublea_sTrain
 
ublas::vector< adoublea_plasticStrain
 
ublas::vector< adoublea_internalVariables
 
ublas::vector< adoublea_sTress
 
ublas::vector< adoublea_internalFluxes
 
adouble a_w
 
adouble a_y
 
adouble a_h
 

Detailed Description

Definition at line 1003 of file ADOLCPlasticityMaterialModels.hpp.

Member Function Documentation

◆ addMatBlockOps()

MoFEMErrorCode ADOLCPlasticity::HeterogeneousParaboloidalPlasticity::addMatBlockOps ( MoFEM::Interface & m_field,
boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > & pip,
std::string block_name,
Sev sev )
inlineoverridevirtual

Get model parameters from blocks.

Reimplemented from ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 1008 of file ADOLCPlasticityMaterialModels.hpp.

1010 {
1012
1013 struct OpMatBlocks : public ForcesAndSourcesCore::UserDataOperator {
1014
1015 using OP = ForcesAndSourcesCore::UserDataOperator;
1016
1017 double poissonRatio = 0.25;
1018 Tag thYoungModulus;
1019 Tag thCompressiveYieldStress;
1020 Tag thTensionYieldStress;
1021 double youngsModulus;
1022 double yieldStrengthC;
1023 double yieldStrengthT;
1024
1025 OpMatBlocks(boost::shared_ptr<ParamsArray> params_array_ptr,
1026 boost::shared_ptr<ParamsArray> def_params_array_ptr,
1027 MoFEM::Interface &m_field, Sev sev,
1028 std::vector<const CubitMeshSets *> meshset_vec_ptr)
1029 : OP(NOSPACE, OP::OPSPACE), paramsArrayPtr(params_array_ptr),
1030 defParamsArray(def_params_array_ptr), mField(m_field) {
1031 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-poisson_ratio",
1032 &poissonRatio, 0);
1033 rval = mField.get_moab().tag_get_handle(
1034 "YOUNGS_MODULUS", 1, MB_TYPE_DOUBLE, thYoungModulus, MB_TAG_SPARSE);
1035 if (rval != MB_SUCCESS) {
1036 MOFEM_LOG("ADOLCPlasticityWold", Sev::warning)
1037 << "YOUNGS_MODULUS tag does not exist using default value";
1038 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-young_modulus",
1039 &youngsModulus, 0);
1040 }
1041
1042 rval = mField.get_moab().tag_get_handle(
1043 "Yield_Strength_C", 1, MB_TYPE_DOUBLE, thCompressiveYieldStress,
1044 MB_TAG_SPARSE);
1045
1046 if (rval != MB_SUCCESS) {
1047 MOFEM_LOG("ADOLCPlasticityWold", Sev::warning)
1048 << "Yield_Strength_C tag does not exist using default value";
1049 yieldStrengthC = defParamsArray->at(SIGMA_YC);
1050 }
1051
1052 rval = mField.get_moab().tag_get_handle(
1053 "Yield_Strength_T", 1, MB_TYPE_DOUBLE, thTensionYieldStress,
1054 MB_TAG_SPARSE);
1055
1056 if (rval != MB_SUCCESS) {
1057 MOFEM_LOG("ADOLCPlasticityWold", Sev::warning)
1058 << "Yield_Strength_T tag does not exist using default value";
1059 yieldStrengthT = defParamsArray->at(SIGMA_YT);
1060 }
1061
1062 CHK_THROW_MESSAGE(extractBlockData(m_field, meshset_vec_ptr, sev),
1063 "Can not get data from block");
1064 }
1065
1066 MoFEMErrorCode doWork(int side, EntityType type,
1067 EntitiesFieldData::EntData &data) {
1069
1070 for (auto &b : blockData) {
1071
1072 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
1073 EntityHandle fe_ent = getFEEntityHandle();
1074
1075 if (thYoungModulus != 0)
1076 CHKERR mField.get_moab().tag_get_data(thYoungModulus, &fe_ent, 1,
1077 &youngsModulus);
1078
1079 if (thCompressiveYieldStress != 0)
1080 CHKERR mField.get_moab().tag_get_data(
1081 thCompressiveYieldStress, &fe_ent, 1, &yieldStrengthC);
1082
1083 if (thTensionYieldStress != 0)
1084 CHKERR mField.get_moab().tag_get_data(
1085 thTensionYieldStress, &fe_ent, 1, &yieldStrengthT);
1086
1087 std::copy(b.paramsArray.begin(), b.paramsArray.end(),
1088 paramsArrayPtr->begin());
1089
1090 (*paramsArrayPtr)[LAMBDA] = LAMBDA(youngsModulus, poissonRatio);
1091 (*paramsArrayPtr)[MU] = MU(youngsModulus, poissonRatio);
1092 (*paramsArrayPtr)[SIGMA_YC] = yieldStrengthC;
1093 (*paramsArrayPtr)[SIGMA_YT] = yieldStrengthT;
1094
1096 }
1097 }
1098
1099 std::copy(defParamsArray->begin(), defParamsArray->end(),
1100 paramsArrayPtr->begin());
1101
1103 }
1104
1105 private:
1106 boost::shared_ptr<ParamsArray> paramsArrayPtr;
1107 boost::shared_ptr<ParamsArray> defParamsArray;
1108 MoFEM::Interface &mField;
1109
1110 struct BlockData {
1111 ParamsArray paramsArray;
1112 Range blockEnts;
1113 };
1114
1115 std::vector<BlockData> blockData;
1116
1118 extractBlockData(MoFEM::Interface &m_field,
1119 std::vector<const CubitMeshSets *> meshset_vec_ptr,
1120 Sev sev) {
1122
1123 for (auto m : meshset_vec_ptr) {
1124 MOFEM_TAG_AND_LOG("ADOLCPlasticityWold", sev,
1125 "Heterogeneous Paraboloidal")
1126 << *m;
1127 std::vector<double> block_data;
1128 CHKERR m->getAttributes(block_data);
1129 if (block_data.size() != 0) {
1130 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1131 "Expected that block has zero attributes");
1132 }
1133 auto get_block_ents = [&]() {
1134 Range ents;
1135 CHKERR
1136 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
1137 return ents;
1138 };
1139
1140 blockData.push_back({*defParamsArray, get_block_ents()});
1141 std::copy(block_data.begin(), block_data.end(),
1142 blockData.back().paramsArray.begin());
1143 const auto young_modulus = blockData.back().paramsArray[LAMBDA];
1144 const auto poisson_ratio = blockData.back().paramsArray[MU];
1145
1146 // It is assumed that user provide young_modulus and poisson_ratio in
1147 // first two arguments of the block
1148
1149 MOFEM_LOG("ADOLCPlasticityWold", sev)
1150 << "Young modulus: " << (blockData.back().paramsArray)[LAMBDA];
1151 MOFEM_LOG("ADOLCPlasticityWold", sev)
1152 << "Poisson ratio: " << (blockData.back().paramsArray)[MU];
1153
1154 blockData.back().paramsArray[LAMBDA] =
1156 blockData.back().paramsArray[MU] = MU(young_modulus, poisson_ratio);
1157
1158 MOFEM_LOG("ADOLCPlasticityWold", sev)
1159 << "LAMBDA: " << (blockData.back().paramsArray)[LAMBDA];
1160 MOFEM_LOG("ADOLCPlasticityWold", sev)
1161 << "MU: " << (blockData.back().paramsArray)[MU];
1162 MOFEM_LOG("ADOLCPlasticityWold", sev)
1163 << "NUP: " << (blockData.back().paramsArray)[NUP];
1164 MOFEM_LOG("ADOLCPlasticityWold", sev)
1165 << "SIGMA_YT: " << (blockData.back().paramsArray)[SIGMA_YT];
1166 MOFEM_LOG("ADOLCPlasticityWold", sev)
1167 << "SIGMA_YC: " << (blockData.back().paramsArray)[SIGMA_YC];
1168 MOFEM_LOG("ADOLCPlasticityWold", sev)
1169 << "HT: " << (blockData.back().paramsArray)[HT];
1170 MOFEM_LOG("ADOLCPlasticityWold", sev)
1171 << "HC: " << (blockData.back().paramsArray)[HC];
1172 MOFEM_LOG("ADOLCPlasticityWold", sev)
1173 << "NT: " << (blockData.back().paramsArray)[NT];
1174 MOFEM_LOG("ADOLCPlasticityWold", sev)
1175 << "NC: " << (blockData.back().paramsArray)[NC];
1176 }
1178 }
1179 };
1180
1181 auto cubit_meshsets_vec_ptr =
1182 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
1183 std::regex((boost::format("%s(.*)") % block_name).str()));
1184
1185 pip.push_back(new OpMatBlocks(paramsArrayPtr, defaultParamsArrayPtr,
1186 m_field, sev, cubit_meshsets_vec_ptr));
1187
1189 }
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#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()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MU(E, NU)
Definition fem_tools.h:23
#define LAMBDA(E, NU)
Definition fem_tools.h:22
#define MOFEM_LOG(channel, severity)
Log.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)
double young_modulus
Young modulus.
Definition plastic.cpp:125
double poisson_ratio
Poisson ratio.
Definition plastic.cpp:126
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ ParaboloidalPlasticity()

ADOLCPlasticity::ParaboloidalPlasticity::ParaboloidalPlasticity ( )
inline

The documentation for this struct was generated from the following file: