v0.14.0
Classes | Typedefs | Functions | Variables
PlasticOps Namespace Reference

Classes

struct  CommonData
 [Common data] More...
 
struct  Monitor
 
struct  OpCalculateConstraintsLhs_dEPImpl
 
struct  OpCalculateConstraintsLhs_dEPImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculateConstraintsLhs_dTAUImpl
 
struct  OpCalculateConstraintsLhs_dTAUImpl< GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculateConstraintsLhs_dUImpl
 
struct  OpCalculateConstraintsLhs_dUImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculateConstraintsLhs_LogStrain_dUImpl
 
struct  OpCalculateConstraintsLhs_LogStrain_dUImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculateConstraintsRhsImpl
 
struct  OpCalculateConstraintsRhsImpl< GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticFlowLhs_dEPImpl
 
struct  OpCalculatePlasticFlowLhs_dEPImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticFlowLhs_dTAUImpl
 
struct  OpCalculatePlasticFlowLhs_dTAUImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticFlowLhs_dUImpl
 
struct  OpCalculatePlasticFlowLhs_dUImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticFlowLhs_LogStrain_dUImpl
 
struct  OpCalculatePlasticFlowLhs_LogStrain_dUImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticFlowRhsImpl
 
struct  OpCalculatePlasticFlowRhsImpl< DIM, GAUSS, AssemblyDomainEleOp >
 [Calculate stress] More...
 
struct  OpCalculatePlasticInternalForceLhs_dEPImpl
 
struct  OpCalculatePlasticInternalForceLhs_dEPImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticInternalForceLhs_LogStrain_dEPImpl
 
struct  OpCalculatePlasticInternalForceLhs_LogStrain_dEPImpl< DIM, GAUSS, AssemblyDomainEleOp >
 
struct  OpCalculatePlasticityImpl
 
struct  OpCalculatePlasticityImpl< DIM, GAUSS, DomainEleOp >
 
struct  OpCalculatePlasticSurfaceImpl
 
struct  OpCalculatePlasticSurfaceImpl< DIM, GAUSS, DomainEleOp >
 [Auxiliary functions functions More...
 
struct  OpPlasticStressImpl
 
struct  OpPlasticStressImpl< DIM, GAUSS, DomainEleOp >
 
struct  PlasticityIntegrators
 

Typedefs

using Pip = boost::ptr_deque< ForcesAndSourcesCore::UserDataOperator >
 
using CommonPlasticPtr = boost::shared_ptr< PlasticOps::CommonData >
 
using CommonHenkyPtr = 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
 

Typedef Documentation

◆ CommonHenkyPtr

using PlasticOps::CommonHenkyPtr = typedef boost::shared_ptr<HenckyOps::CommonData>
Examples
PlasticOps.hpp.

Definition at line 244 of file PlasticOps.hpp.

◆ CommonPlasticPtr

using PlasticOps::CommonPlasticPtr = typedef boost::shared_ptr<PlasticOps::CommonData>
Examples
PlasticOps.hpp.

Definition at line 243 of file PlasticOps.hpp.

◆ Pip

using PlasticOps::Pip = typedef boost::ptr_deque<ForcesAndSourcesCore::UserDataOperator>
Examples
PlasticOps.hpp.

Definition at line 242 of file PlasticOps.hpp.

Function Documentation

◆ addMatBlockOps()

template<int DIM>
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

Parameters
m_field
meshset_vec_ptr
sev
Returns
MoFEMErrorCode

Get elasticity tensor

Calculate elasticity tensor for given material parameters

Parameters
mat_D_ptr
bulk_modulus_K
shear_modulus_G
Returns
MoFEMErrorCode

[Calculate elasticity tensor]

[Calculate elasticity tensor]

Examples
ADOLCPlasticityMaterialModels.hpp, PlasticOps.hpp, seepage.cpp, and thermo_elastic.cpp.

Definition at line 248 of file PlasticOps.hpp.

251  {
253 
254  struct OpMatBlocks : public DomainEleOp {
255  OpMatBlocks(boost::shared_ptr<MatrixDouble> m_D_ptr,
256  boost::shared_ptr<CommonData::BlockParams> mat_params_ptr,
257  double scale_value, MoFEM::Interface &m_field, Sev sev,
258  std::vector<const CubitMeshSets *> meshset_vec_ptr)
259  : DomainEleOp(NOSPACE, DomainEleOp::OPSPACE), matDPtr(m_D_ptr),
260  matParamsPtr(mat_params_ptr), scaleVal(scale_value) {
261  CHK_THROW_MESSAGE(extractBlockData(m_field, meshset_vec_ptr, sev),
262  "Can not get data from block");
263  }
264 
265  MoFEMErrorCode doWork(int side, EntityType type,
268 
269  auto getK = [](auto &p) {
270  auto young_modulus = p[CommonData::YOUNG_MODULUS];
271  auto poisson_ratio = p[CommonData::POISSON_RATIO];
272  return young_modulus / (3 * (1 - 2 * poisson_ratio));
273  };
274 
275  auto getG = [](auto &p) {
276  auto young_modulus = p[CommonData::YOUNG_MODULUS];
277  auto poisson_ratio = p[CommonData::POISSON_RATIO];
278  return young_modulus / (2 * (1 + poisson_ratio));
279  };
280 
281  auto scale_fun = [this](auto &p) {
282  p[CommonData::YOUNG_MODULUS] *= scaleVal;
283  p[CommonData::SIGMA_Y] *= scaleVal;
284  p[CommonData::H] *= scaleVal;
285  p[CommonData::VIS_H] *= scaleVal;
286  p[CommonData::QINF] *= scaleVal;
287  p[CommonData::C1_k] *= scaleVal;
288  };
289 
290  for (auto &b : blockData) {
291  if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
292  *matParamsPtr = b.bParams;
293  scale_fun(*matParamsPtr);
294  CHKERR getMatDPtr(matDPtr, getK(*matParamsPtr), getG(*matParamsPtr));
296  }
297  }
298 
299  (*matParamsPtr) = {young_modulus, poisson_ratio, sigmaY, H,
300  visH, Qinf, b_iso, C1_k};
301  scale_fun(*matParamsPtr);
302  CHKERR getMatDPtr(matDPtr, getK(*matParamsPtr), getG(*matParamsPtr));
303 
305  }
306 
307  private:
308  boost::shared_ptr<MatrixDouble> matDPtr;
309  boost::shared_ptr<CommonData::BlockParams> matParamsPtr;
310  const double scaleVal;
311 
312  struct BlockData {
313  std::array<double, CommonData::LAST_PARAM> bParams;
314  Range blockEnts;
315  };
316  std::vector<BlockData> blockData;
317 
318  /**
319  * @brief Extract block data from meshsets
320  *
321  * @param m_field
322  * @param meshset_vec_ptr
323  * @param sev
324  * @return MoFEMErrorCode
325  */
327  extractBlockData(MoFEM::Interface &m_field,
328  std::vector<const CubitMeshSets *> meshset_vec_ptr,
329  Sev sev) {
331 
332  for (auto m : meshset_vec_ptr) {
333  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock") << *m;
334  std::vector<double> block_data;
335  CHKERR m->getAttributes(block_data);
336  if (block_data.size() != CommonData::LAST_PARAM) {
337  SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
338  "Wrong number of block attribute");
339  }
340  auto get_block_ents = [&]() {
341  Range ents;
342  CHKERR m_field.get_moab().get_entities_by_handle(m->meshset, ents,
343  true);
344  return ents;
345  };
346 
347  CommonData::BlockParams block_params;
348  for (auto i = 0; i != CommonData::LAST_PARAM; ++i) {
349  block_params[i] = block_data[i];
350  }
351 
352  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
353  << "E = " << block_params[CommonData::YOUNG_MODULUS]
354  << " nu = " << block_params[CommonData::POISSON_RATIO];
355  MOFEM_TAG_AND_LOG("WORLD", sev, "MatBlock")
356  << std::endl
357  << "sigma_y = " << block_params[CommonData::SIGMA_Y] << std::endl
358  << "h = " << block_params[CommonData::H] << std::endl
359  << "vis_h = " << block_params[CommonData::VIS_H] << std::endl
360  << "qinf = " << block_params[CommonData::QINF] << std::endl
361  << "biso = " << block_params[CommonData::BISO] << std::endl
362  << "C1_k = " << block_params[CommonData::C1_k] << std::endl;
363 
364  blockData.push_back({block_params, get_block_ents()});
365  }
366  MOFEM_LOG_CHANNEL("WORLD");
368  }
369 
370  /**
371  * @brief Get elasticity tensor
372  *
373  * Calculate elasticity tensor for given material parameters
374  *
375  * @param mat_D_ptr
376  * @param bulk_modulus_K
377  * @param shear_modulus_G
378  * @return MoFEMErrorCode
379  *
380  */
381  MoFEMErrorCode getMatDPtr(boost::shared_ptr<MatrixDouble> mat_D_ptr,
382  double bulk_modulus_K, double shear_modulus_G) {
384  //! [Calculate elasticity tensor]
385  auto set_material_stiffness = [&]() {
391  double A = (DIM == 2)
392  ? 2 * shear_modulus_G /
393  (bulk_modulus_K + (4. / 3.) * shear_modulus_G)
394  : 1;
395  auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*mat_D_ptr);
396  t_D(i, j, k, l) =
397  2 * shear_modulus_G * ((t_kd(i, k) ^ t_kd(j, l)) / 4.) +
398  A * (bulk_modulus_K - (2. / 3.) * shear_modulus_G) * t_kd(i, j) *
399  t_kd(k, l);
400  };
401  //! [Calculate elasticity tensor]
402  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
403  mat_D_ptr->resize(size_symm * size_symm, 1);
404  set_material_stiffness();
406  }
407  };
408 
409  // push operator to calculate material stiffness matrix for each block
410  pip.push_back(new OpMatBlocks(
411  mat_D_Ptr, mat_params_ptr, scale_value, m_field, sev,
412 
413  // Get blockset using regular expression
414  m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(std::regex(
415 
416  (boost::format("%s(.*)") % block_name).str()
417 
418  ))
419 
420  ));
421 
423 }

