|
OOFEM 3.0
|
#include <mplasticmaterial2.h>
Public Member Functions | |
| MPlasticMaterial2 (int n, Domain *d) | |
| virtual | ~MPlasticMaterial2 () |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| const char * | giveClassName () const override |
| LinearElasticMaterial * | giveLinearElasticMaterial () |
| Returns reference to undamaged (bulk) material. | |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &strain, 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< 4 > | giveRealStressVector_PlaneStrain (const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_3d. | |
| FloatArrayF< 3 > | giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| FloatArrayF< 1 > | giveRealStressVector_1d (const FloatArrayF< 1 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual double | computeDamage (GaussPoint *gp, const FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| virtual int | giveSizeOfFullHardeningVarsVector () const |
| virtual int | giveSizeOfReducedHardeningVarsVector (GaussPoint *) const |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const 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 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) |
| 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< 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< 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. | |
| virtual const char * | giveInputRecordName () const =0 |
| 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 Types | |
| enum | ReturnMappingAlgoType { mpm_ClosestPoint , mpm_CuttingPlane } |
| Protected type to determine the return mapping algorithm. More... | |
| enum | functType { yieldFunction , loadFunction } |
| Type that allows to distinguish between yield function and loading function. More... | |
| enum | plastType { associatedPT , nonassociatedPT } |
Protected Member Functions | |
| virtual int | giveMaxNumberOfActiveYieldConds (GaussPoint *gp) const =0 |
| void | closestPointReturn (FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma, GaussPoint *gp, const FloatArray &totalStrain, FloatArray &plasticStrainR, FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const |
| void | cuttingPlaneReturn (FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma, GaussPoint *gp, const FloatArray &totalStrain, FloatArray &plasticStrainR, FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) const |
| void | computeResidualVector (FloatArray &answer, GaussPoint *gp, const FloatArray &gamma, const IntArray &activeConditionMap, const FloatArray &plasticStrainVectorR, const FloatArray &strainSpaceHardeningVariables, std ::vector< FloatArray > &gradVec) const |
| virtual FloatMatrix | giveConsistentStiffnessMatrix (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrix | giveElastoPlasticStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| void | computeAlgorithmicModuli (FloatMatrix &answer, GaussPoint *gp, const FloatMatrix &elasticModuliInverse, const FloatArray &gamma, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const |
| virtual double | computeYieldValueAt (GaussPoint *gp, int isurf, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const =0 |
| Computes the value of yield function. | |
| virtual void | computeStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const =0 |
| Computes the stress gradient of yield/loading function (df/d_sigma). | |
| void | computeReducedStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const |
| virtual void | computeStrainHardeningVarsIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &stress, const FloatArray &dlambda, const FloatArray &dplasticStrain, const IntArray &activeConditionMap) const =0 |
| virtual void | computeKGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const =0 |
| virtual void | computeReducedHardeningVarsSigmaGradient (FloatMatrix &answer, GaussPoint *gp, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &gamma) const =0 |
| virtual void | computeReducedHardeningVarsLamGradient (FloatMatrix &answer, GaussPoint *gp, int actSurf, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &gamma) const =0 |
| computes derivative of \( \kappa \) vector with respect to lambda vector | |
| virtual int | hasHardening () const =0 |
| virtual void | computeReducedSSGradientMatrix (FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const =0 |
| Computes second derivative of loading function with respect to stress. | |
| virtual void | computeReducedSKGradientMatrix (FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const =0 |
| Computes second derivative of loading function with respect to stress and hardening vars. | |
| virtual void | computeTrialStressIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &strainIncrement, TimeStep *tStep) const |
| virtual void | computeReducedElasticModuli (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) const |
| FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 4, 4 > | givePlaneStrainStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 2, 2 > | give2dBeamLayerStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 5, 5 > | givePlateLayerStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 3, 3 > | giveFiberStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override |
| long | getPopulationSignature (IntArray &mask) const |
| int | testPopulation (long pop) const |
| void | clearPopulationSet () const |
| void | addNewPopulation (IntArray &mask) const |
| int | getNewPopulation (IntArray &result, IntArray &candidateMask, int degree, int size) const |
Protected Attributes | |
| LinearElasticMaterial * | linearElasticMaterial |
| Reference to bulk (undamaged) material. | |
| int | nsurf |
| Number of yield surfaces. | |
| enum oofem::MPlasticMaterial2::ReturnMappingAlgoType | rmType |
| enum oofem::MPlasticMaterial2::plastType | plType |
| bool | iterativeUpdateOfActiveConds |
| Flag indicating whether iterative update of a set of active yield conditions takes place. | |
| std ::set< long > | populationSet |
| Set for keeping record of generated populations of active yield conditions during 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 represents a base class for non-associated multisurface plasticity. The Multisurface plasticity is characterized by the following: Let \(\sigma, \varepsilon \), and \(\varepsilon^p \) be the stress, total strain, and plastic strain vectors, respectively. It is assumed that the total strain is decomposed into reversible elastic and irreversible plastic parts \(\varepsilon=\varepsilon^e+\varepsilon^p \). The elastic response is characterized in terms of elastic constitutive matrix \( D^e \) as \( \sigma=D^e(\varepsilon-\varepsilon^e) \) As long as the stress remains inside the elastic domain, the deformation process is purely elastic and the plastic strain does not change.
It is assumed that the elastic domain, denoted as \( IE \) is bounded by a composite yield surface. It is defined as
\[IE=\{(\sigma,\kappa)|f_i(\sigma,\kappa)<0,\; i\in\{1,\cdots,m\}\} \]
where \( f_i(\sigma,\kappa) \) are \( m\ge1 \) yield functions intersecting in a possibly non-smooth fashion. The vector \( \kappa \) contains internal variables controlling the evolution of yield surfaces (amount of hardening or softening). The evolution of plastic strain \( \varepsilon^p \) is expressed in Koiter's form. Assuming the non-associated plasticity, this reads
\[\label{epe} \varepsilon^p=\sum^{m}_{i=1} \lambda^i \partial_{\sigma}g_i(\sigma,\kappa) \]
where \( g_i \) are plastic potential functions. The \( \lambda^i \) are referred as plastic consistency parameters, which satisfy the following Kuhn-Tucker conditions
\[\label{ktc} \lambda^i\ge0,\;f_i\le0,\;{\rm and}\ \lambda^i f_i=0 \]
These conditions imply that in the elastic regime the yield function must remain negative and the rate of the plastic multiplier is zero (plastic strain remains constant) while in the plastic regime the yield function must be equal to zero (stress remains on the surface) and the rate of the plastic multiplier is positive. The evolution of vector of internal hardening/softening variables \( \kappa \) is expressed in terms of a general hardening/softening law of the form
\[\dot{\kappa} = \dot{\kappa}(\sigma, \lambda) \]
where \( \lambda \) is the vector of plastic consistency parameters \( \lambda_i \).
Definition at line 173 of file mplasticmaterial2.h.
|
protected |
Type that allows to distinguish between yield function and loading function.
| Enumerator | |
|---|---|
| yieldFunction | |
| loadFunction | |
Definition at line 183 of file mplasticmaterial2.h.
|
protected |
| Enumerator | |
|---|---|
| associatedPT | |
| nonassociatedPT | |
Definition at line 184 of file mplasticmaterial2.h.
|
protected |
Protected type to determine the return mapping algorithm.
| Enumerator | |
|---|---|
| mpm_ClosestPoint | |
| mpm_CuttingPlane | |
Definition at line 181 of file mplasticmaterial2.h.
| oofem::MPlasticMaterial2::MPlasticMaterial2 | ( | int | n, |
| Domain * | d ) |
Definition at line 50 of file mplasticmaterial2.C.
References associatedPT, iterativeUpdateOfActiveConds, linearElasticMaterial, mpm_ClosestPoint, plType, rmType, and oofem::StructuralMaterial::StructuralMaterial().
Referenced by oofem::DruckerPragerCutMat::DruckerPragerCutMat(), giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveRealStressVector_1d(), giveRealStressVector_3d(), giveRealStressVector_PlaneStrain(), giveRealStressVector_PlaneStress(), oofem::J2Mat::J2Mat(), and oofem::Masonry02::Masonry02().
|
virtual |
Definition at line 63 of file mplasticmaterial2.C.
References linearElasticMaterial.
|
protected |
Definition at line 1912 of file mplasticmaterial2.C.
References getPopulationSignature(), and populationSet.
Referenced by closestPointReturn().
|
protected |
Definition at line 1906 of file mplasticmaterial2.C.
References populationSet.
Referenced by closestPointReturn().
|
protected |
Definition at line 184 of file mplasticmaterial2.C.
References oofem::FloatArray::add(), oofem::FloatMatrix::add(), addNewPopulation(), associatedPT, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatMatrix::beInverseOf(), oofem::FloatArray::beProductOf(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), clearPopulationSet(), computeAlgorithmicModuli(), computeKGradientVector(), oofem::FloatArray::computeNorm(), computeReducedElasticModuli(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedSKGradientMatrix(), computeReducedStressGradientVector(), computeResidualVector(), computeStrainHardeningVarsIncrement(), computeTrialStressIncrement(), computeYieldValueAt(), getNewPopulation(), oofem::GaussPoint::giveElement(), giveMaxNumberOfActiveYieldConds(), oofem::FEMComponent::giveNumber(), oofem::GaussPoint::giveNumber(), oofem::MPlasticMaterial2Status::givePlasticStrainVector(), oofem::FloatArray::giveSize(), oofem::Material::giveStatus(), oofem::MPlasticMaterial2Status::giveStrainSpaceHardeningVars(), hasHardening(), iterativeUpdateOfActiveConds, loadFunction, nonassociatedPT, nsurf, OOFEM_ERROR, OOFEM_WARNING, PLASTIC_MATERIAL_MAX_ITERATIONS, plType, oofem::FloatArray::printYourself(), oofem::IntArray::printYourself(), RES_TOL, oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), oofem::IntArray::resize(), oofem::FloatMatrix::solveForRhs(), oofem::FloatArray::subtract(), oofem::FloatMatrix::times(), YIELD_TOL, YIELD_TOL_SECONDARY, yieldFunction, oofem::FloatArray::zero(), oofem::FloatMatrix::zero(), and oofem::IntArray::zero().
Referenced by giveRealStressVector().
|
protected |
Definition at line 1277 of file mplasticmaterial2.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatMatrix::beInverseOf(), oofem::FloatMatrix::beProductOf(), computeReducedHardeningVarsSigmaGradient(), computeReducedSKGradientMatrix(), computeReducedSSGradientMatrix(), oofem::FloatMatrix::giveNumberOfRows(), hasHardening(), nsurf, oofem::FloatMatrix::resize(), oofem::sgn(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
Referenced by closestPointReturn(), and giveConsistentStiffnessMatrix().
|
virtual |
Definition at line 178 of file mplasticmaterial2.C.
Referenced by giveRealStressVector().
|
protectedpure virtual |
Computes the derivative of yield/loading function with respect to \( \kappa \) vector
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().
|
protectedvirtual |
Reimplemented in oofem::DruckerPragerCutMat, and oofem::Masonry02.
Definition at line 1661 of file mplasticmaterial2.C.
References linearElasticMaterial.
Referenced by closestPointReturn(), computeTrialStressIncrement(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().
|
protectedpure virtual |
computes derivative of \( \kappa \) vector with respect to lambda vector
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().
|
protectedpure virtual |
Computes derivative of \( \kappa \) vector with respect to stress
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), computeAlgorithmicModuli(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().
|
protectedpure virtual |
Computes second derivative of loading function with respect to stress and hardening vars.
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), computeAlgorithmicModuli(), and giveConsistentStiffnessMatrix().
|
protectedpure virtual |
Computes second derivative of loading function with respect to stress.
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by computeAlgorithmicModuli().
|
protected |
Definition at line 1216 of file mplasticmaterial2.C.
References computeStressGradientVector(), and oofem::GaussPoint::giveMaterialMode().
Referenced by closestPointReturn(), oofem::J2Mat::computeKGradientVector(), oofem::J2Mat::computeReducedHardeningVarsLamGradient(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().
|
protected |
Definition at line 1229 of file mplasticmaterial2.C.
References oofem::FloatArray::at(), oofem::IntArray::at(), oofem::MPlasticMaterial2Status::givePlasticStrainVector(), oofem::FloatArray::giveSize(), oofem::Material::giveStatus(), nsurf, and oofem::FloatArray::resize().
Referenced by closestPointReturn().
|
protectedpure virtual |
Computes the increment of strain-space hardening variables.
| answer | Result. |
| gp | Gauss point to compute at. |
| stress | Updated stress (corresponds to newly reached state). |
| dlambda | Increment of consistency parameters. |
| dplasticStrain | Actual plastic strain increment. |
| activeConditionMap | Array of active yield conditions. |
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), and cuttingPlaneReturn().
|
protectedpure virtual |
Computes the stress gradient of yield/loading function (df/d_sigma).
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by computeReducedStressGradientVector().
|
protectedvirtual |
Computes full-space trial stress increment (elastic).
| answer | Contains result in full space. |
| gp | Integration point. |
| strainIncrement | Strain increment vector. |
| tStep | Solution step. |
Definition at line 1260 of file mplasticmaterial2.C.
References oofem::FloatArray::beProductOf(), computeReducedElasticModuli(), and oofem::GaussPoint::giveMaterialMode().
Referenced by closestPointReturn(), and cuttingPlaneReturn().
|
protectedpure virtual |
Computes the value of yield function.
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), and cuttingPlaneReturn().
|
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 88 of file mplasticmaterial2.C.
References giveSizeOfReducedHardeningVarsVector().
|
protected |
Definition at line 832 of file mplasticmaterial2.C.
References oofem::FloatArray::add(), oofem::FloatMatrix::add(), associatedPT, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), computeKGradientVector(), computeReducedElasticModuli(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedStressGradientVector(), computeStrainHardeningVarsIncrement(), computeTrialStressIncrement(), computeYieldValueAt(), oofem::GaussPoint::giveElement(), giveMaxNumberOfActiveYieldConds(), oofem::FEMComponent::giveNumber(), oofem::GaussPoint::giveNumber(), oofem::FloatMatrix::giveNumberOfRows(), oofem::MPlasticMaterial2Status::givePlasticStrainVector(), oofem::Material::giveStatus(), oofem::MPlasticMaterial2Status::giveStrainSpaceHardeningVars(), hasHardening(), iterativeUpdateOfActiveConds, loadFunction, nonassociatedPT, nsurf, OOFEM_WARNING, PLASTIC_MATERIAL_MAX_ITERATIONS, plType, oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), oofem::IntArray::resize(), oofem::FloatMatrix::solveForRhs(), oofem::FloatArray::subtract(), YIELD_TOL, yieldFunction, oofem::FloatArray::zero(), oofem::FloatMatrix::zero(), and oofem::IntArray::zero().
Referenced by giveRealStressVector().
|
protected |
Definition at line 1920 of file mplasticmaterial2.C.
References oofem::IntArray::at(), getPopulationSignature(), oofem::IntArray::giveSize(), populationSet, oofem::IntArray::resize(), testPopulation(), and oofem::IntArray::zero().
Referenced by closestPointReturn().
|
protected |
Definition at line 1878 of file mplasticmaterial2.C.
References oofem::IntArray::at(), and oofem::IntArray::giveSize().
Referenced by addNewPopulation(), and getNewPopulation().
|
overrideprotectedvirtual |
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.
Definition at line 1760 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
overrideprotectedvirtual |
Method for computing 2d beam layer stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d beam layer 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 1782 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
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 1672 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 195 of file mplasticmaterial2.h.
|
protectedvirtual |
Definition at line 1339 of file mplasticmaterial2.C.
References oofem::FloatArray::add(), oofem::FloatMatrix::add(), associatedPT, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatMatrix::beInverseOf(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), computeAlgorithmicModuli(), computeKGradientVector(), computeReducedElasticModuli(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedSKGradientMatrix(), computeReducedStressGradientVector(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatMatrix::giveNumberOfRows(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::MPlasticMaterial2Status::giveTempActiveConditionMap(), oofem::MPlasticMaterial2Status::giveTempGamma(), oofem::MPlasticMaterial2Status::giveTempStateFlag(), oofem::MPlasticMaterial2Status::giveTempStrainSpaceHardeningVarsVector(), oofem::StructuralMaterialStatus::giveTempStressVector(), hasHardening(), loadFunction, MPlasticMaterial2(), oofem::FloatMatrix::negated(), nonassociatedPT, nsurf, plType, oofem::FloatMatrix::resize(), oofem::FloatMatrix::times(), yieldFunction, and oofem::FloatMatrix::zero().
Referenced by give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), oofem::Masonry02::give2dInterfaceMaterialStiffnessMatrix(), give3dMaterialStiffnessMatrix(), giveFiberStiffMtrx(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), and givePlateLayerStiffMtrx().
|
protectedvirtual |
Definition at line 1497 of file mplasticmaterial2.C.
References oofem::FloatArray::add(), oofem::FloatMatrix::add(), associatedPT, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatMatrix::beInverseOf(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), computeKGradientVector(), computeReducedElasticModuli(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedStressGradientVector(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatMatrix::giveNumberOfRows(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::MPlasticMaterial2Status::giveStrainSpaceHardeningVars(), oofem::StructuralMaterialStatus::giveStressVector(), oofem::MPlasticMaterial2Status::giveTempActiveConditionMap(), oofem::MPlasticMaterial2Status::giveTempGamma(), oofem::MPlasticMaterial2Status::giveTempStateFlag(), hasHardening(), loadFunction, MPlasticMaterial2(), oofem::FloatMatrix::negated(), nonassociatedPT, nsurf, plType, oofem::FloatMatrix::resize(), and yieldFunction.
Referenced by give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), oofem::Masonry02::give2dInterfaceMaterialStiffnessMatrix(), give3dMaterialStiffnessMatrix(), giveFiberStiffMtrx(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), and givePlateLayerStiffMtrx().
|
overrideprotectedvirtual |
Method for computing 1d fiber stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d fiber 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 1830 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
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 1854 of file mplasticmaterial2.C.
References oofem::StructuralMaterial::computePrincipalValues(), oofem::GaussPoint::giveMaterialMode(), oofem::MPlasticMaterial2Status::givePlasticStrainVector(), oofem::Material::giveStatus(), and oofem::principal_strain.
|
inline |
Returns reference to undamaged (bulk) material.
Definition at line 198 of file mplasticmaterial2.h.
References linearElasticMaterial.
|
protectedpure virtual |
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), and cuttingPlaneReturn().
|
overrideprotectedvirtual |
Method for computing plane strain stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane strain 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. Note: as already described, if zero strain component is imposed (Plane strain, ..) this condition must be taken into account in geometrical relations, and corresponding component has to be included in reduced vector. (So plane strain conditions are \( \epsilon_z = \gamma_{xz} = \gamma_{yz} = 0 \), but relations for \( \epsilon_z\) and \(\sigma_z\) are included).
| 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 1737 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
overrideprotectedvirtual |
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 1712 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
overrideprotectedvirtual |
Method for computing 2d plate layer stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d plate layer 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 1807 of file mplasticmaterial2.C.
References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), linearElasticMaterial, mpm_ClosestPoint, and rmType.
|
overridevirtual |
Computes the real stress vector for given total strain and integration point. The total strain is defined as strain computed directly from displacement field at given time. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process.
| answer | Stress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress. |
| gp | Integration point. |
| reducedStrain | Strain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 95 of file mplasticmaterial2.C.
References oofem::FloatArray::at(), closestPointReturn(), computeDamage(), cuttingPlaneReturn(), oofem::MPlasticMaterial2Status::giveStateFlag(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::Material::initTempStatus(), oofem::MPlasticMaterial2Status::letTempDamageBe(), oofem::MPlasticMaterial2Status::letTempPlasticStrainVectorBe(), oofem::MPlasticMaterial2Status::letTempStateFlagBe(), oofem::MPlasticMaterial2Status::letTempStrainSpaceHardeningVarsVectorBe(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), mpm_ClosestPoint, nsurf, rmType, oofem::MPlasticMaterial2Status::setTempActiveConditionMap(), oofem::MPlasticMaterial2Status::setTempGamma(), and oofem::FloatArray::times().
Referenced by giveRealStressVector_1d(), giveRealStressVector_3d(), giveRealStressVector_PlaneStrain(), and giveRealStressVector_PlaneStress().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 231 of file mplasticmaterial2.h.
References giveRealStressVector(), and MPlasticMaterial2().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 213 of file mplasticmaterial2.h.
References giveRealStressVector(), and MPlasticMaterial2().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_3d.
Reimplemented from oofem::StructuralMaterial.
Definition at line 219 of file mplasticmaterial2.h.
References giveRealStressVector(), and MPlasticMaterial2().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 225 of file mplasticmaterial2.h.
References giveRealStressVector(), and MPlasticMaterial2().
|
inlinevirtual |
Reimplemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Definition at line 243 of file mplasticmaterial2.h.
|
inlinevirtual |
Reimplemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Definition at line 244 of file mplasticmaterial2.h.
Referenced by CreateStatus().
|
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 202 of file mplasticmaterial2.h.
References oofem::Material::give(), and tAlpha.
|
protectedpure virtual |
Indicates, whether receiver model has hardening/softening behavior or behaves according to perfect plasticity theory.
Implemented in oofem::DruckerPragerCutMat, oofem::J2Mat, and oofem::Masonry02.
Referenced by closestPointReturn(), computeAlgorithmicModuli(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().
|
overridevirtual |
Tests if material supports material mode.
| mode | Required material mode. |
Reimplemented from oofem::Material.
Definition at line 73 of file mplasticmaterial2.C.
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 200 of file mplasticmaterial2.h.
|
protected |
Definition at line 1896 of file mplasticmaterial2.C.
References populationSet.
Referenced by getNewPopulation().
|
protected |
Flag indicating whether iterative update of a set of active yield conditions takes place.
Definition at line 186 of file mplasticmaterial2.h.
Referenced by closestPointReturn(), cuttingPlaneReturn(), and MPlasticMaterial2().
|
protected |
Reference to bulk (undamaged) material.
Definition at line 177 of file mplasticmaterial2.h.
Referenced by oofem::DruckerPragerCutMat::computeReducedElasticModuli(), oofem::Masonry02::computeReducedElasticModuli(), computeReducedElasticModuli(), oofem::DruckerPragerCutMat::DruckerPragerCutMat(), give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), give3dMaterialStiffnessMatrix(), giveFiberStiffMtrx(), oofem::DruckerPragerCutMat::giveLinearElasticMaterial(), giveLinearElasticMaterial(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), givePlateLayerStiffMtrx(), oofem::DruckerPragerCutMat::initializeFrom(), oofem::J2Mat::initializeFrom(), oofem::Masonry02::initializeFrom(), oofem::J2Mat::J2Mat(), oofem::Masonry02::Masonry02(), MPlasticMaterial2(), and ~MPlasticMaterial2().
|
protected |
Number of yield surfaces.
Definition at line 179 of file mplasticmaterial2.h.
Referenced by closestPointReturn(), computeAlgorithmicModuli(), computeResidualVector(), cuttingPlaneReturn(), oofem::DruckerPragerCutMat::DruckerPragerCutMat(), giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveRealStressVector(), oofem::J2Mat::J2Mat(), and oofem::Masonry02::Masonry02().
|
protected |
|
mutableprotected |
Set for keeping record of generated populations of active yield conditions during return.
Definition at line 188 of file mplasticmaterial2.h.
Referenced by addNewPopulation(), clearPopulationSet(), getNewPopulation(), and testPopulation().
|
protected |
Referenced by oofem::DruckerPragerCutMat::DruckerPragerCutMat(), give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), oofem::Masonry02::give2dInterfaceMaterialStiffnessMatrix(), give3dMaterialStiffnessMatrix(), giveFiberStiffMtrx(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), givePlateLayerStiffMtrx(), giveRealStressVector(), oofem::J2Mat::initializeFrom(), oofem::Masonry02::initializeFrom(), oofem::Masonry02::Masonry02(), and MPlasticMaterial2().