|
OOFEM 3.0
|
#include <misesmat.h>
Public Member Functions | |
| MisesMat (int n, Domain *d) | |
| void | performPlasticityReturn (const FloatArray &totalStrain, GaussPoint *gp, TimeStep *tStep) const |
| void | performPlasticityReturn_PlaneStress (const FloatArrayF< 3 > &totalStrain, GaussPoint *gp, TimeStep *tStep) const |
| double | checkYieldStress (double &dKappa, double kappa, GaussPoint *gp, TimeStep *tStep) const |
| double | computeYieldStress (double kappa, GaussPoint *gp, TimeStep *tStep) const |
| double | computeYieldStressPrime (double kappa) const |
| double | computeDamage (GaussPoint *gp, TimeStep *tStep) const |
| double | computeDamageParam (double tempKappa) const |
| double | computeDamageParamPrime (double tempKappa) const |
| virtual double | computeCumPlastStrain (GaussPoint *gp, TimeStep *tStep) const |
| void | initializeFrom (InputRecord &ir) override |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| const char * | giveInputRecordName () const override |
| const char * | giveClassName () const override |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) 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. | |
| FloatArrayF< 3 > | giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &totalStrain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| FloatArrayF< 1 > | giveRealStressVector_1d (const FloatArrayF< 1 > &reducedE, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| double | giveS (GaussPoint *gp, TimeStep *tStep) const |
| double | giveTemperature (GaussPoint *gp, TimeStep *tStep) const |
| 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< 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) |
| virtual FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const |
| 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< 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 > | 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 double | predictRelativeComputationalCost (GaussPoint *gp) |
| virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
| 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 Member Functions | |
| void | computeGLPlasticStrain (const FloatMatrix &F, FloatMatrix &Ep, FloatMatrix b, double J) |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
Protected Attributes | |
| IsotropicLinearElasticMaterial | linearElasticMaterial |
| Reference to the basic elastic material. | |
| double | G = 0. |
| Elastic shear modulus. | |
| double | K = 0. |
| Elastic bulk modulus. | |
| double | H = 0. |
| Hardening modulus. | |
| ScalarFunction | sig0 |
| Initial (uniaxial) yield stress. | |
| int | hType |
| type of hardening function | |
| FloatArray | h_eps |
| user-defined hardening (yield stress - kappa) | |
| FloatArray | h_function_eps |
| double | omega_crit = 0. |
| critical(maximal) damage. | |
| double | a = 0. |
| exponent in damage function. | |
| double | yieldTol = 0. |
| tolerance for the yield function in RRM algorithm. | |
| 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 an isotropic elastoplastic material with Mises yield condition, associated flow rule and linear isotropic hardening.
It differs from other similar materials (such as J2Mat) by implementation - here we use the radial return, which is the most efficient algorithm for this specific model. Also, an extension to large strain will be available.
The model also exemplifies how to implement non-3d material modes, in this case 1D, by overloading the default implementations that iterates over the 3D method.
Definition at line 76 of file misesmat.h.
| oofem::MisesMat::MisesMat | ( | int | n, |
| Domain * | d ) |
Definition at line 54 of file misesmat.C.
References linearElasticMaterial, and oofem::StructuralMaterial::StructuralMaterial().
Referenced by giveRealStressVector_PlaneStress(), oofem::MisesMatGrad::MisesMatGrad(), and oofem::MisesMatNl::MisesMatNl().
| double oofem::MisesMat::checkYieldStress | ( | double & | dKappa, |
| double | kappa, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 379 of file misesmat.C.
References h_eps, h_function_eps, hType, and OOFEM_ERROR.
Referenced by performPlasticityReturn().
|
virtual |
Reimplemented in oofem::MisesMatGrad.
Definition at line 486 of file misesmat.C.
References oofem::Material::giveStatus(), and oofem::MisesMatStatus::giveTempCumulativePlasticStrain().
Referenced by computeDamage().
| double oofem::MisesMat::computeDamage | ( | GaussPoint * | gp, |
| TimeStep * | tStep ) const |
Definition at line 472 of file misesmat.C.
References computeCumPlastStrain(), computeDamageParam(), oofem::MisesMatStatus::giveDamage(), and oofem::Material::giveStatus().
Referenced by giveRealStressVector_1d(), giveRealStressVector_3d(), giveRealStressVector_PlaneStress(), and oofem::MisesMatGrad::giveRealStressVectorGradientDamage().
| double oofem::MisesMat::computeDamageParam | ( | double | tempKappa | ) | const |
Definition at line 451 of file misesmat.C.
References a, and omega_crit.
Referenced by computeDamage(), and oofem::MisesMatNl::computeDamage().
| double oofem::MisesMat::computeDamageParamPrime | ( | double | tempKappa | ) | const |
Definition at line 461 of file misesmat.C.
References a, and omega_crit.
Referenced by oofem::MisesMatGrad::give1dGprime(), give1dStressStiffMtrx(), oofem::MisesMatGrad::give1dStressStiffMtrx(), oofem::MisesMatNl::give1dStressStiffMtrx(), oofem::MisesMatGrad::give3dGprime(), give3dMaterialStiffnessMatrix(), oofem::MisesMatGrad::give3dMaterialStiffnessMatrix(), oofem::MisesMatNl::giveLocalNonlocalStiffnessContribution(), oofem::MisesMatGrad::givePlaneStrainGprime(), and oofem::MisesMatGrad::givePlaneStrainStiffMtrx().
|
protected |
| double oofem::MisesMat::computeYieldStress | ( | double | kappa, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 404 of file misesmat.C.
References giveS(), H, h_eps, h_function_eps, hType, and OOFEM_ERROR.
Referenced by give3dMaterialStiffnessMatrix(), performPlasticityReturn(), and performPlasticityReturn_PlaneStress().
| double oofem::MisesMat::computeYieldStressPrime | ( | double | kappa | ) | const |
Definition at line 426 of file misesmat.C.
References H, h_eps, h_function_eps, hType, and OOFEM_ERROR.
Referenced by give1dStressStiffMtrx(), give3dMaterialStiffnessMatrix(), givePlaneStressStiffMtrx(), performPlasticityReturn(), and performPlasticityReturn_PlaneStress().
|
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.
Reimplemented in oofem::MisesMatGrad, and oofem::MisesMatNl.
Definition at line 103 of file misesmat.C.
|
overridevirtual |
Method for computing 1d stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Reimplemented in oofem::MisesMatGrad, and oofem::MisesMatNl.
Definition at line 616 of file misesmat.C.
References oofem::FloatArray::at(), computeDamageParamPrime(), computeYieldStressPrime(), E, oofem::MisesMatStatus::giveCumulativePlasticStrain(), oofem::Material::giveStatus(), oofem::MisesMatStatus::giveTempCumulativePlasticStrain(), oofem::MisesMatStatus::giveTempDamage(), oofem::MisesMatStatus::giveTempEffectiveStress(), linearElasticMaterial, and oofem::signum().
|
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.
Reimplemented in oofem::MisesMatGrad.
Definition at line 496 of file misesmat.C.
References computeDamageParamPrime(), oofem::StructuralMaterial::computeStressNorm(), computeYieldStress(), computeYieldStressPrime(), oofem::dyad(), G, oofem::MisesMatStatus::giveCumulativePlasticStrain(), oofem::Material::giveStatus(), oofem::I_dev6, and linearElasticMaterial.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Reimplemented in oofem::MisesMatGrad, and oofem::MisesMatNl.
Definition at line 128 of file misesmat.h.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Reimplemented in oofem::MisesMatGrad, and oofem::MisesMatNl.
Definition at line 127 of file misesmat.h.
References _IFT_MisesMat_Name.
|
overrideprotectedvirtual |
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 647 of file misesmat.C.
References oofem::FloatArray::at(), oofem::MisesMatStatus::giveCumulativePlasticStrain(), oofem::MisesMatStatus::giveDamage(), oofem::StructuralMaterial::giveIPValue(), oofem::MisesMatStatus::givePlasDef(), giveS(), oofem::Material::giveStatus(), and oofem::FloatArray::resize().
|
overridevirtual |
Method for computing plane stress stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane stress stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.
| answer | Stiffness matrix. |
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 548 of file misesmat.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrixF< N, M >::at(), oofem::FloatMatrix::beDyadicProductOf(), oofem::FloatArray::beProductOf(), computeYieldStressPrime(), E, G, oofem::StructuralMaterial::giveReducedSymVectorForm(), oofem::Material::giveStatus(), linearElasticMaterial, and oofem::FloatMatrix::times().
|
overridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Reimplemented in oofem::MisesMatNl.
Definition at line 109 of file misesmat.C.
References computeDamage(), E, oofem::StructuralMaterial::giveRealStressVector_1d(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), linearElasticMaterial, and performPlasticityReturn().
|
overridevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Reimplemented in oofem::MisesMatNl.
Definition at line 162 of file misesmat.C.
References computeDamage(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), and performPlasticityReturn().
|
overridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 145 of file misesmat.C.
References computeDamage(), oofem::Material::giveStatus(), oofem::Material::initTempStatus(), MisesMat(), and performPlasticityReturn_PlaneStress().
| double oofem::MisesMat::giveS | ( | GaussPoint * | gp, |
| TimeStep * | tStep ) const |
Definition at line 732 of file misesmat.C.
References oofem::FEMComponent::giveDomain(), oofem::TimeStep::giveIntrinsicTime(), giveTemperature(), and sig0.
Referenced by computeYieldStress(), oofem::MisesMatGrad::give3dMaterialStiffnessMatrix(), giveIPValue(), and oofem::MisesMatGrad::givePlaneStrainStiffMtrx().
| double oofem::MisesMat::giveTemperature | ( | GaussPoint * | gp, |
| TimeStep * | tStep ) const |
Definition at line 739 of file misesmat.C.
References oofem::FloatArray::at(), oofem::Element::computeGlobalCoordinates(), oofem::FEMComponent::domain, oofem::GaussPoint::giveElement(), oofem::FieldManager::giveField(), oofem::GaussPoint::giveNaturalCoordinates(), oofem::FEMComponent::giveNumber(), and OOFEM_ERROR.
Referenced by giveS().
|
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::Material.
Reimplemented in oofem::MisesMatGrad, and oofem::MisesMatNl.
Definition at line 60 of file misesmat.C.
References _IFT_MisesMat_a, _IFT_MisesMat_h, _IFT_MisesMat_h_eps, _IFT_MisesMat_h_function_eps, _IFT_MisesMat_htype, _IFT_MisesMat_omega_crit, _IFT_MisesMat_sig0, _IFT_MisesMat_yieldTol, a, G, H, h_eps, h_function_eps, hType, oofem::StructuralMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, K, linearElasticMaterial, omega_crit, sig0, and yieldTol.
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 125 of file misesmat.h.
| void oofem::MisesMat::performPlasticityReturn | ( | const FloatArray & | totalStrain, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 181 of file misesmat.C.
References oofem::FloatArray::add(), oofem::StructuralMaterial::applyDeviatoricElasticCompliance(), oofem::StructuralMaterial::applyDeviatoricElasticStiffness(), oofem::FloatArray::at(), checkYieldStress(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeDeviatoricVolumetricSum(), oofem::StructuralMaterial::computeStressNorm(), computeYieldStress(), computeYieldStressPrime(), E, G, oofem::FloatArray::giveSize(), oofem::Material::giveStatus(), hType, K, linearElasticMaterial, oofem::FloatArray::resize(), oofem::signum(), oofem::FloatArray::subtract(), and oofem::FloatArray::times().
Referenced by giveRealStressVector_1d(), oofem::MisesMatNl::giveRealStressVector_1d(), giveRealStressVector_3d(), and oofem::MisesMatNl::updateBeforeNonlocAverage().
| void oofem::MisesMat::performPlasticityReturn_PlaneStress | ( | const FloatArrayF< 3 > & | totalStrain, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 266 of file misesmat.C.
References oofem::FloatArray::at(), oofem::FloatArrayF< N >::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beProductOf(), computeYieldStress(), computeYieldStressPrime(), E, G, oofem::MisesMatStatus::giveCumulativePlasticStrain(), oofem::StructuralMaterial::giveFullSymVectorForm(), oofem::MisesMatStatus::givePlasticStrain(), oofem::StructuralMaterial::giveReducedSymVectorForm(), oofem::Material::giveStatus(), K, oofem::MisesMatStatus::letTempEffectiveStressBe(), oofem::MisesMatStatus::letTempPlasticStrainBe(), linearElasticMaterial, OOFEM_WARNING, oofem::MisesMatStatus::setTempCumulativePlasticStrain(), oofem::FloatArray::subtract(), oofem::FloatArray::times(), and yieldTol.
Referenced by giveRealStressVector_PlaneStress().
|
protected |
exponent in damage function.
Definition at line 103 of file misesmat.h.
Referenced by computeDamageParam(), computeDamageParamPrime(), and initializeFrom().
|
protected |
Elastic shear modulus.
Definition at line 83 of file misesmat.h.
Referenced by oofem::MisesMatGrad::give3dKappaMatrix(), give3dMaterialStiffnessMatrix(), oofem::MisesMatGrad::give3dMaterialStiffnessMatrix(), oofem::MisesMatGrad::givePlaneStrainKappaMatrix(), oofem::MisesMatGrad::givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), initializeFrom(), performPlasticityReturn(), and performPlasticityReturn_PlaneStress().
|
protected |
Hardening modulus.
Definition at line 89 of file misesmat.h.
Referenced by computeYieldStress(), computeYieldStressPrime(), oofem::MisesMatGrad::give1dKappaMatrix(), oofem::MisesMatGrad::give1dStressStiffMtrx(), oofem::MisesMatNl::give1dStressStiffMtrx(), oofem::MisesMatGrad::give3dKappaMatrix(), oofem::MisesMatGrad::give3dMaterialStiffnessMatrix(), oofem::MisesMatGrad::givePlaneStrainKappaMatrix(), oofem::MisesMatGrad::givePlaneStrainStiffMtrx(), oofem::MisesMatNl::giveRemoteNonlocalStiffnessContribution(), and initializeFrom().
|
protected |
user-defined hardening (yield stress - kappa)
Definition at line 98 of file misesmat.h.
Referenced by checkYieldStress(), computeYieldStress(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Definition at line 98 of file misesmat.h.
Referenced by checkYieldStress(), computeYieldStress(), computeYieldStressPrime(), and initializeFrom().
|
protected |
type of hardening function
Definition at line 95 of file misesmat.h.
Referenced by checkYieldStress(), computeYieldStress(), computeYieldStressPrime(), initializeFrom(), and performPlasticityReturn().
|
protected |
Elastic bulk modulus.
Definition at line 86 of file misesmat.h.
Referenced by initializeFrom(), performPlasticityReturn(), and performPlasticityReturn_PlaneStress().
|
protected |
Reference to the basic elastic material.
Definition at line 80 of file misesmat.h.
Referenced by oofem::MisesMatGrad::give1dKappaMatrix(), give1dStressStiffMtrx(), oofem::MisesMatGrad::give1dStressStiffMtrx(), oofem::MisesMatNl::give1dStressStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::MisesMatGrad::give3dMaterialStiffnessMatrix(), oofem::MisesMatGrad::givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), giveRealStressVector_1d(), oofem::MisesMatNl::giveRemoteNonlocalStiffnessContribution(), initializeFrom(), MisesMat(), performPlasticityReturn(), and performPlasticityReturn_PlaneStress().
|
protected |
critical(maximal) damage.
Definition at line 101 of file misesmat.h.
Referenced by computeDamageParam(), computeDamageParamPrime(), and initializeFrom().
|
protected |
Initial (uniaxial) yield stress.
Definition at line 92 of file misesmat.h.
Referenced by giveS(), and initializeFrom().
|
protected |
tolerance for the yield function in RRM algorithm.
Definition at line 106 of file misesmat.h.
Referenced by initializeFrom(), and performPlasticityReturn_PlaneStress().