◆ constrain_abs()

double PlasticOps::constrain_abs ( double  x)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 265 of file PlasticOpsGeneric.hpp.

265  {
266  const auto y = -x / zeta;
267  if (y > std::numeric_limits<float>::max_exponent10 ||
268  y < std::numeric_limits<float>::min_exponent10) {
269  return std::abs(x);
270  } else {
271  const double e = std::exp(y);
272  return x + 2 * zeta * std::log1p(e);
273  }
274 };

◆ constraint()

double PlasticOps::constraint ( double  eqiv,
double  dot_tau,
double  f,
double  sigma_y,
double  abs_w,
double  vis_H,
double  sigma_Y 
)
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 \]

Examples
PlasticOpsGeneric.hpp.

Definition at line 292 of file PlasticOpsGeneric.hpp.

293  {
294  return vis_H * dot_tau +
295  (sigma_Y / 2) * ((cn0 * (dot_tau - eqiv) + cn1 * (dot_tau) -
296  (f - sigma_y) / sigma_Y) -
297  abs_w);
298 };

◆ constrian_sign()

double PlasticOps::constrian_sign ( double  x)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 251 of file PlasticOpsGeneric.hpp.

251  {
252  const auto y = x / zeta;
253  if (y > std::numeric_limits<float>::max_exponent10 ||
254  y < std::numeric_limits<float>::min_exponent10) {
255  if (x > 0)
256  return 1.;
257  else
258  return -1.;
259  } else {
260  const auto e = std::exp(y);
261  return (e - 1) / (1 + e);
262  }
263 };

◆ createCommonPlasticOps()

template<int DIM, IntegrationType I, typename DomainEleOp >
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 
)
Examples
PlasticOps.hpp.

Definition at line 426 of file PlasticOps.hpp.

