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

[J2 2D] More...

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

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

Public Types

enum  ModelParams {
  LAMBDA , MU , NUP , SIGMA_YT ,
  SIGMA_YC , HT , HC , NT ,
  NC , LAST_PARAM
}
 
using ParamsArray = std::array<double, LAST_PARAM>
 
- Public Types inherited from ADOLCPlasticity::ClosestPointProjection
enum  TypesTags { W = 0 , Y , H , LAST_TAPE }
 

Public Member Functions

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

Public Attributes

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

Detailed Description

[J2 2D]

Paraboloidal yield criterion

See paper for reference [ullah2017three] and [stier2015comparing]

Linear hardening

\[ \psi = \frac{1}{2} \lambda \textrm{tr}[\varepsilon]^2 + \mu \varepsilon : \varepsilon + \sigma_{yt}\alpha_0 + \frac{1}{2} H_t \alpha_0^2 + \sigma_{yc}\alpha_1 + \frac{1}{2} H_c \alpha_1^2 \]

Exponential hardening

\[ \psi = \frac{1}{2} \lambda \textrm{tr}[\varepsilon]^2 + \mu \varepsilon : \varepsilon\\ + (\sigma_{yt}+H_t)\alpha_0 + \frac{H_t}{n_t} \exp{(-n_t \alpha_0)}\\ + (\sigma_{yc}+H_c)\alpha_0 + \frac{H_c}{n_c} \exp{(-n_c \alpha_1)} \]

\[ I_1 = \textrm{tr} (\boldsymbol{\sigma}) \]

\[ \eta=\textrm{dev}[\boldsymbol{\sigma}] \]

\[ J_2 = \frac{1}{2} \eta:\eta \]

\[ y = 6J_2 + 2I_1\left(\overline{\alpha_1}-\overline{\alpha_0}\right) - 2\overline{\alpha_0} \,\overline{\alpha_1} \]

where

\[ \overline{\alpha_0}=\frac{\partial \psi}{\partial \alpha_0}=\sigma_{yt} + H_t \alpha_0 \]

\[ \overline{\alpha_1}=\frac{\partial \psi}{\partial \alpha_1}=\sigma_{yc} + H_c \alpha_1 \]

\[ \Psi = 6J_2 + 2\alpha I_1 \left(\overline{\alpha_1}-\overline{\alpha_0}\right) - 2\overline{\alpha_0} \,\overline{\alpha_1} \]

\[ \alpha= \frac{1-2\nu_p}{1+\nu_p} \]

Definition at line 596 of file ADOLCPlasticityMaterialModels.hpp.

Member Typedef Documentation

◆ ParamsArray

Member Enumeration Documentation

◆ ModelParams

Constructor & Destructor Documentation

◆ ParaboloidalPlasticity()

ADOLCPlasticity::ParaboloidalPlasticity::ParaboloidalPlasticity ( )
inline

Member Function Documentation

◆ addMatBlockOps()

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

Get model parameters from blocks.

Reimplemented from ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 714 of file ADOLCPlasticityMaterialModels.hpp.

