v0.15.0
Loading...
Searching...
No Matches
Public Member Functions | Protected Attributes | List of all members
PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp > Struct Template Reference

#include "users_modules/multifield-thermoplasticity-private/tutorials/adv-0/src/PlasticOpsGeneric.hpp"

Inheritance diagram for PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >:
[legend]
Collaboration diagram for PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >:
[legend]

Public Member Functions

 OpCalculatePlasticityWithoutRatesImpl (const std::string field_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< MatrixDouble > m_D_ptr, boost::shared_ptr< ThermoPlasticOps::ThermoPlasticBlockedParameters > tp_common_data_ptr=nullptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 

Protected Attributes

boost::shared_ptr< CommonDatacommonDataPtr
 
boost::shared_ptr< ThermoPlasticOps::ThermoPlasticBlockedParametersTPCommonDataPtr
 
boost::shared_ptr< MatrixDouble > mDPtr
 

Detailed Description

template<int DIM, typename DomainEleOp>
struct PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >

Definition at line 882 of file PlasticOpsGeneric.hpp.

Constructor & Destructor Documentation

◆ OpCalculatePlasticityWithoutRatesImpl()

template<int DIM, typename DomainEleOp >
PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >::OpCalculatePlasticityWithoutRatesImpl ( const std::string  field_name,
boost::shared_ptr< CommonData common_data_ptr,
boost::shared_ptr< MatrixDouble >  m_D_ptr,
boost::shared_ptr< ThermoPlasticOps::ThermoPlasticBlockedParameters tp_common_data_ptr = nullptr 
)

Definition at line 899 of file PlasticOpsGeneric.hpp.

904 : DomainEleOp(field_name, DomainEleOp::OPROW),
905 commonDataPtr(common_data_ptr), mDPtr(m_D_ptr),
906 TPCommonDataPtr(tp_common_data_ptr) {
907 // Opetor is only executed for vertices
908 std::fill(&DomainEleOp::doEntities[MBEDGE],
909 &DomainEleOp::doEntities[MBMAXTYPE], false);
910}
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
constexpr auto field_name
boost::shared_ptr< ThermoPlasticOps::ThermoPlasticBlockedParameters > TPCommonDataPtr

Member Function Documentation

◆ doWork()

template<int DIM, typename DomainEleOp >
MoFEMErrorCode PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >::doWork ( int  side,
EntityType  type,
EntData data 
)

< material parameters

Definition at line 913 of file PlasticOpsGeneric.hpp.

914 {
916
917 FTensor::Index<'i', DIM> i;
918 FTensor::Index<'j', DIM> j;
919 FTensor::Index<'k', DIM> k;
920 FTensor::Index<'l', DIM> l;
921 FTensor::Index<'m', DIM> m;
922 FTensor::Index<'n', DIM> n;
923
924 auto &params = commonDataPtr->blockParams; ///< material parameters
925
926 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
927 auto t_w = DomainEleOp::getFTensor0IntegrationWeight();
928 auto t_tau = getFTensor0FromVec(commonDataPtr->plasticTau);
929 auto t_tau_dot = getFTensor0FromVec(commonDataPtr->plasticTauDot);
930 auto t_f = getFTensor0FromVec(commonDataPtr->plasticSurface);
931 auto t_flow = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->plasticFlow);
932 auto t_plastic_strain =
933 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->plasticStrain);
934 auto t_stress =
935 getFTensor2SymmetricFromMat<DIM>(*(commonDataPtr->mStressPtr));
936
937 auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*commonDataPtr->mDPtr);
938 auto t_D_Op = getFTensor4DdgFromMat<DIM, DIM, 0>(*mDPtr);
939
940 auto t_temp = getFTensor0FromVec(TPCommonDataPtr->temperature);
941
942 auto t_diff_plastic_strain = diff_tensor(FTensor::Number<DIM>());
943 auto t_diff_deviator = diff_deviator(diff_tensor(FTensor::Number<DIM>()));
944
945 FTensor::Ddg<double, DIM, DIM> t_flow_dir_dstress;
946 FTensor::Ddg<double, DIM, DIM> t_flow_dir_dstrain;
947 t_flow_dir_dstress(i, j, k, l) =
948 1.5 * (t_diff_deviator(M, N, i, j) * t_diff_deviator(M, N, k, l));
949 t_flow_dir_dstrain(i, j, k, l) =
950 t_flow_dir_dstress(i, j, m, n) * t_D_Op(m, n, k, l);
951
952 auto t_alpha_dir =
953 kinematic_hardening_dplastic_strain<DIM>(params[CommonData::C1_k]);
954
955 commonDataPtr->resC.resize(nb_gauss_pts, false);
956 commonDataPtr->resCdTau.resize(nb_gauss_pts, false);
957 commonDataPtr->resCdStrain.resize(size_symm, nb_gauss_pts, false);
958 commonDataPtr->resCdPlasticStrain.resize(size_symm, nb_gauss_pts, false);
959 TPCommonDataPtr->resCdTemperature.resize(nb_gauss_pts, false);
960 commonDataPtr->resFlow.resize(size_symm, nb_gauss_pts, false);
961 commonDataPtr->resFlowDtau.resize(size_symm, nb_gauss_pts, false);
962 commonDataPtr->resFlowDstrain.resize(size_symm * size_symm, nb_gauss_pts,
963 false);
964 commonDataPtr->resFlowDstrainDot.resize(size_symm * size_symm, nb_gauss_pts,
965 false);
966 TPCommonDataPtr->resFlowDtemp.resize(size_symm, nb_gauss_pts, false);
967
968 commonDataPtr->resC.clear();
969 commonDataPtr->resCdTau.clear();
970 commonDataPtr->resCdStrain.clear();
971 commonDataPtr->resCdPlasticStrain.clear();
972 TPCommonDataPtr->resCdTemperature.clear();
973 commonDataPtr->resFlow.clear();
974 commonDataPtr->resFlowDtau.clear();
975 commonDataPtr->resFlowDstrain.clear();
976 commonDataPtr->resFlowDstrainDot.clear();
977 TPCommonDataPtr->resFlowDtemp.clear();
978
979 auto t_res_c = getFTensor0FromVec(commonDataPtr->resC);
980 auto t_res_c_dtau = getFTensor0FromVec(commonDataPtr->resCdTau);
981 auto t_res_c_dstrain =
982 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resCdStrain);
983 auto t_res_c_plastic_strain =
984 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resCdPlasticStrain);
985 auto t_res_c_temperature =
986 getFTensor0FromVec(TPCommonDataPtr->resCdTemperature);
987 auto t_res_flow = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resFlow);
988 auto t_res_flow_dtau =
989 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resFlowDtau);
990 auto t_res_flow_dstrain =
991 getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->resFlowDstrain);
992 auto t_res_flow_dplastic_strain =
993 getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->resFlowDstrainDot);
994 auto t_res_flow_dtemp =
995 getFTensor2SymmetricFromMat<DIM>(TPCommonDataPtr->resFlowDtemp);
996
997 auto next = [&]() {
998 ++t_tau;
999 ++t_tau_dot;
1000 ++t_f;
1001 ++t_flow;
1002 ++t_plastic_strain;
1003 ++t_stress;
1004 ++t_res_c;
1005 ++t_res_c_dtau;
1006 ++t_res_c_dstrain;
1007 ++t_res_c_plastic_strain;
1008 ++t_res_c_temperature;
1009 ++t_res_flow;
1010 ++t_res_flow_dtau;
1011 ++t_res_flow_dstrain;
1012 ++t_res_flow_dplastic_strain;
1013 ++t_res_flow_dtemp;
1014 ++t_w;
1015 ++t_temp;
1016 };
1017
1018 auto get_avtive_pts = [&]() {
1019 int nb_points_avtive_on_elem = 0;
1020 int nb_points_on_elem = 0;
1021
1022 auto t_tau = getFTensor0FromVec(commonDataPtr->plasticTau);
1023 auto t_tau_dot = getFTensor0FromVec(commonDataPtr->plasticTauDot);
1024 auto t_f = getFTensor0FromVec(commonDataPtr->plasticSurface);
1025 auto t_plastic_strain_dot =
1026 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->plasticStrainDot);
1027 auto t_temp = getFTensor0FromVec(TPCommonDataPtr->temperature);
1028
1029 auto dt = this->getTStimeStep();
1030
1031 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1032 auto eqiv = equivalent_strain_dot(t_plastic_strain_dot);
1033 const auto ww =
1034 w(eqiv, t_tau_dot, t_f,
1035 iso_hardening(t_tau, params[CommonData::H], TPCommonDataPtr->omega0,
1036 params[CommonData::QINF], TPCommonDataPtr->omegaH,
1037 params[CommonData::BISO], params[CommonData::SIGMA_Y],
1038 TPCommonDataPtr->temp0, t_temp),
1039 params[CommonData::SIGMA_Y]);
1040 const auto sign_ww = constrian_sign(ww, dt);
1041
1042 ++nb_points_on_elem;
1043 if (sign_ww > 0) {
1044 ++nb_points_avtive_on_elem;
1045 }
1046
1047 ++t_tau;
1048 ++t_tau_dot;
1049 ++t_f;
1050 ++t_plastic_strain_dot;
1051 ++t_temp;
1052 }
1053
1054 int &active_points = PlasticOps::CommonData::activityData[0];
1055 int &avtive_full_elems = PlasticOps::CommonData::activityData[1];
1056 int &avtive_elems = PlasticOps::CommonData::activityData[2];
1057 int &nb_points = PlasticOps::CommonData::activityData[3];
1058 int &nb_elements = PlasticOps::CommonData::activityData[4];
1059
1060 ++nb_elements;
1061 nb_points += nb_points_on_elem;
1062 if (nb_points_avtive_on_elem > 0) {
1063 ++avtive_elems;
1064 active_points += nb_points_avtive_on_elem;
1065 if (nb_points_avtive_on_elem == nb_points_on_elem) {
1066 ++avtive_full_elems;
1067 }
1068 }
1069
1070 if (nb_points_avtive_on_elem != nb_points_on_elem)
1071 return 1;
1072 else
1073 return 0;
1074 };
1075
1076 if (DomainEleOp::getTSCtx() == TSMethod::TSContext::CTX_TSSETIJACOBIAN) {
1077 get_avtive_pts();
1078 }
1079
1080 for (auto gg = 0; gg != nb_gauss_pts; ++gg) {
1081
1082 double eqiv;
1084
1085 const auto sigma_y =
1086 iso_hardening(t_tau, params[CommonData::H], TPCommonDataPtr->omega0,
1087 params[CommonData::QINF], TPCommonDataPtr->omegaH,
1088 params[CommonData::BISO], params[CommonData::SIGMA_Y],
1089 TPCommonDataPtr->temp0, t_temp);
1090 const auto d_sigma_y = iso_hardening_dtau(
1091 t_tau, params[CommonData::H], TPCommonDataPtr->omega0,
1092 params[CommonData::QINF], TPCommonDataPtr->omegaH,
1093 params[CommonData::BISO], params[CommonData::SIGMA_Y],
1094 TPCommonDataPtr->temp0, t_temp);
1095
1096 double c_dot_tau, c_sigma_y, c, c_f, c_temperature, c_equiv;
1097
1098 eqiv = 0.0;
1099 t_diff_eqiv(i, j) = 0.0;
1100 c = 0.0;
1101 c_dot_tau = 0.0;
1102 c_equiv = 0.0;
1103 c_sigma_y = 0.0;
1104 c_f = 0.0;
1105
1106 auto d_sigma_y_dtemp = iso_hardening_dtemp(
1107 t_tau, params[CommonData::H], TPCommonDataPtr->omega0,
1108 params[CommonData::QINF], TPCommonDataPtr->omegaH,
1109 params[CommonData::BISO], params[CommonData::SIGMA_Y],
1110 TPCommonDataPtr->temp0, t_temp);
1111
1112 auto t_dev_stress = deviator(
1113
1114 t_stress, trace(t_stress),
1115
1116 kinematic_hardening(t_plastic_strain, params[CommonData::C1_k])
1117
1118 );
1119
1120 next();
1121 }
1122
1124}
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
FTensor::Index< 'i', SPACE_DIM > i
double dt
const double c
speed of light (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
FTensor::Index< 'N', 3 > N
double constrian_sign(double x, double dt)
FTensor::Index< 'M', 3 > M
auto diff_tensor(FTensor::Number< DIM >)
[Lambda functions]
auto diff_deviator(FTensor::Ddg< double, DIM, DIM > &&t_diff_stress, FTensor::Number< DIM >)
double trace(FTensor::Tensor2_symmetric< T, 2 > &t_stress)
auto deviator(FTensor::Tensor2_symmetric< T, DIM > &t_stress, double trace, FTensor::Tensor2_symmetric< double, DIM > &t_alpha, FTensor::Number< DIM >)
double w(double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y)
auto equivalent_strain_dot(FTensor::Tensor2_symmetric< T, DIM > &t_plastic_strain_dot)
FTensor::Index< 'm', 3 > m
static std::array< int, 5 > activityData
double iso_hardening_dtemp(double tau, double H, double omega_0, double Qinf, double omega_h, double b_iso, double sigmaY, double temp_0, double temp)
auto kinematic_hardening(FTensor::Tensor2_symmetric< T, DIM > &t_plastic_strain, double C1_k)
Definition plastic.cpp:92
double iso_hardening_dtau(double tau, double H, double Qinf, double b_iso)
Definition plastic.cpp:78
constexpr auto size_symm
Definition plastic.cpp:42
double iso_hardening(double tau, double H, double Qinf, double b_iso, double sigmaY)
Definition plastic.cpp:73

Member Data Documentation

◆ commonDataPtr

template<int DIM, typename DomainEleOp >
boost::shared_ptr<CommonData> PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >::commonDataPtr
protected

Definition at line 892 of file PlasticOpsGeneric.hpp.

◆ mDPtr

template<int DIM, typename DomainEleOp >
boost::shared_ptr<MatrixDouble> PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >::mDPtr
protected

Definition at line 895 of file PlasticOpsGeneric.hpp.

◆ TPCommonDataPtr

template<int DIM, typename DomainEleOp >
boost::shared_ptr<ThermoPlasticOps::ThermoPlasticBlockedParameters> PlasticOps::OpCalculatePlasticityWithoutRatesImpl< DIM, GAUSS, DomainEleOp >::TPCommonDataPtr
protected

Definition at line 894 of file PlasticOpsGeneric.hpp.


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