429  {
430 
431  using P = PlasticityIntegrators<DomainEleOp>;
432 
433  auto common_plastic_ptr = boost::make_shared<PlasticOps::CommonData>();
434  common_plastic_ptr = boost::make_shared<PlasticOps::CommonData>();
435 
436  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
437  auto make_d_mat = []() {
438  return boost::make_shared<MatrixDouble>(size_symm * size_symm, 1);
439  };
440 
441  common_plastic_ptr->mDPtr = make_d_mat();
442  common_plastic_ptr->mGradPtr = boost::make_shared<MatrixDouble>();
443  common_plastic_ptr->mStrainPtr = boost::make_shared<MatrixDouble>();
444  common_plastic_ptr->mStressPtr = boost::make_shared<MatrixDouble>();
445 
446  auto m_D_ptr = common_plastic_ptr->mDPtr;
447 
448  CHK_THROW_MESSAGE(addMatBlockOps<DIM>(m_field, block_name, pip, m_D_ptr,
449  common_plastic_ptr->getParamsPtr(),
450  scale, sev),
451  "add mat block ops");
452 
453  pip.push_back(new OpCalculateScalarFieldValues(
454  tau, common_plastic_ptr->getPlasticTauPtr()));
455  pip.push_back(new OpCalculateTensor2SymmetricFieldValues<DIM>(
456  ep, common_plastic_ptr->getPlasticStrainPtr()));
458  u, common_plastic_ptr->mGradPtr));
459 
460  CommonHenkyPtr common_henky_ptr;
461 
462  if (is_large_strains) {
463  common_henky_ptr = boost::make_shared<HenckyOps::CommonData>();
464  common_henky_ptr->matGradPtr = common_plastic_ptr->mGradPtr;
465  common_henky_ptr->matDPtr = common_plastic_ptr->mDPtr;
466  common_henky_ptr->matLogCPlastic =
467  common_plastic_ptr->getPlasticStrainPtr();
468  common_plastic_ptr->mStrainPtr = common_henky_ptr->getMatLogC();
469  common_plastic_ptr->mStressPtr = common_henky_ptr->getMatHenckyStress();
470 
472 
473  pip.push_back(new typename H::template OpCalculateEigenVals<DIM, I>(
474  u, common_henky_ptr));
475  pip.push_back(
476  new typename H::template OpCalculateLogC<DIM, I>(u, common_henky_ptr));
477  pip.push_back(new typename H::template OpCalculateLogC_dC<DIM, I>(
478  u, common_henky_ptr));
479  pip.push_back(
480  new typename H::template OpCalculateHenckyPlasticStress<DIM, I, 0>(
481  u, common_henky_ptr, m_D_ptr));
482  pip.push_back(new typename H::template OpCalculatePiolaStress<DIM, I, 0>(
483  u, common_henky_ptr));
484  } else {
485 
486  pip.push_back(new OpSymmetrizeTensor<SPACE_DIM>(
487  common_plastic_ptr->mGradPtr, common_plastic_ptr->mStrainPtr));
488  pip.push_back(new typename P::template OpPlasticStress<DIM, I>(
489  common_plastic_ptr, m_D_ptr));
490  }
491 
492  pip.push_back(new typename P::template OpCalculatePlasticSurface<DIM, I>(
493  u, common_plastic_ptr));
494 
495  return std::make_tuple(common_plastic_ptr, common_henky_ptr);
496 }

◆ deviator() [1/3]

template<typename T >
auto PlasticOps::deviator ( FTensor::Tensor2_symmetric< T, 2 > &  t_stress,
double  trace,
FTensor::Tensor2_symmetric< double, 2 > &&  t_alpha 
)
inline

Definition at line 113 of file PlasticOpsGeneric.hpp.

114  {
115  return deviator(t_stress, trace, t_alpha, FTensor::Number<2>());
116 };

◆ deviator() [2/3]

template<typename T >
auto PlasticOps::deviator ( FTensor::Tensor2_symmetric< T, 3 > &  t_stress,
double  trace,
FTensor::Tensor2_symmetric< double, 3 > &&  t_alpha 
)
inline

Definition at line 119 of file PlasticOpsGeneric.hpp.

120  {
121  return deviator(t_stress, trace, t_alpha, FTensor::Number<3>());
122 };

◆ deviator() [3/3]

