v0.14.0 |
Typedefs | |
using | Pip = boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > |
using | CommonPlasticPtr = boost::shared_ptr< PlasticOps::CommonData > |
using | CommonHenckyPtr = boost::shared_ptr< HenckyOps::CommonData > |
Functions | |
template<int DIM> | |
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) |
template<int DIM, IntegrationType I, typename DomainEleOp > | |
auto | createCommonPlasticOps (MoFEM::Interface &m_field, std::string block_name, boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > &pip, std::string u, std::string ep, std::string tau, double scale, Sev sev) |
template<int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp > | |
MoFEMErrorCode | opFactoryDomainRhs (MoFEM::Interface &m_field, std::string block_name, Pip &pip, std::string u, std::string ep, std::string tau) |
template<int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp > | |
MoFEMErrorCode | opFactoryDomainLhs (MoFEM::Interface &m_field, std::string block_name, Pip &pip, std::string u, std::string ep, std::string tau) |
template<int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp > | |
MoFEMErrorCode | opFactoryDomainReactions (MoFEM::Interface &m_field, std::string block_name, Pip &pip, std::string u, std::string ep, std::string tau) |
template<int DIM> | |
auto | diff_tensor (FTensor::Number< DIM >) |
[Lambda functions] More... | |
template<int DIM> | |
auto | symm_L_tensor (FTensor::Number< DIM >) |
template<int DIM> | |
auto | diff_symmetrize (FTensor::Number< DIM >) |
template<typename T > | |
double | trace (FTensor::Tensor2_symmetric< T, 2 > &t_stress) |
template<typename T > | |
double | trace (FTensor::Tensor2_symmetric< T, 3 > &t_stress) |
template<typename T , int DIM> | |
auto | deviator (FTensor::Tensor2_symmetric< T, DIM > &t_stress, double trace, FTensor::Tensor2_symmetric< double, DIM > &t_alpha, FTensor::Number< DIM >) |
template<typename T > | |
auto | deviator (FTensor::Tensor2_symmetric< T, 2 > &t_stress, double trace, FTensor::Tensor2_symmetric< double, 2 > &&t_alpha) |
template<typename T > | |
auto | deviator (FTensor::Tensor2_symmetric< T, 3 > &t_stress, double trace, FTensor::Tensor2_symmetric< double, 3 > &&t_alpha) |
template<int DIM> | |
auto | diff_deviator (FTensor::Ddg< double, DIM, DIM > &&t_diff_stress, FTensor::Number< DIM >) |
auto | diff_deviator (FTensor::Ddg< double, 2, 2 > &&t_diff_stress) |
auto | diff_deviator (FTensor::Ddg< double, 3, 3 > &&t_diff_stress) |
double | platsic_surface (FTensor::Tensor2_symmetric< double, 3 > &&t_stress_deviator) |
template<int DIM> | |
auto | plastic_flow (long double f, FTensor::Tensor2_symmetric< double, 3 > &&t_dev_stress, FTensor::Ddg< double, 3, DIM > &&t_diff_deviator) |
template<typename T , int DIM> | |
auto | diff_plastic_flow_dstress (long double f, FTensor::Tensor2_symmetric< T, DIM > &t_flow, FTensor::Ddg< double, 3, DIM > &&t_diff_deviator) |
template<typename T , int DIM> | |
auto | diff_plastic_flow_dstrain (FTensor::Ddg< T, DIM, DIM > &t_D, FTensor::Ddg< double, DIM, DIM > &&t_diff_plastic_flow_dstress) |
double | constrian_sign (double x) |
double | constrain_abs (double x) |
double | w (double eqiv, double dot_tau, double f, double sigma_y, double sigma_Y) |
double | constraint (double eqiv, double dot_tau, double f, double sigma_y, double abs_w, double vis_H, double sigma_Y) |
double | diff_constrain_ddot_tau (double sign, double eqiv, double dot_tau, double vis_H, double sigma_Y) |
double | diff_constrain_deqiv (double sign, double eqiv, double dot_tau, double sigma_Y) |
auto | diff_constrain_df (double sign) |
auto | diff_constrain_dsigma_y (double sign) |
template<typename T , int DIM> | |
auto | diff_constrain_dstress (double diff_constrain_df, FTensor::Tensor2_symmetric< T, DIM > &t_plastic_flow) |
template<typename T1 , typename T2 , int DIM> | |
auto | diff_constrain_dstrain (T1 &t_D, T2 &&t_diff_constrain_dstress) |
template<typename T , int DIM> | |
auto | equivalent_strain_dot (FTensor::Tensor2_symmetric< T, DIM > &t_plastic_strain_dot) |
template<typename T1 , typename T2 , typename T3 , int DIM> | |
auto | diff_equivalent_strain_dot (const T1 eqiv, T2 &t_plastic_strain_dot, T3 &t_diff_plastic_strain, FTensor::Number< DIM >) |
static auto | get_mat_tensor_sym_dvector (size_t rr, MatrixDouble &mat, FTensor::Number< 2 >) |
[Lambda functions] More... | |
static auto | get_mat_tensor_sym_dvector (size_t rr, MatrixDouble &mat, FTensor::Number< 3 >) |
static auto | get_mat_tensor_sym_dtensor_sym (size_t rr, MatrixDouble &mat, FTensor::Number< 2 >) |
static auto | get_mat_tensor_sym_dtensor_sym (size_t rr, MatrixDouble &mat, FTensor::Number< 3 >) |
static auto | get_mat_tensor_sym_dscalar (size_t rr, MatrixDouble &mat, FTensor::Number< 2 >) |
static auto | get_mat_tensor_sym_dscalar (size_t rr, MatrixDouble &mat, FTensor::Number< 3 >) |
auto | get_mat_scalar_dtensor_sym (MatrixDouble &mat, FTensor::Number< 2 >) |
auto | get_mat_scalar_dtensor_sym (MatrixDouble &mat, FTensor::Number< 3 >) |
static FTensor::Tensor2< FTensor::PackPtr< double *, 3 >, 2, 3 > | get_mat_vector_dtensor_sym (size_t rr, MatrixDouble &mat, FTensor::Number< 2 >) |
static FTensor::Tensor2< FTensor::PackPtr< double *, 6 >, 3, 6 > | get_mat_vector_dtensor_sym (size_t rr, MatrixDouble &mat, FTensor::Number< 3 >) |
Variables | |
FTensor::Index< 'I', 3 > | I |
[Common data] More... | |
FTensor::Index< 'J', 3 > | J |
FTensor::Index< 'M', 3 > | M |
FTensor::Index< 'N', 3 > | N |
using PlasticOps::CommonHenckyPtr = typedef boost::shared_ptr<HenckyOps::CommonData> |
Definition at line 244 of file PlasticOps.hpp.
using PlasticOps::CommonPlasticPtr = typedef boost::shared_ptr<PlasticOps::CommonData> |
Definition at line 243 of file PlasticOps.hpp.
using PlasticOps::Pip = typedef boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator> |
Definition at line 242 of file PlasticOps.hpp.
MoFEMErrorCode PlasticOps::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 | ||
) |
Extract block data from meshsets
m_field | |
meshset_vec_ptr | |
sev |
Get elasticity tensor
Calculate elasticity tensor for given material parameters
mat_D_ptr | |
bulk_modulus_K | |
shear_modulus_G |
[Calculate elasticity tensor]
[Calculate elasticity tensor]
Definition at line 248 of file PlasticOps.hpp.
Definition at line 265 of file PlasticOpsGeneric.hpp.
|
inline |
\[ \dot{\tau} - \frac{1}{2}\left\{\dot{\tau} + (f(\pmb\sigma) - \sigma_y) + \| \dot{\tau} + (f(\pmb\sigma) - \sigma_y) \|\right\} = 0 \\ c_n \sigma_y \dot{\tau} - \frac{1}{2}\left\{c_n\sigma_y \dot{\tau} + (f(\pmb\sigma) - \sigma_y) + \| c_n \sigma_y \dot{\tau} + (f(\pmb\sigma) - \sigma_y) \|\right\} = 0 \]
Definition at line 292 of file PlasticOpsGeneric.hpp.
auto PlasticOps::createCommonPlasticOps | ( | MoFEM::Interface & | m_field, |
std::string | block_name, | ||
boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator > & | pip, | ||
std::string | u, | ||
std::string | ep, | ||
std::string | tau, | ||
double | scale, | ||
Sev | sev | ||
) |
|
inline |
Definition at line 113 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 119 of file PlasticOpsGeneric.hpp.
|
inline |
|
inline |
Definition at line 300 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 305 of file PlasticOpsGeneric.hpp.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Definition at line 157 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 161 of file PlasticOpsGeneric.hpp.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
Definition at line 43 of file PlasticOpsGeneric.hpp.
|
inline |
[Lambda functions]
Definition at line 10 of file PlasticOpsGeneric.hpp.
|
inline |
auto PlasticOps::get_mat_scalar_dtensor_sym | ( | MatrixDouble & | mat, |
FTensor::Number< 2 > | |||
) |
auto PlasticOps::get_mat_scalar_dtensor_sym | ( | MatrixDouble & | mat, |
FTensor::Number< 3 > | |||
) |
Definition at line 1186 of file PlasticOpsGeneric.hpp.
|
inlinestatic |
|
inlinestatic |
Definition at line 1096 of file PlasticOpsGeneric.hpp.
|
inlinestatic |
|
inlinestatic |
Definition at line 966 of file PlasticOpsGeneric.hpp.
|
inlinestatic |
[Lambda functions]
[Auxiliary functions functions
Definition at line 366 of file PlasticOpsGeneric.hpp.
|
inlinestatic |
Definition at line 373 of file PlasticOpsGeneric.hpp.
|
inlinestatic |
Definition at line 38 of file PlasticOpsSmallStrains.hpp.
|
inlinestatic |
Definition at line 45 of file PlasticOpsSmallStrains.hpp.
MoFEMErrorCode PlasticOps::opFactoryDomainLhs | ( | MoFEM::Interface & | m_field, |
std::string | block_name, | ||
Pip & | pip, | ||
std::string | u, | ||
std::string | ep, | ||
std::string | tau | ||
) |
MoFEMErrorCode PlasticOps::opFactoryDomainReactions | ( | MoFEM::Interface & | m_field, |
std::string | block_name, | ||
Pip & | pip, | ||
std::string | u, | ||
std::string | ep, | ||
std::string | tau | ||
) |
MoFEMErrorCode PlasticOps::opFactoryDomainRhs | ( | MoFEM::Interface & | m_field, |
std::string | block_name, | ||
Pip & | pip, | ||
std::string | u, | ||
std::string | ep, | ||
std::string | tau | ||
) |
|
inline |
|
inline |
\[ \begin{split} f&=\sqrt{s_{ij}s_{ij}}\\ A_{ij}&=\frac{\partial f}{\partial \sigma_{ij}}= \frac{1}{f} s_{kl} \frac{\partial s_{kl}}{\partial \sigma_{ij}}\\ \frac{\partial A_{ij}}{\partial \sigma_{kl}}&= \frac{\partial^2 f}{\partial \sigma_{ij}\partial\sigma_{mn}}= \frac{1}{f} \left( \frac{\partial s_{kl}}{\partial \sigma_{mn}}\frac{\partial s_{kl}}{\partial \sigma_{ij}} -A_{mn}A_{ij} \right)\\ \frac{\partial f}{\partial \varepsilon_{ij}}&=A_{mn}D_{mnij} \\ \frac{\partial A_{ij}}{\partial \varepsilon_{kl}}&= \frac{\partial A_{ij}}{\partial \sigma_{mn}} \frac{\partial \sigma_{mn}}{\partial \varepsilon_{kl}}= \frac{\partial A_{ij}}{\partial \sigma_{mn}} D_{mnkl} \end{split} \]
Definition at line 189 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 21 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 80 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 86 of file PlasticOpsGeneric.hpp.
|
inline |
Definition at line 276 of file PlasticOpsGeneric.hpp.
FTensor::Index<'I', 3> PlasticOps::I |
[Common data]
Definition at line 115 of file PlasticOps.hpp.
FTensor::Index<'J', 3> PlasticOps::J |
Definition at line 116 of file PlasticOps.hpp.
FTensor::Index<'M', 3> PlasticOps::M |
Definition at line 117 of file PlasticOps.hpp.
FTensor::Index<'N', 3> PlasticOps::N |
Definition at line 118 of file PlasticOps.hpp.