716 {
718
719 struct OpMatBlocks : public ForcesAndSourcesCore::UserDataOperator {
720
721 using OP = ForcesAndSourcesCore::UserDataOperator;
722
723 OpMatBlocks(boost::shared_ptr<ParamsArray> params_array_ptr,
724 boost::shared_ptr<ParamsArray> def_params_array_ptr,
725 MoFEM::Interface &m_field, Sev sev,
726 std::vector<const CubitMeshSets *> meshset_vec_ptr)
727 : OP(NOSPACE, OP::OPSPACE), paramsArrayPtr(params_array_ptr),
728 defParamsArray(def_params_array_ptr) {
729 CHK_THROW_MESSAGE(extractBlockData(m_field, meshset_vec_ptr, sev),
730 "Can not get data from block");
731 }
732
733 MoFEMErrorCode doWork(int side, EntityType type,
734 EntitiesFieldData::EntData &data) {
736
737 for (auto &b : blockData) {
738
739 if (b.blockEnts.find(getFEEntityHandle()) != b.blockEnts.end()) {
740 std::copy(b.paramsArray.begin(), b.paramsArray.end(),
741 paramsArrayPtr->begin());
743 }
744 }
745
746 std::copy(defParamsArray->begin(), defParamsArray->end(),
747 paramsArrayPtr->begin());
748
750 }
751
752 private:
753 boost::shared_ptr<ParamsArray> paramsArrayPtr;
754 boost::shared_ptr<ParamsArray> defParamsArray;
755
756 struct BlockData {
757 ParamsArray paramsArray;
758 Range blockEnts;
759 };
760
761 std::vector<BlockData> blockData;
762
764 extractBlockData(MoFEM::Interface &m_field,
765 std::vector<const CubitMeshSets *> meshset_vec_ptr,
766 Sev sev) {
768
769 for (auto m : meshset_vec_ptr) {
770 MOFEM_TAG_AND_LOG("ADOLCPlasticityWold", sev, "Paraboloidal") << *m;
771 std::vector<double> block_data;
772 CHKERR m->getAttributes(block_data);
773 if (block_data.size() != 9) {
774 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
775 "Expected that block has nine attributes");
776 }
777 auto get_block_ents = [&]() {
778 Range ents;
779 CHKERR
780 m_field.get_moab().get_entities_by_handle(m->meshset, ents, true);
781 return ents;
782 };
783
784 blockData.push_back({*defParamsArray, get_block_ents()});
785 std::copy(block_data.begin(), block_data.end(),
786 blockData.back().paramsArray.begin());
787
788 const auto young_modulus = blockData.back().paramsArray[LAMBDA];
789 const auto poisson_ratio = blockData.back().paramsArray[MU];
790
791 // It is assumed that user provide young_modulus and poisson_ratio in
792 // first two arguments of the block
793
794 MOFEM_LOG("ADOLCPlasticityWold", sev)
795 << "Young modulus: " << (blockData.back().paramsArray)[LAMBDA];
796 MOFEM_LOG("ADOLCPlasticityWold", sev)
797 << "Poisson ratio: " << (blockData.back().paramsArray)[MU];
798
799 blockData.back().paramsArray[LAMBDA] =
801 blockData.back().paramsArray[MU] = MU(young_modulus, poisson_ratio);
802
803 MOFEM_LOG("ADOLCPlasticityWold", sev)
804 << "LAMBDA: " << (blockData.back().paramsArray)[LAMBDA];
805 MOFEM_LOG("ADOLCPlasticityWold", sev)
806 << "MU: " << (blockData.back().paramsArray)[MU];
807 MOFEM_LOG("ADOLCPlasticityWold", sev)
808 << "NUP: " << (blockData.back().paramsArray)[NUP];
809 MOFEM_LOG("ADOLCPlasticityWold", sev)
810 << "SIGMA_YT: " << (blockData.back().paramsArray)[SIGMA_YT];
811 MOFEM_LOG("ADOLCPlasticityWold", sev)
812 << "SIGMA_YC: " << (blockData.back().paramsArray)[SIGMA_YC];
813 MOFEM_LOG("ADOLCPlasticityWold", sev)
814 << "HT: " << (blockData.back().paramsArray)[HT];
815 MOFEM_LOG("ADOLCPlasticityWold", sev)
816 << "HC: " << (blockData.back().paramsArray)[HC];
817 MOFEM_LOG("ADOLCPlasticityWold", sev)
818 << "NT: " << (blockData.back().paramsArray)[NT];
819 MOFEM_LOG("ADOLCPlasticityWold", sev)
820 << "NC: " << (blockData.back().paramsArray)[NC];
821 }
823 }
824 };
825
826 auto cubit_meshsets_vec_ptr =
827 m_field.getInterface<MeshsetsManager>()->getCubitMeshsetPtr(
828 std::regex((boost::format("%s(.*)") % block_name).str()));
829
830 pip.push_back(new OpMatBlocks(paramsArrayPtr, defaultParamsArrayPtr,
831 m_field, sev, cubit_meshsets_vec_ptr));
832
834 }
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MU(E, NU)
Definition fem_tools.h:23
#define LAMBDA(E, NU)
Definition fem_tools.h:22
#define MOFEM_LOG(channel, severity)
Log.
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
double young_modulus
Young modulus.
Definition plastic.cpp:125
double poisson_ratio
Poisson ratio.
Definition plastic.cpp:126
FTensor::Index< 'm', 3 > m
virtual moab::Interface & get_moab()=0
Deprecated interface functions.
MoFEMErrorCode getInterface(IFACE *&iface) const
Get interface reference to pointer of interface.

◆ codedHessianW()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::codedHessianW ( vector< double * > hessian)
inlinevirtual

[Paraboloidal free energy]

Reimplemented from ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 900 of file ADOLCPlasticityMaterialModels.hpp.