template<typename T , int DIM>
auto PlasticOps::deviator ( FTensor::Tensor2_symmetric< T, DIM > &  t_stress,
double  trace,
FTensor::Tensor2_symmetric< double, DIM > &  t_alpha,
FTensor::Number< DIM >   
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 92 of file PlasticOpsGeneric.hpp.

97  {
99  t_dev(I, J) = 0;
100  for (int ii = 0; ii != DIM; ++ii)
101  for (int jj = ii; jj != DIM; ++jj)
102  t_dev(ii, jj) = t_stress(ii, jj);
103  t_dev(0, 0) -= trace;
104  t_dev(1, 1) -= trace;
105  t_dev(2, 2) -= trace;
106  for (int ii = 0; ii != DIM; ++ii)
107  for (int jj = ii; jj != DIM; ++jj)
108  t_dev(ii, jj) -= t_alpha(ii, jj);
109  return t_dev;
110 };

◆ diff_constrain_ddot_tau()

double PlasticOps::diff_constrain_ddot_tau ( double  sign,
double  eqiv,
double  dot_tau,
double  vis_H,
double  sigma_Y 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 300 of file PlasticOpsGeneric.hpp.

301  {
302  return vis_H + (sigma_Y / 2) * (cn0 + cn1 * (1 - sign));
303 };

◆ diff_constrain_deqiv()

double PlasticOps::diff_constrain_deqiv ( double  sign,
double  eqiv,
double  dot_tau,
double  sigma_Y 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 305 of file PlasticOpsGeneric.hpp.

306  {
307  return (sigma_Y / 2) * (-cn0);
308 };

◆ diff_constrain_df()

auto PlasticOps::diff_constrain_df ( double  sign)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 310 of file PlasticOpsGeneric.hpp.

310 { return (-1 - sign) / 2; };

◆ diff_constrain_dsigma_y()

auto PlasticOps::diff_constrain_dsigma_y ( double  sign)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 312 of file PlasticOpsGeneric.hpp.

312 { return (1 + sign) / 2; }

◆ diff_constrain_dstrain()

template<typename T1 , typename T2 , int DIM>
auto PlasticOps::diff_constrain_dstrain ( T1 &  t_D,
T2 &&  t_diff_constrain_dstress 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 326 of file PlasticOpsGeneric.hpp.

326  {
331  FTensor::Tensor2_symmetric<double, DIM> t_diff_constrain_dstrain;
332  t_diff_constrain_dstrain(k, l) =
333  t_diff_constrain_dstress(i, j) * t_D(i, j, k, l);
334  return t_diff_constrain_dstrain;
335 };

◆ diff_constrain_dstress()

template<typename T , int DIM>
auto PlasticOps::diff_constrain_dstress ( double  diff_constrain_df,
FTensor::Tensor2_symmetric< T, DIM > &  t_plastic_flow 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 316 of file PlasticOpsGeneric.hpp.

317  {
320  FTensor::Tensor2_symmetric<double, DIM> t_diff_constrain_dstress;
321  t_diff_constrain_dstress(i, j) = diff_constrain_df * t_plastic_flow(i, j);
322  return t_diff_constrain_dstress;
323 };

◆ diff_deviator() [1/3]

auto PlasticOps::diff_deviator ( FTensor::Ddg< double, 2, 2 > &&  t_diff_stress)
inline

Definition at line 157 of file PlasticOpsGeneric.hpp.

157  {
158  return diff_deviator(std::move(t_diff_stress), FTensor::Number<2>());
159 }

◆ diff_deviator() [2/3]

auto PlasticOps::diff_deviator ( FTensor::Ddg< double, 3, 3 > &&  t_diff_stress)
inline

Definition at line 161 of file PlasticOpsGeneric.hpp.

161  {
162  return diff_deviator(std::move(t_diff_stress), FTensor::Number<3>());
163 }

◆ diff_deviator() [3/3]

template<int DIM>
auto PlasticOps::diff_deviator ( FTensor::Ddg< double, DIM, DIM > &&  t_diff_stress,
FTensor::Number< DIM >   
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 125 of file PlasticOpsGeneric.hpp.

126  {
129  FTensor::Ddg<double, 3, DIM> t_diff_deviator;
130  t_diff_deviator(I, J, k, l) = 0;
131  for (int ii = 0; ii != DIM; ++ii)
132  for (int jj = ii; jj != DIM; ++jj)
133  for (int kk = 0; kk != DIM; ++kk)
134  for (int ll = kk; ll != DIM; ++ll)
135  t_diff_deviator(ii, jj, kk, ll) = t_diff_stress(ii, jj, kk, ll);
136 
137  constexpr double third = boost::math::constants::third<double>();
138 
139  t_diff_deviator(0, 0, 0, 0) -= third;
140  t_diff_deviator(0, 0, 1, 1) -= third;
141 
142  t_diff_deviator(1, 1, 0, 0) -= third;
143  t_diff_deviator(1, 1, 1, 1) -= third;
144 
145  t_diff_deviator(2, 2, 0, 0) -= third;
146  t_diff_deviator(2, 2, 1, 1) -= third;
147 
148  if constexpr (DIM == 3) {
149  t_diff_deviator(0, 0, 2, 2) -= third;
150  t_diff_deviator(1, 1, 2, 2) -= third;
151  t_diff_deviator(2, 2, 2, 2) -= third;
152  }
153 
154  return t_diff_deviator;
155 };

◆ diff_equivalent_strain_dot()

template<typename T1 , typename T2 , typename T3 , int DIM>
auto PlasticOps::diff_equivalent_strain_dot ( const T1  eqiv,
T2 &  t_plastic_strain_dot,
T3 &  t_diff_plastic_strain,
FTensor::Number< DIM >   
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 349 of file PlasticOpsGeneric.hpp.

351  {
356  constexpr double A = 2. / 3;
358  t_diff_eqiv(i, j) = A * (t_plastic_strain_dot(k, l) / eqiv) *
359  t_diff_plastic_strain(k, l, i, j);
360  return t_diff_eqiv;
361 };

◆ diff_plastic_flow_dstrain()

template<typename T , int DIM>
auto PlasticOps::diff_plastic_flow_dstrain ( FTensor::Ddg< T, DIM, DIM > &  t_D,
FTensor::Ddg< double, DIM, DIM > &&  t_diff_plastic_flow_dstress 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 224 of file PlasticOpsGeneric.hpp.

226  {
233  FTensor::Ddg<double, DIM, DIM> t_diff_flow;
234  t_diff_flow(i, j, k, l) =
235  t_diff_plastic_flow_dstress(i, j, m, n) * t_D(m, n, k, l);
236  return t_diff_flow;
237 };

◆ diff_plastic_flow_dstress()

template<typename T , int DIM>
auto PlasticOps::diff_plastic_flow_dstress ( long double  f,
FTensor::Tensor2_symmetric< T, DIM > &  t_flow,
FTensor::Ddg< double, 3, DIM > &&  t_diff_deviator 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 208 of file PlasticOpsGeneric.hpp.

210  {
215  FTensor::Ddg<double, DIM, DIM> t_diff_flow;
216  t_diff_flow(i, j, k, l) =
217  (1.5 * (t_diff_deviator(M, N, i, j) * t_diff_deviator(M, N, k, l) -
218  (2. / 3.) * t_flow(i, j) * t_flow(k, l))) /
219  f;
220  return t_diff_flow;
221 };

◆ diff_symmetrize()

template<int DIM>
auto PlasticOps::diff_symmetrize ( FTensor::Number< DIM >  )
inline
Examples
PlasticOpsGeneric.hpp, PlasticOpsLargeStrains.hpp, and PlasticOpsSmallStrains.hpp.

Definition at line 43 of file PlasticOpsGeneric.hpp.

43  {
44 
49 
51 
52  t_diff(i, j, k, l) = 0;
53  t_diff(0, 0, 0, 0) = 1;
54  t_diff(1, 1, 1, 1) = 1;
55 
56  t_diff(1, 0, 1, 0) = 0.5;
57  t_diff(1, 0, 0, 1) = 0.5;
58 
59  t_diff(0, 1, 0, 1) = 0.5;
60  t_diff(0, 1, 1, 0) = 0.5;
61 
62  if constexpr (DIM == 3) {
63  t_diff(2, 2, 2, 2) = 1;
64 
65  t_diff(2, 0, 2, 0) = 0.5;
66  t_diff(2, 0, 0, 2) = 0.5;
67  t_diff(0, 2, 0, 2) = 0.5;
68  t_diff(0, 2, 2, 0) = 0.5;
69 
70  t_diff(2, 1, 2, 1) = 0.5;
71  t_diff(2, 1, 1, 2) = 0.5;
72  t_diff(1, 2, 1, 2) = 0.5;
73  t_diff(1, 2, 2, 1) = 0.5;
74  }
75 
76  return t_diff;
77 };

◆ diff_tensor()

template<int DIM>
auto PlasticOps::diff_tensor ( FTensor::Number< DIM >  )
inline

[Lambda functions]

Examples
PlasticOpsGeneric.hpp.

Definition at line 10 of file PlasticOpsGeneric.hpp.

10  {
17  t_diff(i, j, k, l) = (t_kd(i, k) ^ t_kd(j, l)) / 4.;
18  return t_diff;
19 };

◆ equivalent_strain_dot()

template<typename T , int DIM>
auto PlasticOps::equivalent_strain_dot ( FTensor::Tensor2_symmetric< T, DIM > &  t_plastic_strain_dot)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 338 of file PlasticOpsGeneric.hpp.

339  {
342  constexpr double A = 2. / 3;
343  return std::sqrt(A * t_plastic_strain_dot(i, j) *
344  t_plastic_strain_dot(i, j)) +
345  std::numeric_limits<double>::epsilon();
346 };

◆ get_mat_scalar_dtensor_sym() [1/2]

auto PlasticOps::get_mat_scalar_dtensor_sym ( MatrixDouble &  mat,
FTensor::Number< 2 >   
)
Examples
PlasticOpsGeneric.hpp.

Definition at line 1184 of file PlasticOpsGeneric.hpp.

1184  {
1186  &mat(0, 0), &mat(0, 1), &mat(0, 2)};
1187 }

◆ get_mat_scalar_dtensor_sym() [2/2]

auto PlasticOps::get_mat_scalar_dtensor_sym ( MatrixDouble &  mat,
FTensor::Number< 3 >   
)

Definition at line 1189 of file PlasticOpsGeneric.hpp.

1189  {
1191  &mat(0, 0), &mat(0, 1), &mat(0, 2), &mat(0, 3), &mat(0, 4), &mat(0, 5)};
1192 }

◆ get_mat_tensor_sym_dscalar() [1/2]

static auto PlasticOps::get_mat_tensor_sym_dscalar ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 2 >   
)
inlinestatic
Examples
PlasticOpsGeneric.hpp.

Definition at line 1090 of file PlasticOpsGeneric.hpp.

1091  {
1093  &mat(3 * rr + 0, 0), &mat(3 * rr + 1, 0), &mat(3 * rr + 2, 0)};
1094 }

◆ get_mat_tensor_sym_dscalar() [2/2]

static auto PlasticOps::get_mat_tensor_sym_dscalar ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 3 >   
)
inlinestatic

Definition at line 1096 of file PlasticOpsGeneric.hpp.

1097  {
1099  &mat(6 * rr + 0, 0), &mat(6 * rr + 1, 0), &mat(6 * rr + 2, 0),
1100  &mat(6 * rr + 3, 0), &mat(6 * rr + 4, 0), &mat(6 * rr + 5, 0)};
1101 }

◆ get_mat_tensor_sym_dtensor_sym() [1/2]

static auto PlasticOps::get_mat_tensor_sym_dtensor_sym ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 2 >   
)
inlinestatic
Examples
PlasticOpsGeneric.hpp.

Definition at line 958 of file PlasticOpsGeneric.hpp.

959  {
961  &mat(3 * rr + 0, 0), &mat(3 * rr + 0, 1), &mat(3 * rr + 0, 2),
962  &mat(3 * rr + 1, 0), &mat(3 * rr + 1, 1), &mat(3 * rr + 1, 2),
963  &mat(3 * rr + 2, 0), &mat(3 * rr + 2, 1), &mat(3 * rr + 2, 2)};
964 }

◆ get_mat_tensor_sym_dtensor_sym() [2/2]

static auto PlasticOps::get_mat_tensor_sym_dtensor_sym ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 3 >   
)
inlinestatic

