|
OOFEM 3.0
|
#include <druckerpragercutmat.h>
Public Member Functions | |
| DruckerPragerCutMat (int n, Domain *d) | |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| bool | hasCastingTimeSupport () const override |
| void | initializeFrom (InputRecord &ir) override |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| LinearElasticMaterial * | giveLinearElasticMaterial () |
| Returns a reference to the basic elastic material. | |
| int | giveSizeOfFullHardeningVarsVector () const override |
| int | giveSizeOfReducedHardeningVarsVector (GaussPoint *) const override |
| Public Member Functions inherited from oofem::MPlasticMaterial2 | |
| 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 |
| 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 |
| 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 | |
| int | giveMaxNumberOfActiveYieldConds (GaussPoint *gp) const override |
| double | computeYieldValueAt (GaussPoint *gp, int isurf, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const override |
| Computes the value of yield function. | |
| void | computeStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &stressSpaceHardeningVars) const override |
| Computes the stress gradient of yield/loading function (df/d_sigma). | |
| void | computeReducedSSGradientMatrix (FloatMatrix &gradientMatrix, int isurf, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const override |
| Computes second derivative of yield/loading function with respect to stress. | |
| void | computeReducedElasticModuli (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) const override |
| int | hasHardening () const override |
| Functions related to hardening. | |
| void | computeStrainHardeningVarsIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &stress, const FloatArray &dlambda, const FloatArray &dplasticStrain, const IntArray &activeConditionMap) const override |
| Compute dot(kappa_1), dot(kappa_2) etc. | |
| void | computeKGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const override |
| Computes the derivative of yield/loading function with respect to kappa_1, kappa_2 etc. | |
| void | computeReducedSKGradientMatrix (FloatMatrix &gradientMatrix, int isurf, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) const override |
| computes mixed derivative of load function with respect to stress and hardening variables | |
| void | computeReducedHardeningVarsSigmaGradient (FloatMatrix &answer, GaussPoint *gp, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &dlambda) const override |
| computes dk(i)/dsig(j) gradient matrix | |
| void | computeReducedHardeningVarsLamGradient (FloatMatrix &answer, GaussPoint *gp, int actSurf, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &dlambda) const override |
| computes dKappa_i/dLambda_j | |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| Protected Member Functions inherited from oofem::MPlasticMaterial2 | |
| 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 |
| void | computeReducedStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) const |
| virtual void | computeTrialStressIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &strainIncrement, 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 | |
| double | G = 0. |
| Elastic shear modulus. | |
| double | K = 0. |
| Elastic bulk modulus. | |
| double | H = 0. |
| Hardening modulus. | |
| double | sigT = 0. |
| Uniaxial tensile strength for cut-off. | |
| double | tau0 = 0. |
| Initial yield stress under pure shear. | |
| double | alpha = 0. |
| Friction coefficient. | |
| double | alphaPsi = 0. |
| Dilatancy coefficient (allowing non-associated plasticity). | |
| double | yieldTol = 0. |
| Tolerance of the error in the yield criterion. | |
| int | newtonIter = 30 |
| Maximum number of iterations in lambda search. | |
| double | omegaCrit = 0. |
| Maximum damage value. | |
| double | a = 0. |
| Parameter for damage computation from cumulative plastic strain. | |
| Protected Attributes inherited from oofem::MPlasticMaterial2 | |
| 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. | |
| Protected Types inherited from oofem::MPlasticMaterial2 | |
| 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 } |
This class implements an isotropic elasto-plasto-damage material with Drucker-Prager yield condition, tension cut-off, non-associated flow rule, linear isotropic hardening and isotropic damage.
Definition at line 63 of file druckerpragercutmat.h.
| oofem::DruckerPragerCutMat::DruckerPragerCutMat | ( | int | n, |
| Domain * | d ) |
Definition at line 49 of file druckerpragercutmat.C.
References oofem::MPlasticMaterial2::linearElasticMaterial, oofem::MPlasticMaterial2::MPlasticMaterial2(), oofem::MPlasticMaterial2::mpm_CuttingPlane, oofem::MPlasticMaterial2::nonassociatedPT, oofem::MPlasticMaterial2::nsurf, oofem::MPlasticMaterial2::plType, and oofem::MPlasticMaterial2::rmType.
|
overrideprotectedvirtual |
Computes the derivative of yield/loading function with respect to kappa_1, kappa_2 etc.
Implements oofem::MPlasticMaterial2.
Definition at line 231 of file druckerpragercutmat.C.
References oofem::FloatArray::at(), H, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overrideprotectedvirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 212 of file druckerpragercutmat.C.
References oofem::MPlasticMaterial2::linearElasticMaterial.
|
overrideprotectedvirtual |
computes dKappa_i/dLambda_j
Implements oofem::MPlasticMaterial2.
Definition at line 259 of file druckerpragercutmat.C.
References alphaPsi, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
computes dk(i)/dsig(j) gradient matrix
Implements oofem::MPlasticMaterial2.
Definition at line 251 of file druckerpragercutmat.C.
References oofem::GaussPoint::giveMaterialMode(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
computes mixed derivative of load function with respect to stress and hardening variables
Implements oofem::MPlasticMaterial2.
Definition at line 243 of file druckerpragercutmat.C.
References oofem::GaussPoint::giveMaterialMode(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
Computes second derivative of yield/loading function with respect to stress.
Implements oofem::MPlasticMaterial2.
Definition at line 149 of file druckerpragercutmat.C.
References oofem::FloatMatrix::at(), oofem::StructuralMaterial::computeDeviator(), oofem::StructuralMaterial::computeSecondStressInvariant(), oofem::GaussPoint::giveMaterialMode(), OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), and oofem::FloatMatrix::zero().
|
overrideprotectedvirtual |
Compute dot(kappa_1), dot(kappa_2) etc.
Implements oofem::MPlasticMaterial2.
Definition at line 220 of file druckerpragercutmat.C.
References alphaPsi, oofem::FloatArray::at(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overrideprotectedvirtual |
Computes the stress gradient of yield/loading function (df/d_sigma).
Implements oofem::MPlasticMaterial2.
Definition at line 115 of file druckerpragercutmat.C.
References alphaPsi, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::StructuralMaterial::computeDeviator(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondStressInvariant(), oofem::principal_stress, and oofem::FloatArray::resize().
|
overrideprotectedvirtual |
Computes the value of yield function.
Implements oofem::MPlasticMaterial2.
Definition at line 95 of file druckerpragercutmat.C.
References alpha, oofem::FloatArray::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computePrincipalValues(), oofem::StructuralMaterial::computeSecondStressInvariant(), H, oofem::principal_stress, sigT, and tau0.
|
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 89 of file druckerpragercutmat.C.
References giveSizeOfReducedHardeningVarsVector().
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 115 of file druckerpragercutmat.h.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 116 of file druckerpragercutmat.h.
References _IFT_DruckerPragerCutMat_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 274 of file druckerpragercutmat.C.
References oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
inline |
Returns a reference to the basic elastic material.
Definition at line 119 of file druckerpragercutmat.h.
References oofem::MPlasticMaterial2::linearElasticMaterial.
|
inlineoverrideprotectedvirtual |
Implements oofem::MPlasticMaterial2.
Definition at line 125 of file druckerpragercutmat.h.
|
inlineoverridevirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 121 of file druckerpragercutmat.h.
|
inlineoverridevirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 122 of file druckerpragercutmat.h.
Referenced by CreateStatus().
|
inlineoverridevirtual |
Tests if material supports casting time
Reimplemented from oofem::Material.
Definition at line 107 of file druckerpragercutmat.h.
|
inlineoverrideprotectedvirtual |
Functions related to hardening.
Implements oofem::MPlasticMaterial2.
Definition at line 137 of file druckerpragercutmat.h.
|
overridevirtual |
Tests if material supports material mode.
| mode | Required material mode. |
Reimplemented from oofem::Material.
Definition at line 60 of file druckerpragercutmat.C.
|
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 67 of file druckerpragercutmat.C.
References _IFT_DruckerPragerCutMat_a, _IFT_DruckerPragerCutMat_alpha, _IFT_DruckerPragerCutMat_alphapsi, _IFT_DruckerPragerCutMat_h, _IFT_DruckerPragerCutMat_newtonIter, _IFT_DruckerPragerCutMat_omegaCrit, _IFT_DruckerPragerCutMat_sigT, _IFT_DruckerPragerCutMat_tau0, _IFT_DruckerPragerCutMat_yieldTol, a, alpha, alphaPsi, G, H, IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, K, oofem::MPlasticMaterial2::linearElasticMaterial, newtonIter, omegaCrit, sigT, tau0, and yieldTol.
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 113 of file druckerpragercutmat.h.
|
protected |
Parameter for damage computation from cumulative plastic strain.
Definition at line 100 of file druckerpragercutmat.h.
Referenced by initializeFrom().
|
protected |
Friction coefficient.
Definition at line 85 of file druckerpragercutmat.h.
Referenced by computeYieldValueAt(), and initializeFrom().
|
protected |
Dilatancy coefficient (allowing non-associated plasticity).
Definition at line 88 of file druckerpragercutmat.h.
Referenced by computeReducedHardeningVarsLamGradient(), computeStrainHardeningVarsIncrement(), computeStressGradientVector(), and initializeFrom().
|
protected |
Elastic shear modulus.
Definition at line 70 of file druckerpragercutmat.h.
Referenced by initializeFrom().
|
protected |
Hardening modulus.
Definition at line 76 of file druckerpragercutmat.h.
Referenced by computeKGradientVector(), computeYieldValueAt(), and initializeFrom().
|
protected |
Elastic bulk modulus.
Definition at line 73 of file druckerpragercutmat.h.
Referenced by initializeFrom().
|
protected |
Maximum number of iterations in lambda search.
Definition at line 94 of file druckerpragercutmat.h.
Referenced by initializeFrom().
|
protected |
Maximum damage value.
Definition at line 97 of file druckerpragercutmat.h.
Referenced by initializeFrom().
|
protected |
Uniaxial tensile strength for cut-off.
Definition at line 79 of file druckerpragercutmat.h.
Referenced by computeYieldValueAt(), and initializeFrom().
|
protected |
Initial yield stress under pure shear.
Definition at line 82 of file druckerpragercutmat.h.
Referenced by computeYieldValueAt(), and initializeFrom().
|
protected |
Tolerance of the error in the yield criterion.
Definition at line 91 of file druckerpragercutmat.h.
Referenced by initializeFrom().