901 {
903 // first 6x6 is \frac{\partial^2 \psi}{\partial \varepsilon^p^2}
904 // lower symmetric part
905
906 lambda = (*paramsArrayPtr)[0];
907 mu = (*paramsArrayPtr)[1];
908
909 Ht = (*paramsArrayPtr)[5];
910 Hc = (*paramsArrayPtr)[6];
911
912 double lambda_plus_2mu = lambda + 2 * mu;
913
914 // Diagonal elements
915 hessian[0][0] = lambda_plus_2mu;
916 hessian[1][1] = lambda_plus_2mu;
917 hessian[2][2] = lambda_plus_2mu;
918
919 // Off-diagonal elements
920 hessian[1][0] = lambda;
921 hessian[2][0] = lambda;
922 hessian[2][1] = lambda;
923
924 // Shear components
925 for (int i = 3; i <= 5; ++i) {
926 hessian[i][i] = mu;
927 }
928
929 //next 6x6 is \frac{\partial^2 \psi}{\partial \varepsilon \partial \varepsilon^p}
930 double neg_lambda_plus_2mu = -(lambda + 2 * mu);
931 double neg_lambda = -lambda;
932 double neg_mu = -mu;
933
934 for (int i = 6; i <= 8; ++i) {
935 for (int j = 0; j <= 2; ++j) {
936 hessian[i][j] = (i - 6 == j) ? neg_lambda_plus_2mu : neg_lambda;
937 }
938 }
939
940 hessian[9][3] = neg_mu;
941 hessian[10][4] = neg_mu;
942 hessian[11][5] = neg_mu;
943
944 //next 6x6 is \frac{\partial^2 \psi}{\partial \varepsilon \partial \varepsilon}
945 hessian[6][6] = lambda_plus_2mu;
946 hessian[7][6] = lambda;
947 hessian[7][7] = lambda_plus_2mu;
948 hessian[8][6] = lambda;
949 hessian[8][7] = lambda;
950 hessian[8][8] = lambda_plus_2mu;
951 hessian[9][9] = mu;
952 hessian[10][10] = mu;
953 hessian[11][11] = mu;
954
955 // internalVariables is \frac{\partial^2 \psi}{\partial \alpha^2}
956 hessian[12][12] = Ht;
957 hessian[13][13] = Hc;
958
960 }

◆ evalF()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::evalF ( )
inline

[Paraboloidal yield function]

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 965 of file ADOLCPlasticityMaterialModels.hpp.

965 {
967 auto t_voight_stress_op = voight_to_stress_op();
968 auto t_vioght_stress = getFTensor1FromPtr<6>(&a_sTress[0]);
969 tStress(i, j) = t_voight_stress_op(i, j, Z) * t_vioght_stress(Z);
970 tR = tStress(i, i);
971 I1 = tR;
973 tR /= 3.;
974 tDeviator(i, j) = tStress(i, j) - tR * t_kd(i, j);
975 J2 = 0.5 * (tDeviator(i, j) * tDeviator(i, j));
977 }
Kronecker Delta class symmetric.
constexpr auto t_kd
FTensor::Dg< double, 3, 6 > voight_to_stress_op()
Op convert Vight stress vector to stress tensor.
FTensor::Tensor2_symmetric< adouble, 3 > tDeviator
FTensor::Tensor2_symmetric< adouble, 3 > tStress

◆ flowPotential()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::flowPotential ( )
inlinevirtual

[Paraboloidal yield function]

[Paraboloidal flow potential]

Implements ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 989 of file ADOLCPlasticityMaterialModels.hpp.

989 {
991 CHKERR evalF();
992 double alpha =
993 (1 - 2 * nup) /
994 (1 + nup); // relation between alpha and plastic poisson's ratio
995 a_h = 6 * J2 +
996 2 * alpha * I1 * (a_internalFluxes[1] - a_internalFluxes[0]) -
999 }
MoFEMErrorCode evalF()
[Paraboloidal yield function]

◆ freeHemholtzEnergy()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::freeHemholtzEnergy ( )
inlinevirtual

[Paraboloidal free energy]

Implements ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 864 of file ADOLCPlasticityMaterialModels.hpp.