Definition at line 966 of file PlasticOpsGeneric.hpp.

967  {
969  &mat(6 * rr + 0, 0), &mat(6 * rr + 0, 1), &mat(6 * rr + 0, 2),
970  &mat(6 * rr + 0, 3), &mat(6 * rr + 0, 4), &mat(6 * rr + 0, 5),
971  &mat(6 * rr + 1, 0), &mat(6 * rr + 1, 1), &mat(6 * rr + 1, 2),
972  &mat(6 * rr + 1, 3), &mat(6 * rr + 1, 4), &mat(6 * rr + 1, 5),
973  &mat(6 * rr + 2, 0), &mat(6 * rr + 2, 1), &mat(6 * rr + 2, 2),
974  &mat(6 * rr + 2, 3), &mat(6 * rr + 2, 4), &mat(6 * rr + 2, 5),
975  &mat(6 * rr + 3, 0), &mat(6 * rr + 3, 1), &mat(6 * rr + 3, 2),
976  &mat(6 * rr + 3, 3), &mat(6 * rr + 3, 4), &mat(6 * rr + 3, 5),
977  &mat(6 * rr + 4, 0), &mat(6 * rr + 4, 1), &mat(6 * rr + 4, 2),
978  &mat(6 * rr + 4, 3), &mat(6 * rr + 4, 4), &mat(6 * rr + 4, 5),
979  &mat(6 * rr + 5, 0), &mat(6 * rr + 5, 1), &mat(6 * rr + 5, 2),
980  &mat(6 * rr + 5, 3), &mat(6 * rr + 5, 4), &mat(6 * rr + 5, 5)};
981 }

◆ get_mat_tensor_sym_dvector() [1/2]

static auto PlasticOps::get_mat_tensor_sym_dvector ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 2 >   
)
inlinestatic

[Lambda functions]

[Auxiliary functions functions

Examples
PlasticOpsGeneric.hpp, PlasticOpsLargeStrains.hpp, and PlasticOpsSmallStrains.hpp.

Definition at line 366 of file PlasticOpsGeneric.hpp.

367  {
369  &mat(3 * rr + 0, 0), &mat(3 * rr + 0, 1), &mat(3 * rr + 1, 0),
370  &mat(3 * rr + 1, 1), &mat(3 * rr + 2, 0), &mat(3 * rr + 2, 1)};
371 }

◆ get_mat_tensor_sym_dvector() [2/2]

static auto PlasticOps::get_mat_tensor_sym_dvector ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 3 >   
)
inlinestatic

Definition at line 373 of file PlasticOpsGeneric.hpp.

374  {
376  &mat(6 * rr + 0, 0), &mat(6 * rr + 0, 1), &mat(6 * rr + 0, 2),
377  &mat(6 * rr + 1, 0), &mat(6 * rr + 1, 1), &mat(6 * rr + 1, 2),
378  &mat(6 * rr + 2, 0), &mat(6 * rr + 2, 1), &mat(6 * rr + 2, 2),
379  &mat(6 * rr + 3, 0), &mat(6 * rr + 3, 1), &mat(6 * rr + 3, 2),
380  &mat(6 * rr + 4, 0), &mat(6 * rr + 4, 1), &mat(6 * rr + 4, 2),
381  &mat(6 * rr + 5, 0), &mat(6 * rr + 5, 1), &mat(6 * rr + 5, 2)};
382 }

