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

#include <users_modules/tutorials/adv-0/src/PlasticOpsGeneric.hpp>

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

Public Member Functions

 OpCalculatePlasticityImpl (const std::string field_name, boost::shared_ptr< CommonData > common_data_ptr, boost::shared_ptr< MatrixDouble > m_D_ptr)
 
MoFEMErrorCode doWork (int side, EntityType type, EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
- Public Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
int getNumNodes ()
 get element number of nodes More...
 
const EntityHandlegetConn ()
 get element connectivity More...
 
double getVolume () const
 element volume (linear geometry) More...
 
doublegetVolume ()
 element volume (linear geometry) More...
 
FTensor::Tensor2< double *, 3, 3 > & getJac ()
 get element Jacobian More...
 
FTensor::Tensor2< double *, 3, 3 > & getInvJac ()
 get element inverse Jacobian More...
 
VectorDoublegetCoords ()
 nodal coordinates More...
 
VolumeElementForcesAndSourcesCoregetVolumeFE () const
 return pointer to Generic Volume Finite Element object More...
 
- Public Member Functions inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
 UserDataOperator (const FieldSpace space, const char type=OPSPACE, const bool symm=true)
 
 UserDataOperator (const std::string field_name, const char type, const bool symm=true)
 
 UserDataOperator (const std::string row_field_name, const std::string col_field_name, const char type, const bool symm=true)
 
boost::shared_ptr< const NumeredEntFiniteElementgetNumeredEntFiniteElementPtr () const
 Return raw pointer to NumeredEntFiniteElement. More...
 
EntityHandle getFEEntityHandle () const
 Return finite element entity handle. More...
 
int getFEDim () const
 Get dimension of finite element. More...
 
EntityType getFEType () const
 Get dimension of finite element. More...
 
boost::weak_ptr< SideNumbergetSideNumberPtr (const int side_number, const EntityType type)
 Get the side number pointer. More...
 
EntityHandle getSideEntity (const int side_number, const EntityType type)
 Get the side entity. More...
 
int getNumberOfNodesOnElement () const
 Get the number of nodes on finite element. More...
 
MoFEMErrorCode getProblemRowIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get row indices. More...
 
MoFEMErrorCode getProblemColIndices (const std::string filed_name, const EntityType type, const int side, VectorInt &indices) const
 Get col indices. More...
 
const FEMethodgetFEMethod () const
 Return raw pointer to Finite Element Method object. More...
 
int getOpType () const
 Get operator types. More...
 
void setOpType (const OpType type)
 Set operator type. More...
 
void addOpType (const OpType type)
 Add operator type. More...
 
int getNinTheLoop () const
 get number of finite element in the loop More...
 
int getLoopSize () const
 get size of elements in the loop More...
 
std::string getFEName () const
 Get name of the element. More...
 
ForcesAndSourcesCoregetPtrFE () const
 
ForcesAndSourcesCoregetSidePtrFE () const
 
ForcesAndSourcesCoregetRefinePtrFE () const
 
const PetscData::SwitchesgetDataCtx () const
 
const KspMethod::KSPContext getKSPCtx () const
 
const SnesMethod::SNESContext getSNESCtx () const
 
const TSMethod::TSContext getTSCtx () const
 
Vec getKSPf () const
 
Mat getKSPA () const
 
Mat getKSPB () const
 
Vec getSNESf () const
 
Vec getSNESx () const
 
Mat getSNESA () const
 
Mat getSNESB () const
 
Vec getTSu () const
 
Vec getTSu_t () const
 
Vec getTSu_tt () const
 
Vec getTSf () const
 
Mat getTSA () const
 
Mat getTSB () const
 
int getTSstep () const
 
double getTStime () const
 
double getTStimeStep () const
 
double getTSa () const
 
double getTSaa () const
 
MatrixDoublegetGaussPts ()
 matrix of integration (Gauss) points for Volume Element More...
 
auto getFTensor0IntegrationWeight ()
 Get integration weights. More...
 
MatrixDoublegetCoordsAtGaussPts ()
 Gauss points and weight, matrix (nb. of points x 3) More...
 
auto getFTensor1CoordsAtGaussPts ()
 Get coordinates at integration points assuming linear geometry. More...
 
double getMeasure () const
 get measure of element More...
 
doublegetMeasure ()
 get measure of element More...
 
MoFEMErrorCode loopSide (const string &fe_name, ForcesAndSourcesCore *side_fe, const size_t dim, const EntityHandle ent_for_side=0, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy, AdjCache *adj_cache=nullptr)
 User calls this function to loop over elements on the side of face. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopThis (const string &fe_name, ForcesAndSourcesCore *this_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over the same element using a different set of integration points. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopParent (const string &fe_name, ForcesAndSourcesCore *parent_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
MoFEMErrorCode loopChildren (const string &fe_name, ForcesAndSourcesCore *child_fe, const int verb=QUIET, const LogManager::SeverityLevel sev=Sev::noisy)
 User calls this function to loop over parent elements. This function calls finite element with its operator to do calculations. More...
 
- Public Member Functions inherited from MoFEM::DataOperator
 DataOperator (const bool symm=true)
 
virtual ~DataOperator ()=default
 
virtual MoFEMErrorCode doWork (int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)
 Operator for bi-linear form, usually to calculate values on left hand side. More...
 
virtual MoFEMErrorCode opLhs (EntitiesFieldData &row_data, EntitiesFieldData &col_data)
 
virtual MoFEMErrorCode doWork (int side, EntityType type, EntitiesFieldData::EntData &data)
 Operator for linear form, usually to calculate values on right hand side. More...
 
virtual MoFEMErrorCode opRhs (EntitiesFieldData &data, const bool error_if_no_base=false)
 
bool getSymm () const
 Get if operator uses symmetry of DOFs or not. More...
 
void setSymm ()
 set if operator is executed taking in account symmetry More...
 
void unSetSymm ()
 unset if operator is executed for non symmetric problem More...
 

Protected Attributes

boost::shared_ptr< CommonDatacommonDataPtr
 
boost::shared_ptr< MatrixDouble > mDPtr
 
- Protected Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
ForcesAndSourcesCoreptrFE
 

Additional Inherited Members

- Public Types inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
enum  OpType {
  OPROW = 1 << 0 , OPCOL = 1 << 1 , OPROWCOL = 1 << 2 , OPSPACE = 1 << 3 ,
  OPLAST = 1 << 3
}
 Controls loop over entities on element. More...
 
using AdjCache = std::map< EntityHandle, std::vector< boost::weak_ptr< NumeredEntFiniteElement > > >
 
- Public Types inherited from MoFEM::DataOperator
using DoWorkLhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int row_side, int col_side, EntityType row_type, EntityType col_type, EntitiesFieldData::EntData &row_data, EntitiesFieldData::EntData &col_data)>
 
using DoWorkRhsHookFunType = boost::function< MoFEMErrorCode(DataOperator *op_ptr, int side, EntityType type, EntitiesFieldData::EntData &data)>
 
- Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
char opType
 
std::string rowFieldName
 
std::string colFieldName
 
FieldSpace sPace
 
- Public Attributes inherited from MoFEM::DataOperator
DoWorkLhsHookFunType doWorkLhsHook
 
DoWorkRhsHookFunType doWorkRhsHook
 
bool sYmm
 If true assume that matrix is symmetric structure. More...
 
std::array< bool, MBMAXTYPE > doEntities
 If true operator is executed for entity. More...
 
booldoVertices
 \deprectaed If false skip vertices More...
 
booldoEdges
 \deprectaed If false skip edges More...
 
booldoQuads
 \deprectaed More...
 
booldoTris
 \deprectaed More...
 
booldoTets
 \deprectaed More...
 
booldoPrisms
 \deprectaed More...
 
- Static Public Attributes inherited from MoFEM::ForcesAndSourcesCore::UserDataOperator
static const char *const OpTypeNames []
 
- Protected Member Functions inherited from MoFEM::VolumeElementForcesAndSourcesCore::UserDataOperator
MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 
virtual MoFEMErrorCode setPtrFE (ForcesAndSourcesCore *ptr)
 

Detailed Description

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

Definition at line 457 of file PlasticOpsGeneric.hpp.

Constructor & Destructor Documentation

◆ OpCalculatePlasticityImpl()

template<int DIM, typename DomainEleOp >
PlasticOps::OpCalculatePlasticityImpl< DIM, GAUSS, DomainEleOp >::OpCalculatePlasticityImpl ( const std::string  field_name,
boost::shared_ptr< CommonData common_data_ptr,
boost::shared_ptr< MatrixDouble >  m_D_ptr 
)

Definition at line 469 of file PlasticOpsGeneric.hpp.

473 commonDataPtr(common_data_ptr), mDPtr(m_D_ptr) {
474 // Opetor is only executed for vertices
475 std::fill(&DomainEleOp::doEntities[MBEDGE],
476 &DomainEleOp::doEntities[MBMAXTYPE], false);
477}
DomainEle::UserDataOperator DomainEleOp
Finire element operator type.
constexpr auto field_name
std::array< bool, MBMAXTYPE > doEntities
If true operator is executed for entity.
@ OPROW
operator doWork function is executed on FE rows

Member Function Documentation

◆ doWork()

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

Operator for linear form, usually to calculate values on right hand side.

< material parameters

Reimplemented from MoFEM::DataOperator.

Definition at line 480 of file PlasticOpsGeneric.hpp.

481 {
483
490
491 auto &params = commonDataPtr->blockParams; ///< material parameters
492
493 const size_t nb_gauss_pts = DomainEleOp::getGaussPts().size2();
495 auto t_tau = getFTensor0FromVec(commonDataPtr->plasticTau);
496 auto t_tau_dot = getFTensor0FromVec(commonDataPtr->plasticTauDot);
497 auto t_f = getFTensor0FromVec(commonDataPtr->plasticSurface);
498 auto t_flow = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->plasticFlow);
499 auto t_plastic_strain =
500 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->plasticStrain);
501 auto t_plastic_strain_dot =
502 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->plasticStrainDot);
503 auto t_stress =
504 getFTensor2SymmetricFromMat<DIM>(*(commonDataPtr->mStressPtr));
505
506 auto t_D = getFTensor4DdgFromMat<DIM, DIM, 0>(*commonDataPtr->mDPtr);
507 auto t_D_Op = getFTensor4DdgFromMat<DIM, DIM, 0>(*mDPtr);
508
509 auto t_diff_plastic_strain = diff_tensor(FTensor::Number<DIM>());
510 auto t_diff_deviator = diff_deviator(diff_tensor(FTensor::Number<DIM>()));
511
512 FTensor::Ddg<double, DIM, DIM> t_flow_dir_dstress;
513 FTensor::Ddg<double, DIM, DIM> t_flow_dir_dstrain;
514 t_flow_dir_dstress(i, j, k, l) =
515 1.5 * (t_diff_deviator(M, N, i, j) * t_diff_deviator(M, N, k, l));
516 t_flow_dir_dstrain(i, j, k, l) =
517 t_flow_dir_dstress(i, j, m, n) * t_D_Op(m, n, k, l);
518
519
520 auto t_alpha_dir =
521 kinematic_hardening_dplastic_strain<DIM>(params[CommonData::C1_k]);
522
523 commonDataPtr->resC.resize(nb_gauss_pts, false);
524 commonDataPtr->resCdTau.resize(nb_gauss_pts, false);
525 commonDataPtr->resCdStrain.resize(size_symm, nb_gauss_pts, false);
526 commonDataPtr->resCdPlasticStrain.resize(size_symm, nb_gauss_pts, false);
527 commonDataPtr->resFlow.resize(size_symm, nb_gauss_pts, false);
528 commonDataPtr->resFlowDtau.resize(size_symm, nb_gauss_pts, false);
529 commonDataPtr->resFlowDstrain.resize(size_symm * size_symm, nb_gauss_pts,
530 false);
531 commonDataPtr->resFlowDstrainDot.resize(size_symm * size_symm, nb_gauss_pts,
532 false);
533
534 commonDataPtr->resC.clear();
535 commonDataPtr->resCdTau.clear();
536 commonDataPtr->resCdStrain.clear();
537 commonDataPtr->resCdPlasticStrain.clear();
538 commonDataPtr->resFlow.clear();
539 commonDataPtr->resFlowDtau.clear();
540 commonDataPtr->resFlowDstrain.clear();
541 commonDataPtr->resFlowDstrainDot.clear();
542
543 auto t_res_c = getFTensor0FromVec(commonDataPtr->resC);
544 auto t_res_c_dtau = getFTensor0FromVec(commonDataPtr->resCdTau);
545 auto t_res_c_dstrain =
546 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resCdStrain);
547 auto t_res_c_plastic_strain =
548 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resCdPlasticStrain);
549 auto t_res_flow = getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resFlow);
550 auto t_res_flow_dtau =
551 getFTensor2SymmetricFromMat<DIM>(commonDataPtr->resFlowDtau);
552 auto t_res_flow_dstrain =
553 getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->resFlowDstrain);
554 auto t_res_flow_dplastic_strain =
555 getFTensor4DdgFromMat<DIM, DIM>(commonDataPtr->resFlowDstrainDot);
556
557 auto next = [&]() {
558 ++t_tau;
559 ++t_tau_dot;
560 ++t_f;
561 ++t_flow;
562 ++t_plastic_strain;
563 ++t_plastic_strain_dot;
564 ++t_stress;
565 ++t_res_c;
566 ++t_res_c_dtau;
567 ++t_res_c_dstrain;
568 ++t_res_c_plastic_strain;
569 ++t_res_flow;
570 ++t_res_flow_dtau;
571 ++t_res_flow_dstrain;
572 ++t_res_flow_dplastic_strain;
573 ++t_w;
574 };
575
576 auto get_avtive_pts = [&]() {
577 int nb_points_avtive_on_elem = 0;
578 int nb_points_on_elem = 0;
579
580 auto t_tau = getFTensor0FromVec(commonDataPtr->plasticTau);
581 auto t_tau_dot = getFTensor0FromVec(commonDataPtr->plasticTauDot);
582 auto t_f = getFTensor0FromVec(commonDataPtr->plasticSurface);
583 auto t_plastic_strain_dot =
584 getFTensor2SymmetricFromMat<SPACE_DIM>(commonDataPtr->plasticStrainDot);
585
586 for (auto &f : commonDataPtr->plasticSurface) {
587 auto eqiv = equivalent_strain_dot(t_plastic_strain_dot);
588 const auto ww = w(
589 eqiv, t_tau_dot, t_f,
590 iso_hardening(t_tau, params[CommonData::H], params[CommonData::QINF],
591 params[CommonData::BISO], params[CommonData::SIGMA_Y]),
592 params[CommonData::SIGMA_Y]);
593 const auto sign_ww = constrian_sign(ww);
594
595 ++nb_points_on_elem;
596 if (sign_ww > 0) {
597 ++nb_points_avtive_on_elem;
598 }
599
600 ++t_tau;
601 ++t_tau_dot;
602 ++t_f;
603 ++t_plastic_strain_dot;
604 }
605
606 int &active_points = PlasticOps::CommonData::activityData[0];
607 int &avtive_full_elems = PlasticOps::CommonData::activityData[1];
608 int &avtive_elems = PlasticOps::CommonData::activityData[2];
609 int &nb_points = PlasticOps::CommonData::activityData[3];
610 int &nb_elements = PlasticOps::CommonData::activityData[4];
611
612 ++nb_elements;
613 nb_points += nb_points_on_elem;
614 if (nb_points_avtive_on_elem > 0) {
615 ++avtive_elems;
616 active_points += nb_points_avtive_on_elem;
617 if (nb_points_avtive_on_elem == nb_points_on_elem) {
618 ++avtive_full_elems;
619 }
620 }
621
622 if (nb_points_avtive_on_elem != nb_points_on_elem)
623 return 1;
624 else
625 return 0;
626 };
627
628 if (DomainEleOp::getTSCtx() == TSMethod::TSContext::CTX_TSSETIJACOBIAN) {
629 get_avtive_pts();
630 }
631
632 for (auto &f : commonDataPtr->plasticSurface) {
633
634 auto eqiv = equivalent_strain_dot(t_plastic_strain_dot);
635 auto t_diff_eqiv = diff_equivalent_strain_dot(eqiv, t_plastic_strain_dot,
636 t_diff_plastic_strain,
638
639 const auto sigma_y =
640 iso_hardening(t_tau, params[CommonData::H], params[CommonData::QINF],
641 params[CommonData::BISO], params[CommonData::SIGMA_Y]);
642 const auto d_sigma_y =
643 iso_hardening_dtau(t_tau, params[CommonData::H],
644 params[CommonData::QINF], params[CommonData::BISO]);
645
646 auto ww = w(eqiv, t_tau_dot, t_f, sigma_y, params[CommonData::SIGMA_Y]);
647 auto abs_ww = constrain_abs(ww);
648 auto sign_ww = constrian_sign(ww);
649
650 auto c = constraint(eqiv, t_tau_dot, t_f, sigma_y, abs_ww,
651 params[CommonData::VIS_H], params[CommonData::SIGMA_Y]);
652 auto c_dot_tau = diff_constrain_ddot_tau(sign_ww, eqiv, t_tau_dot,
653 params[CommonData::VIS_H],
654 params[CommonData::SIGMA_Y]);
655 auto c_equiv = diff_constrain_deqiv(sign_ww, eqiv, t_tau_dot,
656 params[CommonData::SIGMA_Y]);
657 auto c_sigma_y = diff_constrain_dsigma_y(sign_ww);
658 auto c_f = diff_constrain_df(sign_ww);
659
660 auto t_dev_stress = deviator(
661
662 t_stress, trace(t_stress),
663
664 kinematic_hardening(t_plastic_strain, params[CommonData::C1_k])
665
666 );
667
669 t_flow_dir(k, l) = 1.5 * (t_dev_stress(I, J) * t_diff_deviator(I, J, k, l));
671 t_flow_dstrain(i, j) = t_flow(k, l) * t_D_Op(k, l, i, j);
672
673 auto get_res_c = [&]() { return c; };
674
675 auto get_res_c_dstrain = [&](auto &t_diff_res) {
676 t_diff_res(i, j) = c_f * t_flow_dstrain(i, j);
677 };
678
679 auto get_res_c_dplastic_strain = [&](auto &t_diff_res) {
680 t_diff_res(i, j) = (this->getTSa() * c_equiv) * t_diff_eqiv(i, j);
681 t_diff_res(k, l) -= c_f * t_flow(i, j) * t_alpha_dir(i, j, k, l);
682 };
683
684 auto get_res_c_dtau = [&]() {
685 return this->getTSa() * c_dot_tau + c_sigma_y * d_sigma_y;
686 };
687
688 auto get_res_c_plastic_strain = [&](auto &t_diff_res) {
689 t_diff_res(k, l) = -c_f * t_flow(i, j) * t_alpha_dir(i, j, k, l);
690 };
691
692 auto get_res_flow = [&](auto &t_res_flow) {
693 const auto a = sigma_y;
694 const auto b = t_tau_dot;
695 t_res_flow(k, l) = a * t_plastic_strain_dot(k, l) - b * t_flow_dir(k, l);
696 };
697
698 auto get_res_flow_dtau = [&](auto &t_res_flow_dtau) {
699 const auto da = d_sigma_y;
700 const auto db = this->getTSa();
701 t_res_flow_dtau(k, l) =
702 da * t_plastic_strain_dot(k, l) - db * t_flow_dir(k, l);
703 };
704
705 auto get_res_flow_dstrain = [&](auto &t_res_flow_dstrain) {
706 const auto b = t_tau_dot;
707 t_res_flow_dstrain(m, n, k, l) = -t_flow_dir_dstrain(m, n, k, l) * b;
708 };
709
710 auto get_res_flow_dplastic_strain = [&](auto &t_res_flow_dplastic_strain) {
711 const auto a = sigma_y;
712 t_res_flow_dplastic_strain(m, n, k, l) =
713 (a * this->getTSa()) * t_diff_plastic_strain(m, n, k, l);
714 const auto b = t_tau_dot;
715 t_res_flow_dplastic_strain(m, n, i, j) +=
716 (t_flow_dir_dstrain(m, n, k, l) * t_alpha_dir(k, l, i, j)) * b;
717 };
718
719 t_res_c = get_res_c();
720 get_res_flow(t_res_flow);
721
722 if (this->getTSCtx() == TSMethod::TSContext::CTX_TSSETIJACOBIAN) {
723 t_res_c_dtau = get_res_c_dtau();
724 get_res_c_dstrain(t_res_c_dstrain);
725 get_res_c_dplastic_strain(t_res_c_plastic_strain);
726 get_res_flow_dtau(t_res_flow_dtau);
727 get_res_flow_dstrain(t_res_flow_dstrain);
728 get_res_flow_dplastic_strain(t_res_flow_dplastic_strain);
729 }
730
731 next();
732 }
733
735}
constexpr double a
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
Definition: definitions.h:346
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
Definition: definitions.h:416
FTensor::Index< 'n', SPACE_DIM > n
FTensor::Index< 'm', SPACE_DIM > m
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
FTensor::Index< 'l', 3 > l
FTensor::Index< 'j', 3 > j
FTensor::Index< 'k', 3 > k
Params params
double constrian_sign(double x)
double diff_constrain_ddot_tau(double sign, double eqiv, double dot_tau, double vis_H, double sigma_Y)
FTensor::Index< 'M', 3 > M
Definition: PlasticOps.hpp:117
double diff_constrain_deqiv(double sign, double eqiv, double dot_tau, double sigma_Y)
auto diff_equivalent_strain_dot(const T1 eqiv, T2 &t_plastic_strain_dot, T3 &t_diff_plastic_strain, FTensor::Number< DIM >)
double constraint(double eqiv, double dot_tau, double f, double sigma_y, double abs_w, double vis_H, double sigma_Y)
auto diff_constrain_df(double sign)
FTensor::Index< 'N', 3 > N
Definition: PlasticOps.hpp:118
auto diff_constrain_dsigma_y(double sign)
FTensor::Index< 'I', 3 > I
[Common data]
Definition: PlasticOps.hpp:115
double constrain_abs(double x)
FTensor::Index< 'J', 3 > J
Definition: PlasticOps.hpp:116
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)
auto kinematic_hardening(FTensor::Tensor2_symmetric< T, DIM > &t_plastic_strain, double C1_k)
Definition: plastic.cpp:142
double iso_hardening_dtau(double tau, double H, double Qinf, double b_iso)
Definition: plastic.cpp:128
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:123
auto getFTensor0IntegrationWeight()
Get integration weights.
MatrixDouble & getGaussPts()
matrix of integration (Gauss) points for Volume Element
static std::array< int, 5 > activityData
Definition: PlasticOps.hpp:107

Member Data Documentation

◆ commonDataPtr

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

Definition at line 464 of file PlasticOpsGeneric.hpp.

◆ mDPtr

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

Definition at line 465 of file PlasticOpsGeneric.hpp.


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