864 {
866
867 auto p_lambda = mkparam(lambda); // 0
868 auto p_mu = mkparam(mu); // 1
869 auto p_nup = mkparam(nup); // 2
870 auto p_sigma_yt = mkparam(sigmaYt); // 3
871 auto p_sigma_yc = mkparam(sigmaYc); // 4
872 auto p_Ht = mkparam(Ht); // 5
873 auto p_Hc = mkparam(Hc); // 6
874 auto p_nt = mkparam(nt); // 7
875 auto p_nc = mkparam(nc); // 8
876
877 auto t_voight_strain_op = voight_to_strain_op();
878 auto t_voight_total_strain = getFTensor1FromPtr<6>(&a_sTrain[0]);
879 auto t_voight_plastic_strain = getFTensor1FromPtr<6>(&a_plasticStrain[0]);
880
881 tTotalStrain(i, j) = t_voight_strain_op(i, j, Z) * t_voight_total_strain(Z);
882 tPlasticStrain(i, j) =
883 t_voight_strain_op(i, j, Z) * t_voight_plastic_strain(Z);
885 tR = tElasticStrain(i, i);
886
887 a_w = 0.5 * p_lambda * tR * tR +
888 p_mu * (tElasticStrain(i, j) * tElasticStrain(i, j));
889
890 a_w += p_sigma_yt * a_internalVariables[0] +
891 0.5 * p_Ht * a_internalVariables[0] * a_internalVariables[0];
892 a_w += p_sigma_yc * a_internalVariables[1] +
893 0.5 * p_Hc * a_internalVariables[1] * a_internalVariables[1];
894
896 }
FTensor::Dg< double, 3, 6 > voight_to_strain_op()
Op convert Vight strain vector to strain tensor.
FTensor::Tensor2_symmetric< adouble, 3 > tPlasticStrain
FTensor::Tensor2_symmetric< adouble, 3 > tTotalStrain
FTensor::Tensor2_symmetric< adouble, 3 > tElasticStrain

◆ getDefaultMaterialParameters()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::getDefaultMaterialParameters ( )
inline
Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 640 of file ADOLCPlasticityMaterialModels.hpp.

640 {
642
643 double young_modulus = 1;
644 double poisson_ratio = 0.25;
645
646 sigmaYt = 1;
647 sigmaYc = 1;
648 Ht = 0.1;
649 Hc = 0.1;
650 nc = 1.;
651 nt = 1.;
652 nup = 0.;
653
654 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-young_modulus", &young_modulus,
655 PETSC_NULLPTR);
656 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-poisson_ratio", &poisson_ratio,
657 0);
658 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-sigma_yt", &sigmaYt, PETSC_NULLPTR);
659 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-sigma_yc", &sigmaYc, PETSC_NULLPTR);
660 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-Ht", &Ht, PETSC_NULLPTR);
661 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-Hc", &Hc, PETSC_NULLPTR);
662 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-nt", &nt, PETSC_NULLPTR);
663 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-nc", &nc, PETSC_NULLPTR);
664 CHKERR PetscOptionsGetScalar(PETSC_NULLPTR, "-nup", &nup, PETSC_NULLPTR);
665
666 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
667 << "Young modulus: " << young_modulus;
668 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
669 << "Poisson ratio: " << poisson_ratio;
670
673
674 defaultParamsArrayPtr = boost::make_shared<ParamsArray>();
675 (*defaultParamsArrayPtr)[LAMBDA] = lambda;
676 (*defaultParamsArrayPtr)[MU] = mu;
677 (*defaultParamsArrayPtr)[NUP] = nup;
678 (*defaultParamsArrayPtr)[SIGMA_YT] = sigmaYt;
679 (*defaultParamsArrayPtr)[SIGMA_YC] = sigmaYc;
680 (*defaultParamsArrayPtr)[HT] = Ht;
681 (*defaultParamsArrayPtr)[HC] = Hc;
682 (*defaultParamsArrayPtr)[NT] = nt;
683 (*defaultParamsArrayPtr)[NC] = nc;
684
685 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
686 << "LAMBDA: " << (*defaultParamsArrayPtr)[LAMBDA];
687 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
688 << "MU: " << (*defaultParamsArrayPtr)[MU];
689 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
690 << "NUP: " << (*defaultParamsArrayPtr)[NUP];
691 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
692 << "SIGMA_YT: " << (*defaultParamsArrayPtr)[SIGMA_YT];
693 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
694 << "SIGMA_YC: " << (*defaultParamsArrayPtr)[SIGMA_YC];
695 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
696 << "HT: " << (*defaultParamsArrayPtr)[HT];
697 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
698 << "HC: " << (*defaultParamsArrayPtr)[HC];
699 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
700 << "NT: " << (*defaultParamsArrayPtr)[NT];
701 MOFEM_LOG("ADOLCPlasticityWold", Sev::inform)
702 << "NC: " << (*defaultParamsArrayPtr)[NC];
703
704 paramsArrayPtr = boost::make_shared<ParamsArray>();
705 std::copy(defaultParamsArrayPtr->begin(), defaultParamsArrayPtr->end(),
706 paramsArrayPtr->begin());
707 oldParamsArrayPtr = boost::make_shared<ParamsArray>();
708 std::fill(oldParamsArrayPtr->begin(), oldParamsArrayPtr->end(), 0);
709
711 };
PetscErrorCode PetscOptionsGetScalar(PetscOptions *, const char pre[], const char name[], PetscScalar *dval, PetscBool *set)