◆ get_mat_vector_dtensor_sym() [1/2]

static FTensor::Tensor2<FTensor::PackPtr<double *, 3>, 2, 3> PlasticOps::get_mat_vector_dtensor_sym ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 2 >   
)
inlinestatic
Examples
PlasticOpsLargeStrains.hpp, and PlasticOpsSmallStrains.hpp.

Definition at line 38 of file PlasticOpsSmallStrains.hpp.

38  {
40  &mat(2 * rr + 0, 0), &mat(2 * rr + 0, 1), &mat(2 * rr + 0, 2),
41  &mat(2 * rr + 1, 0), &mat(2 * rr + 1, 1), &mat(2 * rr + 1, 2)};
42 }

◆ get_mat_vector_dtensor_sym() [2/2]

static FTensor::Tensor2<FTensor::PackPtr<double *, 6>, 3, 6> PlasticOps::get_mat_vector_dtensor_sym ( size_t  rr,
MatrixDouble &  mat,
FTensor::Number< 3 >   
)
inlinestatic

Definition at line 45 of file PlasticOpsSmallStrains.hpp.

45  {
47  &mat(3 * rr + 0, 0), &mat(3 * rr + 0, 1), &mat(3 * rr + 0, 2),
48  &mat(3 * rr + 0, 3), &mat(3 * rr + 0, 4), &mat(3 * rr + 0, 5),
49 
50  &mat(3 * rr + 1, 0), &mat(3 * rr + 1, 1), &mat(3 * rr + 1, 2),
51  &mat(3 * rr + 1, 3), &mat(3 * rr + 1, 4), &mat(3 * rr + 1, 5),
52 
53  &mat(3 * rr + 2, 0), &mat(3 * rr + 2, 1), &mat(3 * rr + 2, 2),
54  &mat(3 * rr + 2, 3), &mat(3 * rr + 2, 4), &mat(3 * rr + 2, 5)};
55 }

◆ opFactoryDomainLhs()

template<int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp >
MoFEMErrorCode PlasticOps::opFactoryDomainLhs ( MoFEM::Interface m_field,
std::string  block_name,
Pip pip,
std::string  u,
std::string  ep,
std::string  tau 
)
Examples
PlasticOps.hpp.

Definition at line 547 of file PlasticOps.hpp.

548  {
550 
551  using namespace HenckyOps;
552 
553  using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
554  A>::template BiLinearForm<I>;
555  using OpKPiola = typename B::template OpGradTensorGrad<1, DIM, DIM, 1>;
556  using OpKCauchy = typename B::template OpGradSymTensorGrad<1, DIM, DIM, 0>;
557 
558  using P = PlasticityIntegrators<DomainEleOp>;
559 
560  auto [common_plastic_ptr, common_henky_ptr] =
561  createCommonPlasticOps<DIM, I, DomainEleOp>(m_field, block_name, pip, u,
562  ep, tau, scale, Sev::verbose);
563 
564  auto m_D_ptr = common_plastic_ptr->mDPtr;
565 
566  pip.push_back(new OpCalculateTensor2SymmetricFieldValuesDot<DIM>(
567  ep, common_plastic_ptr->getPlasticStrainDotPtr()));
568  pip.push_back(new OpCalculateScalarFieldValuesDot(
569  tau, common_plastic_ptr->getPlasticTauDotPtr()));
570  pip.push_back(new typename P::template OpCalculatePlasticity<DIM, I>(
571  u, common_plastic_ptr, m_D_ptr));
572 
573  if (common_henky_ptr) {
575  pip.push_back(new typename H::template OpHenckyTangent<DIM, I, 0>(
576  u, common_henky_ptr, m_D_ptr));
577  pip.push_back(new OpKPiola(u, u, common_henky_ptr->getMatTangent()));
578  pip.push_back(
579  new typename P::template Assembly<A>::
580  template OpCalculatePlasticInternalForceLhs_LogStrain_dEP<DIM, I>(
581  u, ep, common_plastic_ptr, common_henky_ptr, m_D_ptr));
582  } else {
583  pip.push_back(new OpKCauchy(u, u, m_D_ptr));
584  pip.push_back(new typename P::template Assembly<A>::
585  template OpCalculatePlasticInternalForceLhs_dEP<DIM, I>(
586  u, ep, common_plastic_ptr, m_D_ptr));
587  }
588 
589  if (common_henky_ptr) {
590  pip.push_back(
591  new typename P::template Assembly<A>::
592  template OpCalculateConstraintsLhs_LogStrain_dU<DIM, I>(
593  tau, u, common_plastic_ptr, common_henky_ptr, m_D_ptr));
594  pip.push_back(
595  new typename P::template Assembly<A>::
596  template OpCalculatePlasticFlowLhs_LogStrain_dU<DIM, I>(
597  ep, u, common_plastic_ptr, common_henky_ptr, m_D_ptr));
598  } else {
599  pip.push_back(
600  new
601  typename P::template Assembly<A>::template OpCalculateConstraintsLhs_dU<
602  DIM, I>(tau, u, common_plastic_ptr, m_D_ptr));
603  pip.push_back(
604  new
605  typename P::template Assembly<A>::template OpCalculatePlasticFlowLhs_dU<
606  DIM, I>(ep, u, common_plastic_ptr, m_D_ptr));
607  }
608 
609  pip.push_back(
610  new
611  typename P::template Assembly<A>::template OpCalculatePlasticFlowLhs_dEP<
612  DIM, I>(ep, ep, common_plastic_ptr, m_D_ptr));
613  pip.push_back(
614  new
615  typename P::template Assembly<A>::template OpCalculatePlasticFlowLhs_dTAU<
616  DIM, I>(ep, tau, common_plastic_ptr, m_D_ptr));
617  pip.push_back(
618  new
619  typename P::template Assembly<A>::template OpCalculateConstraintsLhs_dEP<
620  DIM, I>(tau, ep, common_plastic_ptr, m_D_ptr));
621  pip.push_back(
622  new
623  typename P::template Assembly<A>::template OpCalculateConstraintsLhs_dTAU<
624  I>(tau, tau, common_plastic_ptr));
625 
627 }

◆ opFactoryDomainReactions()

