|
OOFEM 3.0
|
#include <druckerPragerPlasticitySM.h>
Public Types | |
| enum | state_flag_values { DP_Elastic , DP_Unloading , DP_Yielding , DP_Vertex } |
| Values of history variable state_flag. More... | |
Public Member Functions | |
| DruckerPragerPlasticitySM (int n, Domain *d) | |
| Constructor. | |
| void | initializeFrom (InputRecord &ir) override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| FloatArrayF< 6 > | giveRealStressVector_3d (const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| void | performLocalStressReturn (GaussPoint *gp, const FloatArrayF< 6 > &strain) const |
| bool | checkForVertexCase (double eM, double gM, double kM, double trialStressJTwo, double volumetricStress, double tempKappa) const |
| void | performRegularReturn (double eM, double gM, double kM, double trialStressJTwo, FloatArrayF< 6 > &stressDeviator, double &volumetricStress, double &tempKappa) const |
| void | performVertexReturn (double eM, double gM, double kM, double trialStressJTwo, FloatArrayF< 6 > &stressDeviator, double &volumetricStress, double &tempKappa, double volumetricElasticTrialStrain, double kappa) const |
| double | computeYieldValue (double meanStress, double JTwo, double kappa, double eM) const |
| virtual double | computeYieldStressInShear (double kappa, double eM) const |
| virtual double | computeYieldStressPrime (double kappa, double eM) const |
| FloatMatrixF< 6, 6 > | giveRegAlgorithmicStiffMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| FloatMatrixF< 6, 6 > | giveVertexAlgorithmicStiffMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const override |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
| double | predictRelativeComputationalCost (GaussPoint *gp) override |
| double | predictRelativeRedistributionCost (GaussPoint *gp) override |
| Public Member Functions inherited from oofem::StructuralMaterial | |
| StructuralMaterial (int n, Domain *d) | |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| const char * | giveClassName () const override |
| void | initializeFrom (InputRecord &ir) override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| void | giveCharacteristicMatrix (FloatMatrix &answer, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override |
| Returns characteristic matrix of the receiver. | |
| void | giveCharacteristicVector (FloatArray &answer, FloatArray &flux, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override |
| Returns characteristic vector of the receiver. | |
| virtual void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) const |
| virtual FloatArrayF< 4 > | giveRealStressVector_PlaneStrain (const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_3d. | |
| virtual FloatArray | giveRealStressVector_StressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| virtual FloatArray | giveRealStressVector_ShellStressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 3 > | giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 1 > | giveRealStressVector_1d (const FloatArrayF< 1 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 2 > | giveRealStressVector_Warping (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 2 > | giveRealStressVector_2dBeamLayer (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 5 > | giveRealStressVector_PlateLayer (const FloatArrayF< 5 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 3 > | giveRealStressVector_Fiber (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 3 > | giveRealStressVector_2dPlateSubSoil (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation is not provided. | |
| virtual FloatArrayF< 6 > | giveRealStressVector_3dBeamSubSoil (const FloatArrayF< 6 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 9 > | giveFirstPKStressVector_3d (const FloatArrayF< 9 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| virtual FloatArrayF< 5 > | giveFirstPKStressVector_PlaneStrain (const FloatArrayF< 5 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_3d. | |
| virtual FloatArray | giveFirstPKStressVector_StressControl (const FloatArray &reducedvF, const IntArray &FControl, GaussPoint *gp, TimeStep *tStep) const |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| virtual FloatArrayF< 4 > | giveFirstPKStressVector_PlaneStress (const FloatArrayF< 4 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_StressControl. | |
| virtual FloatArrayF< 1 > | giveFirstPKStressVector_1d (const FloatArrayF< 1 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_StressControl. | |
| virtual void | giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| double | giveReferenceTemperature () |
| virtual FloatArray | computeStressIndependentStrainVector (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const |
| FloatArrayF< 6 > | computeStressIndependentStrainVector_3d (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const |
| virtual void | giveStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 9, 9 > | give3dMaterialStiffnessMatrix_dPdF (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| void | giveStressDependentPartOfStrainVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) const |
| int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| virtual FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 4, 4 > | givePlaneStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 4, 4 > | givePlaneStrainStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | givePlaneStrainStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 1, 1 > | give1dStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 2, 2 > | give2dBeamLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | givePlateLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 3, 3 > | giveFiberStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 3, 3 > | give2dPlateSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 6, 6 > | give3dBeamSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| Public Member Functions inherited from oofem::Material | |
| Material (int n, Domain *d) | |
| virtual | ~Material ()=default |
| Destructor. | |
| virtual double | giveCharacteristicValue (MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const |
| Returns characteristic value of the receiver. | |
| virtual double | give (int aProperty, GaussPoint *gp) const |
| virtual bool | hasProperty (int aProperty, GaussPoint *gp) const |
| virtual void | modifyProperty (int aProperty, double value, GaussPoint *gp) |
| double | giveCastingTime () const |
| virtual bool | isActivated (TimeStep *tStep) const |
| virtual bool | hasCastingTimeSupport () const |
| void | printYourself () override |
| Prints receiver state on stdout. Useful for debugging. | |
| virtual void | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| virtual void | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| int | checkConsistency () override |
| virtual void | restoreConsistency (GaussPoint *gp) |
| virtual int | initMaterial (Element *element) |
| virtual MaterialStatus * | giveStatus (GaussPoint *gp) const |
| virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
| virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
| virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip) |
| virtual void | initTempStatus (GaussPoint *gp) const |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| Public Member Functions inherited from oofem::FEMComponent | |
| FEMComponent (int n, Domain *d) | |
| virtual | ~FEMComponent ()=default |
| Virtual destructor. | |
| Domain * | giveDomain () const |
| virtual void | setDomain (Domain *d) |
| int | giveNumber () const |
| void | setNumber (int num) |
| virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
| virtual void | initializeFrom (InputRecord &ir, int priority) |
| virtual void | initializeFinish () |
| virtual void | postInitialize () |
| Performs post initialization steps. Called after all components are created and initialized. | |
| virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
| virtual Interface * | giveInterface (InterfaceType t) |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
Protected Attributes | |
| int | hardeningType = 1 |
| double | kappaC = 0. |
| Parameter of the exponential laws. | |
| double | hardeningModulus = 0. |
| Hardening modulus normalized with the elastic modulus, parameter of the linear hardening/softening law. | |
| double | limitYieldStress = 0. |
| Parameter of the exponential hardening law. | |
| double | initialYieldStress = 0. |
| Parameter of all three laws, this is the initial value of the yield stress in pure shear. | |
| double | alpha = 0. |
| Friction coefficient, parameter of the yield criterion. | |
| double | alphaPsi = 0. |
| Dilatancy coefficient, parameter of the flow rule. | |
| double | kFactor = 0. |
| Scalar factor between rate of plastic multiplier and rate of hardening variable. | |
| IsotropicLinearElasticMaterial | LEMaterial |
| Associated linear elastic material. | |
| double | yieldTol = 0. |
| Yield tolerance. | |
| int | newtonIter = 0 |
| Maximum number of iterations for stress return. | |
| Protected Attributes inherited from oofem::StructuralMaterial | |
| double | referenceTemperature = 0. |
| Reference temperature (temperature, when material has been built into structure). | |
| MatResponseMode | SCStiffMode = TangentStiffness |
| stifness mode used in stress control | |
| double | SCRelTol = 1.e-3 |
| relative tolerance for stress control | |
| double | SCAbsTol = 1.e-12 |
| absolute stress tolerance for stress control | |
| int | SCMaxiter = 100000 |
| maximum iterations for stress-control | |
| Protected Attributes inherited from oofem::Material | |
| Dictionary | propertyDictionary |
| double | castingTime |
| int | preCastingTimeMat |
| Material existing before casting time - optional parameter, zero by default. | |
| Protected Attributes inherited from oofem::FEMComponent | |
| int | number |
| Component number. | |
| Domain * | domain |
| Link to domain object, useful for communicating with other FEM components. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from oofem::StructuralMaterial | |
| static int | giveSymVI (int ind1, int ind2) |
| static int | giveVI (int ind1, int ind2) |
| static FloatMatrixF< 9, 9 > | convert_dSdE_2_dPdF_3D (const FloatMatrixF< 6, 6 > &dSdE, const FloatArrayF< 6 > &S, const FloatArrayF< 9 > &F) |
| static FloatMatrixF< 5, 5 > | convert_dSdE_2_dPdF_PlaneStrain (const FloatMatrixF< 4, 4 > &dSdE, const FloatArrayF< 4 > &S, const FloatArrayF< 5 > &F) |
| static FloatMatrixF< 4, 4 > | convert_dSdE_2_dPdF_PlaneStress (const FloatMatrixF< 3, 3 > &dSdE, const FloatArrayF< 3 > &S, const FloatArrayF< 4 > &F) |
| static FloatMatrixF< 1, 1 > | convert_dSdE_2_dPdF_1D (const FloatMatrixF< 1, 1 > &dSdE, const FloatArrayF< 1 > &S, const FloatArrayF< 1 > &F) |
| static void | computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode) |
| Common functions for convenience. | |
| static FloatArrayF< 3 > | computePrincipalValues (const FloatMatrixF< 3, 3 > &s) |
| static FloatArrayF< 3 > | computePrincipalValues (double I1, double I2, double I3) |
| static void | computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode) |
| static std::pair< FloatArrayF< 3 >, FloatMatrixF< 3, 3 > > | computePrincipalValDir (const FloatMatrixF< 3, 3 > &s) |
| static FloatArrayF< 6 > | computeDeviator (const FloatArrayF< 6 > &s) |
| static std::pair< FloatArrayF< 6 >, double > | computeDeviatoricVolumetricSplit (const FloatArrayF< 6 > &s) |
| static FloatArrayF< 6 > | computeDeviatoricVolumetricSum (const FloatArrayF< 6 > &dev, double mean) |
| static FloatArrayF< 6 > | applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double GModulus) |
| static FloatArrayF< 6 > | applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double GModulus) |
| static FloatArrayF< 6 > | applyElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu) |
| static double | computeStressNorm (const FloatArrayF< 6 > &stress) |
| static double | computeFirstInvariant (const FloatArrayF< 6 > &s) |
| static double | computeSecondStressInvariant (const FloatArrayF< 6 > &s) |
| static double | computeThirdStressInvariant (const FloatArrayF< 6 > &s) |
| static double | computeFirstCoordinate (const FloatArrayF< 6 > &s) |
| static double | computeSecondCoordinate (const FloatArrayF< 6 > &s) |
| static double | computeThirdCoordinate (const FloatArrayF< 6 > &s) |
| static int | giveVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| static int | giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode) |
| static void | giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| static int | giveSizeOfVoigtVector (MaterialMode mmode) |
| static int | giveSizeOfVoigtSymVector (MaterialMode mmode) |
| static void | giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. | |
| static void | giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced deformation gradient Voigt vector (2nd order tensor). | |
| static void | giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. | |
| static void | giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. | |
| static void | giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. | |
| static void | giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. | |
| static void | giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. | |
| static void | giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. | |
| static FloatArrayF< 6 > | transformStrainVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &strain, bool transpose=false) |
| static FloatArrayF< 6 > | transformStressVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &stress, bool transpose=false) |
| static double | computeVonMisesStress (const FloatArray ¤tStress) |
| static double | computeVonMisesStress_3D (const FloatArrayF< 6 > &stress) |
| static double | computeVonMisesStress_PlaneStress (const FloatArrayF< 3 > &stress) |
| static FloatMatrixF< 6, 6 > | giveStrainVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false) |
| static FloatMatrixF< 3, 3 > | give2DStrainVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false) |
| static FloatMatrixF< 6, 6 > | giveStressVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false) |
| static FloatMatrixF< 3, 3 > | givePlaneStressVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false) |
| static void | sortPrincDirAndValCloseTo (FloatArray &pVal, FloatMatrix &pDir, const FloatMatrix &toPDir) |
| Static Public Attributes inherited from oofem::StructuralMaterial | |
| static std::array< std::array< int, 3 >, 3 > | vIindex |
| Voigt index map. | |
| static std::array< std::array< int, 3 >, 3 > | svIndex |
| Symmetric Voigt index map. | |
This class implements a (local) nonassociated plasticity model based on the Drucker-Prager yield criterion with hardening and softening.
Definition at line 193 of file druckerPragerPlasticitySM.h.
Values of history variable state_flag.
| Enumerator | |
|---|---|
| DP_Elastic | |
| DP_Unloading | |
| DP_Yielding | |
| DP_Vertex | |
Definition at line 197 of file druckerPragerPlasticitySM.h.
| oofem::DruckerPragerPlasticitySM::DruckerPragerPlasticitySM | ( | int | n, |
| Domain * | d ) |
Constructor.
Definition at line 177 of file druckerPragerPlasticitySM.C.
References LEMaterial, and oofem::StructuralMaterial::StructuralMaterial().
Referenced by giveRealStressVector_3d().
| bool oofem::DruckerPragerPlasticitySM::checkForVertexCase | ( | double | eM, |
| double | gM, | ||
| double | kM, | ||
| double | trialStressJTwo, | ||
| double | volumetricStress, | ||
| double | tempKappa ) const |
Check if the trial stress state falls within the vertex region.
| eM | Elasticity modulus. |
| gM | Shear modulus. |
| kM | Bulk modulus. |
Definition at line 318 of file druckerPragerPlasticitySM.C.
References alphaPsi, computeYieldValue(), kFactor, and yieldTol.
Referenced by performLocalStressReturn().
|
virtual |
Compute the current yield stress in pure shear of the Drucker-Prager model according to the used hardening law. The yield stress is tauY in f(sigma, kappa) = F(sigma) - tauY(kappa).
| kappa | Hardening variable. |
| eM | Elasticity modulus. |
Definition at line 465 of file druckerPragerPlasticitySM.C.
References hardeningModulus, hardeningType, initialYieldStress, kappaC, limitYieldStress, and OOFEM_ERROR.
Referenced by computeYieldValue().
|
virtual |
Compute derivative of yield stress with respect to the hardening variable kappa.
| kappa | Hardening variable. |
| eM | Elasticity modulus. |
Definition at line 485 of file druckerPragerPlasticitySM.C.
References hardeningModulus, hardeningType, initialYieldStress, kappaC, limitYieldStress, and OOFEM_ERROR.
Referenced by giveRegAlgorithmicStiffMatrix(), giveVertexAlgorithmicStiffMatrix(), performRegularReturn(), and performVertexReturn().
| double oofem::DruckerPragerPlasticitySM::computeYieldValue | ( | double | meanStress, |
| double | JTwo, | ||
| double | kappa, | ||
| double | eM ) const |
Compute the yield value based on stress and hardening variable.
| meanStress | 1/3 of trace of sigma. |
| JTwo | Second deviatoric invariant. |
| kappa | Hardening variable. |
| eM | Elasticity modulus. |
Definition at line 456 of file druckerPragerPlasticitySM.C.
References alpha, and computeYieldStressInShear().
Referenced by checkForVertexCase(), performLocalStressReturn(), performRegularReturn(), and performVertexReturn().
|
overridevirtual |
Creates new copy of associated status and inserts it into given integration point.
| gp | Integration point where newly created status will be stored. |
Reimplemented from oofem::Material.
Definition at line 709 of file druckerPragerPlasticitySM.C.
|
overridevirtual |
Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.
| answer | Computed results. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 505 of file druckerPragerPlasticitySM.C.
References giveRegAlgorithmicStiffMatrix(), oofem::Material::giveStatus(), giveVertexAlgorithmicStiffMatrix(), LEMaterial, and OOFEM_ERROR.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 235 of file druckerPragerPlasticitySM.h.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 236 of file druckerPragerPlasticitySM.h.
References _IFT_DruckerPragerPlasticitySM_Name.
|
overridevirtual |
Returns the integration point corresponding value in Reduced form.
| answer | Contain corresponding ip value, zero sized if not available. |
| gp | Integration point to which the value refers. |
| type | Determines the type of internal variable. |
| tStep | Determines the time step. |
Reimplemented from oofem::Material.
Definition at line 685 of file druckerPragerPlasticitySM.C.
References oofem::FloatArray::at(), oofem::Material::giveStatus(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 226 of file druckerPragerPlasticitySM.C.
References oofem::StructuralMaterial::computeStressIndependentStrainVector_3d(), DruckerPragerPlasticitySM(), oofem::Material::giveStatus(), oofem::Material::initTempStatus(), and performLocalStressReturn().
| FloatMatrixF< 6, 6 > oofem::DruckerPragerPlasticitySM::giveRegAlgorithmicStiffMatrix | ( | MatResponseMode | mode, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Compute and give back algorithmic stiffness matrix for the regular case (no vertex).
| mode | Material reponse mode. |
| gp | Gauss point. |
| tStep | Time step. |
Definition at line 543 of file druckerPragerPlasticitySM.C.
References alpha, alphaPsi, oofem::FloatMatrixF< N, M >::at(), oofem::StructuralMaterial::computeDeviator(), oofem::StructuralMaterial::computeSecondStressInvariant(), computeYieldStressPrime(), oofem::dot(), Ex, oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveTempStressVector(), oofem::inv(), kFactor, LEMaterial, NYxz, and OOFEM_ERROR.
Referenced by give3dMaterialStiffnessMatrix().
|
inlineoverridevirtual |
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis.
| answer | Vector of thermal dilatation coefficients. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 324 of file druckerPragerPlasticitySM.h.
References LEMaterial.
| FloatMatrixF< 6, 6 > oofem::DruckerPragerPlasticitySM::giveVertexAlgorithmicStiffMatrix | ( | MatResponseMode | mode, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Compute consistent stiffness matrix for the vertex case.
| mode | Material reponse mode. |
| gp | Gauss point. |
| tStep | Time step. |
Definition at line 625 of file druckerPragerPlasticitySM.C.
References alpha, oofem::FloatMatrixF< N, M >::at(), oofem::StructuralMaterial::computeDeviator(), computeYieldStressPrime(), Ex, oofem::FEMComponent::giveNumber(), oofem::Material::giveStatus(), oofem::DruckerPragerPlasticitySMStatus::giveTempKappa(), LEMaterial, NYxz, and OOFEM_ERROR.
Referenced by give3dMaterialStiffnessMatrix().
|
overridevirtual |
Initializes receiver according to object description stored in input record. This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record. Note that initializeFrom may be called mutiple times.
| ir | Input record to initialize from. |
| priority | Priority of the input record. This is used to determine the order of initialization |
Reimplemented from oofem::FEMComponent.
Definition at line 184 of file druckerPragerPlasticitySM.C.
References _IFT_DruckerPragerPlasticitySM_alpha, _IFT_DruckerPragerPlasticitySM_alphapsi, _IFT_DruckerPragerPlasticitySM_hm, _IFT_DruckerPragerPlasticitySM_ht, _IFT_DruckerPragerPlasticitySM_iys, _IFT_DruckerPragerPlasticitySM_kc, _IFT_DruckerPragerPlasticitySM_lys, _IFT_DruckerPragerPlasticitySM_newtoniter, _IFT_DruckerPragerPlasticitySM_yieldtol, alpha, alphaPsi, hardeningModulus, hardeningType, initialYieldStress, IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, kappaC, kFactor, LEMaterial, limitYieldStress, newtonIter, and yieldTol.
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 322 of file druckerPragerPlasticitySM.h.
| void oofem::DruckerPragerPlasticitySM::performLocalStressReturn | ( | GaussPoint * | gp, |
| const FloatArrayF< 6 > & | strain ) const |
Perform a standard local stress return using the function computeYieldValue at the specified Gauss point. This function computes the history variables, i.e. the temp variable of plastic strain, hardening variable, state flag, and the temp stress, and stores them in the temp-status.
| gp | Gauss point. |
| strain | Strain vector of this Gauss point. |
Definition at line 251 of file druckerPragerPlasticitySM.C.
References oofem::StructuralMaterial::applyDeviatoricElasticCompliance(), oofem::StructuralMaterial::applyDeviatoricElasticStiffness(), checkForVertexCase(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeDeviatoricVolumetricSum(), oofem::StructuralMaterial::computeSecondStressInvariant(), computeYieldValue(), Ex, oofem::Material::giveStatus(), LEMaterial, NYxz, performRegularReturn(), performVertexReturn(), and yieldTol.
Referenced by giveRealStressVector_3d().
| void oofem::DruckerPragerPlasticitySM::performRegularReturn | ( | double | eM, |
| double | gM, | ||
| double | kM, | ||
| double | trialStressJTwo, | ||
| FloatArrayF< 6 > & | stressDeviator, | ||
| double & | volumetricStress, | ||
| double & | tempKappa ) const |
Perform stress return for regular case, i.e. if the trial stress state does not lie within the vertex region.
| eM | Elasticity modulus. |
| gM | Shear modulus. |
| kM | Bulk modulus. |
Definition at line 340 of file druckerPragerPlasticitySM.C.
References alpha, alphaPsi, oofem::StructuralMaterial::computeSecondStressInvariant(), computeYieldStressPrime(), computeYieldValue(), kFactor, newtonIter, OOFEM_ERROR, OOFEM_LOG_DEBUG, and yieldTol.
Referenced by performLocalStressReturn().
| void oofem::DruckerPragerPlasticitySM::performVertexReturn | ( | double | eM, |
| double | gM, | ||
| double | kM, | ||
| double | trialStressJTwo, | ||
| FloatArrayF< 6 > & | stressDeviator, | ||
| double & | volumetricStress, | ||
| double & | tempKappa, | ||
| double | volumetricElasticTrialStrain, | ||
| double | kappa ) const |
Perform stress return for vertex case, i.e. if the trial stress state lies within the vertex region.
| eM | Elasticity modulus. |
| gM | Shear modulus. |
| kM | Bulk modulus. |
Definition at line 400 of file druckerPragerPlasticitySM.C.
References alpha, computeYieldStressPrime(), computeYieldValue(), newtonIter, OOFEM_ERROR, OOFEM_LOG_DEBUG, yieldTol, and oofem::zeros().
Referenced by performLocalStressReturn().
|
overridevirtual |
Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model.
Reimplemented from oofem::Material.
Definition at line 716 of file druckerPragerPlasticitySM.C.
References oofem::DruckerPragerPlasticitySMStatus::giveStateFlag(), and oofem::Material::giveStatus().
|
inlineoverridevirtual |
Returns the relative redistribution cost of the receiver
Reimplemented from oofem::Material.
Definition at line 332 of file druckerPragerPlasticitySM.h.
|
protected |
Friction coefficient, parameter of the yield criterion.
Definition at line 215 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldValue(), giveRegAlgorithmicStiffMatrix(), giveVertexAlgorithmicStiffMatrix(), initializeFrom(), performRegularReturn(), and performVertexReturn().
|
protected |
Dilatancy coefficient, parameter of the flow rule.
Definition at line 217 of file druckerPragerPlasticitySM.h.
Referenced by checkForVertexCase(), giveRegAlgorithmicStiffMatrix(), initializeFrom(), and performRegularReturn().
|
protected |
Hardening modulus normalized with the elastic modulus, parameter of the linear hardening/softening law.
Definition at line 209 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Controls the hardening function in the yield stress: 1: linear hardening/softening with cutoff at zero stress. 2: exponential hardening/softening to limitYieldStress.
Definition at line 205 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Parameter of all three laws, this is the initial value of the yield stress in pure shear.
Definition at line 213 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Parameter of the exponential laws.
Definition at line 207 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Scalar factor between rate of plastic multiplier and rate of hardening variable.
Definition at line 219 of file druckerPragerPlasticitySM.h.
Referenced by checkForVertexCase(), giveRegAlgorithmicStiffMatrix(), initializeFrom(), and performRegularReturn().
|
protected |
Associated linear elastic material.
Definition at line 222 of file druckerPragerPlasticitySM.h.
Referenced by DruckerPragerPlasticitySM(), give3dMaterialStiffnessMatrix(), giveRegAlgorithmicStiffMatrix(), giveThermalDilatationVector(), giveVertexAlgorithmicStiffMatrix(), initializeFrom(), and performLocalStressReturn().
|
protected |
Parameter of the exponential hardening law.
Definition at line 211 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Maximum number of iterations for stress return.
Definition at line 227 of file druckerPragerPlasticitySM.h.
Referenced by initializeFrom(), performRegularReturn(), and performVertexReturn().
|
protected |
Yield tolerance.
Definition at line 225 of file druckerPragerPlasticitySM.h.
Referenced by checkForVertexCase(), initializeFrom(), performLocalStressReturn(), performRegularReturn(), and performVertexReturn().