◆ setParams()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::setParams ( short tag,
bool & recalculate_elastic_tangent )
inlinevirtual

Set parameters for ADOL-C of constitutive relations.

Parameters
tag[in] - tag of the tape
recalculate_elastic_tangent[out] - if setParam set it to true, tangent matrix for elastic domain should be recalculated

Reimplemented from ADOLCPlasticity::ClosestPointProjection.

Examples
ADOLCPlasticityMaterialModels.hpp.

Definition at line 836 of file ADOLCPlasticityMaterialModels.hpp.

836 {
838 switch (tag) {
839 case TypesTags::W: {
840 set_param_vec(tapesTags[tag], paramsArrayPtr->size(),
841 paramsArrayPtr->data());
842 for (auto p = 0; p != LAST_PARAM; p++) {
843 if ((*paramsArrayPtr)[p] != (*oldParamsArrayPtr)[p]) {
844 recalculate_elastic_tangent = true;
845 std::copy(paramsArrayPtr->begin(), paramsArrayPtr->end(),
846 oldParamsArrayPtr->begin());
847 break;
848 }
849 }
850
851 } break;
852 case TypesTags::Y:
853 case TypesTags::H:
854 break;
855 default:
856 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
857 "Unknown tag for setParams");
858 }
860 return 0;
861 }
std::array< int, LAST_TAPE > tapesTags

◆ yieldFunction()

MoFEMErrorCode ADOLCPlasticity::ParaboloidalPlasticity::yieldFunction ( )
inlinevirtual

Member Data Documentation

◆ defaultParamsArrayPtr

boost::shared_ptr<ParamsArray> ADOLCPlasticity::ParaboloidalPlasticity::defaultParamsArrayPtr = nullptr

◆ Hc

double ADOLCPlasticity::ParaboloidalPlasticity::Hc

◆ Ht

double ADOLCPlasticity::ParaboloidalPlasticity::Ht

◆ i

FTensor::Index<'i', 3> ADOLCPlasticity::ParaboloidalPlasticity::i

◆ I1

adouble ADOLCPlasticity::ParaboloidalPlasticity::I1

◆ j

FTensor::Index<'j', 3> ADOLCPlasticity::ParaboloidalPlasticity::j

◆ J2

adouble ADOLCPlasticity::ParaboloidalPlasticity::J2

◆ lambda

double ADOLCPlasticity::ParaboloidalPlasticity::lambda

◆ mu

double ADOLCPlasticity::ParaboloidalPlasticity::mu

◆ nc

double ADOLCPlasticity::ParaboloidalPlasticity::nc

◆ nt

double ADOLCPlasticity::ParaboloidalPlasticity::nt

◆ nup

double ADOLCPlasticity::ParaboloidalPlasticity::nup

◆ oldParamsArrayPtr

boost::shared_ptr<ParamsArray> ADOLCPlasticity::ParaboloidalPlasticity::oldParamsArrayPtr = nullptr

◆ paramsArrayPtr

boost::shared_ptr<ParamsArray> ADOLCPlasticity::ParaboloidalPlasticity::paramsArrayPtr = nullptr

◆ sigmaYc

double ADOLCPlasticity::ParaboloidalPlasticity::sigmaYc

◆ sigmaYt

double ADOLCPlasticity::ParaboloidalPlasticity::sigmaYt

◆ tDeviator

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tDeviator

◆ tElasticStrain

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tElasticStrain

◆ tPlasticStrain

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tPlasticStrain

◆ tR

adouble ADOLCPlasticity::ParaboloidalPlasticity::tR

◆ tStress

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tStress

◆ tTotalStrain

FTensor::Tensor2_symmetric<adouble, 3> ADOLCPlasticity::ParaboloidalPlasticity::tTotalStrain

◆ Z

FTensor::Index<'Z', 6> ADOLCPlasticity::ParaboloidalPlasticity::Z

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