template<int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp >
MoFEMErrorCode PlasticOps::opFactoryDomainReactions ( MoFEM::Interface m_field,
std::string  block_name,
Pip pip,
std::string  u,
std::string  ep,
std::string  tau 
)
Examples
PlasticOps.hpp.

Definition at line 630 of file PlasticOps.hpp.

633  {
635 
636  using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
637  A>::template LinearForm<I>;
638  using OpInternalForceCauchy =
639  typename B::template OpGradTimesSymTensor<1, DIM, DIM>;
640  using OpInternalForcePiola =
641  typename B::template OpGradTimesTensor<1, DIM, DIM>;
642 
643  auto [common_plastic_ptr, common_henky_ptr] =
644  createCommonPlasticOps<DIM, I, DomainEleOp>(m_field, block_name, pip, u,
645  ep, tau, 1., Sev::inform);
646 
647  // Calculate internal forces
648  if (common_henky_ptr) {
649  pip.push_back(new OpInternalForcePiola(
650  u, common_henky_ptr->getMatFirstPiolaStress()));
651  } else {
652  pip.push_back(new OpInternalForceCauchy(u, common_plastic_ptr->mStressPtr));
653  }
654 
656 }

◆ opFactoryDomainRhs()

template<int DIM, AssemblyType A, IntegrationType I, typename DomainEleOp >
MoFEMErrorCode PlasticOps::opFactoryDomainRhs ( MoFEM::Interface m_field,
std::string  block_name,
Pip pip,
std::string  u,
std::string  ep,
std::string  tau 
)
Examples
PlasticOps.hpp.

Definition at line 500 of file PlasticOps.hpp.

501  {
503 
504  using B = typename FormsIntegrators<DomainEleOp>::template Assembly<
505  A>::template LinearForm<I>;
506  using OpInternalForceCauchy =
507  typename B::template OpGradTimesSymTensor<1, DIM, DIM>;
508  using OpInternalForcePiola =
509  typename B::template OpGradTimesTensor<1, DIM, DIM>;
510 
511  using P = PlasticityIntegrators<DomainEleOp>;
512 
513  auto [common_plastic_ptr, common_henky_ptr] =
514  createCommonPlasticOps<DIM, I, DomainEleOp>(m_field, block_name, pip, u,
515  ep, tau, scale, Sev::inform);
516 
517  auto m_D_ptr = common_plastic_ptr->mDPtr;
518 
519  pip.push_back(new OpCalculateTensor2SymmetricFieldValuesDot<DIM>(
520  ep, common_plastic_ptr->getPlasticStrainDotPtr()));
521  pip.push_back(new OpCalculateScalarFieldValuesDot(
522  tau, common_plastic_ptr->getPlasticTauDotPtr()));
523  pip.push_back(new typename P::template OpCalculatePlasticity<DIM, I>(
524  u, common_plastic_ptr, m_D_ptr));
525 
526  // Calculate internal forces
527  if (common_henky_ptr) {
528  pip.push_back(new OpInternalForcePiola(
529  u, common_henky_ptr->getMatFirstPiolaStress()));
530  } else {
531  pip.push_back(new OpInternalForceCauchy(u, common_plastic_ptr->mStressPtr));
532  }
533 
534  pip.push_back(
535  new
536  typename P::template Assembly<A>::template OpCalculateConstraintsRhs<I>(
537  tau, common_plastic_ptr, m_D_ptr));
538  pip.push_back(
539  new typename P::template Assembly<A>::template OpCalculatePlasticFlowRhs<
540  DIM, I>(ep, common_plastic_ptr, m_D_ptr));
541 
543 }

◆ plastic_flow()

template<int DIM>
auto PlasticOps::plastic_flow ( long double  f,
FTensor::Tensor2_symmetric< double, 3 > &&  t_dev_stress,
FTensor::Ddg< double, 3, DIM > &&  t_diff_deviator 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 195 of file PlasticOpsGeneric.hpp.

197  {
201  t_diff_f(k, l) =
202  (1.5 * (t_dev_stress(I, J) * t_diff_deviator(I, J, k, l))) / f;
203  return t_diff_f;
204 };

◆ platsic_surface()

double PlasticOps::platsic_surface ( FTensor::Tensor2_symmetric< double, 3 > &&  t_stress_deviator)
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} \]

Examples
PlasticOpsGeneric.hpp.

Definition at line 189 of file PlasticOpsGeneric.hpp.

189  {
190  return std::sqrt(1.5 * t_stress_deviator(I, J) * t_stress_deviator(I, J)) +
191  std::numeric_limits<double>::epsilon();
192 };

◆ symm_L_tensor()

template<int DIM>
auto PlasticOps::symm_L_tensor ( FTensor::Number< DIM >  )
inline
Examples
EshelbianOperators.cpp, PlasticOpsGeneric.hpp, PlasticOpsLargeStrains.hpp, and PlasticOpsSmallStrains.hpp.

Definition at line 21 of file PlasticOpsGeneric.hpp.

21  {
24  constexpr auto size_symm = (DIM * (DIM + 1)) / 2;
27  t_L(i, j, L) = 0;
28  if constexpr (DIM == 2) {
29  t_L(0, 0, 0) = 1;
30  t_L(1, 0, 1) = 1;
31  t_L(1, 1, 2) = 1;
32  } else {
33  t_L(0, 0, 0) = 1;
34  t_L(1, 0, 1) = 1;
35  t_L(2, 0, 2) = 1;
36  t_L(1, 1, 3) = 1;
37  t_L(2, 1, 4) = 1;
38  t_L(2, 2, 5) = 1;
39  }
40  return t_L;
41 }

◆ trace() [1/2]

template<typename T >
double PlasticOps::trace ( FTensor::Tensor2_symmetric< T, 2 > &  t_stress)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 80 of file PlasticOpsGeneric.hpp.

80  {
81  constexpr double third = boost::math::constants::third<double>();
82  return (t_stress(0, 0) + t_stress(1, 1)) * third;
83 };

◆ trace() [2/2]

template<typename T >
double PlasticOps::trace ( FTensor::Tensor2_symmetric< T, 3 > &  t_stress)
inline

Definition at line 86 of file PlasticOpsGeneric.hpp.

86  {
87  constexpr double third = boost::math::constants::third<double>();
88  return (t_stress(0, 0) + t_stress(1, 1) + t_stress(2, 2)) * third;
89 };

◆ w()

double PlasticOps::w ( double  eqiv,
double  dot_tau,
double  f,
double  sigma_y,
double  sigma_Y 
)
inline
Examples
PlasticOpsGeneric.hpp.

Definition at line 276 of file PlasticOpsGeneric.hpp.

277  {
278  return (f - sigma_y) / sigma_Y + cn1 * (dot_tau);
279 };

Variable Documentation

◆ I

FTensor::Index<'I', 3> PlasticOps::I

[Common data]

Examples
PlasticOps.hpp, and PlasticOpsGeneric.hpp.

Definition at line 115 of file PlasticOps.hpp.

◆ J

FTensor::Index<'J', 3> PlasticOps::J
Examples
PlasticOps.hpp, and PlasticOpsGeneric.hpp.

Definition at line 116 of file PlasticOps.hpp.

◆ M

FTensor::Index<'M', 3> PlasticOps::M

◆ N

FTensor::Index<'N', 3> PlasticOps::N
Examples
PlasticOps.hpp, and PlasticOpsGeneric.hpp.

Definition at line 118 of file PlasticOps.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
ContactOps::sign
double sign(double x)
Definition: ContactOps.hpp:580
PlasticOps::diff_deviator
auto diff_deviator(FTensor::Ddg< double, 3, 3 > &&t_diff_stress)
Definition: PlasticOpsGeneric.hpp:161
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
sigmaY
double sigmaY
Yield stress.
Definition: plastic.cpp:123
PlasticOps::diff_constrain_df
auto diff_constrain_df(double sign)
Definition: PlasticOpsGeneric.hpp:310
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
PlasticOps::CommonHenkyPtr
boost::shared_ptr< HenckyOps::CommonData > CommonHenkyPtr
Definition: PlasticOps.hpp:244
MoFEM::OpCalculateScalarFieldValuesDot
OpCalculateScalarFieldValuesFromPetscVecImpl< PetscData::CTX_SET_X_T > OpCalculateScalarFieldValuesDot
Definition: UserDataOperators.hpp:273
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
J
FTensor::Index< 'J', DIM1 > J
Definition: level_set.cpp:30
HenckyOps
Definition: HenckyOps.hpp:11
FTensor::Tensor2_symmetric
Definition: Tensor2_symmetric_value.hpp:13
scale
double scale
Definition: plastic.cpp:119
zeta
double zeta
Viscous hardening.
Definition: plastic.cpp:126
FTensor::Tensor2
Definition: Tensor2_value.hpp:16
MoFEM::DeprecatedCoreInterface
Deprecated interface functions.
Definition: DeprecatedCoreInterface.hpp:16
I
constexpr IntegrationType I
Definition: operators_tests.cpp:31
MoFEM::Sev
MoFEM::LogManager::SeverityLevel Sev
Definition: CoreTemplates.hpp:17
PlasticOps::M
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
FTensor::Number< 2 >
CHKERR
#define CHKERR
Inline error check.
Definition: definitions.h:548
MoFEM::CoreInterface::get_moab
virtual moab::Interface & get_moab()=0
OpInternalForcePiola
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForcePiola
Definition: seepage.cpp:66
BlockData
Definition: ElasticityMixedFormulation.hpp:12
PlasticOps::trace
double trace(FTensor::Tensor2_symmetric< T, 3 > &t_stress)
Definition: PlasticOpsGeneric.hpp:86
convert.type
type
Definition: convert.py:64
cn1
double cn1
Definition: plastic.cpp:132
EshelbianPlasticity::P
@ P
Definition: EshelbianContact.cpp:197
C1_k
double C1_k
Kinematic hardening.
Definition: plastic.cpp:129
poisson_ratio
double poisson_ratio
Poisson ratio.
Definition: plastic.cpp:122
cn0
double cn0
Definition: plastic.cpp:131
OpKCauchy
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::BiLinearForm< GAUSS >::OpGradSymTensorGrad< 1, SPACE_DIM, SPACE_DIM, 0 > OpKCauchy
[Linear elastic problem]
Definition: thermo_elastic.cpp:34
MoFEM::L
VectorDouble L
Definition: Projection10NodeCoordsOnField.cpp:124
HenckyOps::HenkyIntegrators
Definition: HenckyOps.hpp:573
visH
double visH
Viscous hardening.
Definition: plastic.cpp:125
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
OpInternalForceCauchy
FormsIntegrators< DomainEleOp >::Assembly< PETSC >::LinearForm< GAUSS >::OpGradTimesSymTensor< 1, SPACE_DIM, SPACE_DIM > OpInternalForceCauchy
Definition: thermo_elastic.cpp:38
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
FTensor::Index< 'i', DIM >
convert.n
n
Definition: convert.py:82
N
const int N
Definition: speed_test.cpp:3
b_iso
double b_iso
Saturation exponent.
Definition: plastic.cpp:128
Range
DomainEleOp
MOFEM_TAG_AND_LOG
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
Definition: LogManager.hpp:362
shear_modulus_G
double shear_modulus_G
Definition: dynamic_first_order_con_law.cpp:97
FTensor::Dg
Definition: Dg_value.hpp:9
PlasticOps::deviator
auto deviator(FTensor::Tensor2_symmetric< T, 3 > &t_stress, double trace, FTensor::Tensor2_symmetric< double, 3 > &&t_alpha)
Definition: PlasticOpsGeneric.hpp:119
HenckyOps::f
auto f
Definition: HenckyOps.hpp:15
is_large_strains
PetscBool is_large_strains
Large strains.
Definition: plastic.cpp:116
j
FTensor::Index< 'j', 3 > j
Definition: matrix_function.cpp:19
FTensor::Ddg
Definition: Ddg_value.hpp:7
MOFEM_DATA_INCONSISTENCY
@ MOFEM_DATA_INCONSISTENCY
Definition: definitions.h:31
Qinf
double Qinf
Saturation yield stress.
Definition: plastic.cpp:127
m
FTensor::Index< 'm', 3 > m
Definition: shallow_wave.cpp:80
FTensor::Kronecker_Delta_symmetric
Kronecker Delta class symmetric.
Definition: Kronecker_Delta.hpp:49
third
constexpr double third
Definition: EshelbianADOL-C.cpp:17
PlasticOps::I
FTensor::Index< 'I', 3 > I
[Common data]
Definition: PlasticOps.hpp:115
k
FTensor::Index< 'k', 3 > k
Definition: matrix_function.cpp:20
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
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