|
OOFEM 3.0
|
#include <latticeplasticitydamage.h>
Public Member Functions | |
| LatticePlasticityDamage (int n, Domain *d) | |
| Constructor. | |
| const char * | giveInputRecordName () const override |
| const char * | giveClassName () const override |
| void | initializeFrom (InputRecord &ir) override |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| double | give (int aProperty, GaussPoint *gp) const override |
| FloatMatrixF< 6, 6 > | give3dLatticeStiffnessMatrix (MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| FloatArrayF< 3 > | computeFVector (const FloatArrayF< 3 > &sigma, const double deltaLambda, GaussPoint *gp, TimeStep *tStep) const |
| FloatArrayF< 3 > | computeMVector (const FloatArrayF< 3 > &sigma, const double deltaLambda, GaussPoint *gp, TimeStep *tStep) const |
| FloatMatrixF< 3, 3 > | computeDMMatrix (const FloatArrayF< 3 > &sigma, const double deltaLambda, GaussPoint *gp, TimeStep *tStep) const |
| FloatMatrixF< 4, 4 > | computeJacobian (const FloatArrayF< 3 > &sigma, const double tempKappa, const double deltaLambda, GaussPoint *gp, TimeStep *tStep) const |
| virtual double | computeDamageParam (double kappaOne, double kappaTwo, GaussPoint *gp, TimeStep *tStep) const |
| FloatArrayF< 6 > | giveLatticeStress3d (const FloatArrayF< 6 > &jump, GaussPoint *gp, TimeStep *tStep) override |
| FloatArrayF< 6 > | performPlasticityReturn (GaussPoint *gp, const FloatArrayF< 6 > &reducedStrain, TimeStep *tStep) const |
| void | performDamageEvaluation (GaussPoint *gp, FloatArrayF< 6 > &reducedStrain, TimeStep *tStep) const |
| double | performRegularReturn (FloatArrayF< 3 > &stress, LatticePlasticityDamage_ReturnResult &returnResult, double yieldValue, GaussPoint *gp, TimeStep *tStep) const |
| double | computeYieldValue (const FloatArrayF< 3 > &sigma, const double tempKappa, GaussPoint *gp, TimeStep *tStep) const |
| double | computeHardening (const double kappa, GaussPoint *gp) const |
| double | computeDHardeningDKappa (const double kappa, GaussPoint *gp) const |
| double | computeDDHardeningDDKappa (const double kappa, GaussPoint *gp) const |
| double | computeDuctilityMeasure (FloatArray &stress, double ductilityParameter) const |
| double | computeYieldStress (double kappaP, GaussPoint *gp) |
| const double | computeEquivalentStress (const FloatArray &tempSigma) const |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
| virtual FloatArrayF< 6 > | giveReducedStrain (GaussPoint *gp, TimeStep *tStep) const |
| Public Member Functions inherited from oofem::LatticeLinearElastic | |
| LatticeLinearElastic (int n, Domain *d) | |
| LatticeLinearElastic (int n, Domain *d, double eNormalMean, double alphaOne, double alphaTwo) | |
| FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 3, 3 > | give2dLatticeStiffnessMatrix (MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override |
| Interface * | giveInterface (InterfaceType) override |
| virtual void | giveRandomParameters (FloatArray ¶m) |
| MaterialStatus * | giveStatus (GaussPoint *gp) const override |
| Public Member Functions inherited from oofem::LatticeStructuralMaterial | |
| LatticeStructuralMaterial (int n, Domain *d) | |
| void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override |
| virtual bool | hasAnalyticalTangentStiffness () const |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| const char * | giveClassName () const override |
| virtual double | giveLatticeStress1d (const double strain, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatArrayF< 3 > | giveLatticeStress2d (const FloatArrayF< 3 > &strain, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 1, 1 > | give1dLatticeStiffnessMatrix (MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *atTime) 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 | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) const |
| virtual FloatArrayF< 6 > | giveRealStressVector_3d (const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| virtual FloatArrayF< 4 > | giveRealStressVector_PlaneStrain (const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_3d. | |
| virtual FloatArray | giveRealStressVector_StressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| virtual FloatArray | giveRealStressVector_ShellStressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 3 > | giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 1 > | giveRealStressVector_1d (const FloatArrayF< 1 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 2 > | giveRealStressVector_Warping (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 2 > | giveRealStressVector_2dBeamLayer (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 5 > | giveRealStressVector_PlateLayer (const FloatArrayF< 5 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 3 > | giveRealStressVector_Fiber (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 3 > | giveRealStressVector_2dPlateSubSoil (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation is not provided. | |
| virtual FloatArrayF< 6 > | giveRealStressVector_3dBeamSubSoil (const FloatArrayF< 6 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 9 > | giveFirstPKStressVector_3d (const FloatArrayF< 9 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| virtual FloatArrayF< 5 > | giveFirstPKStressVector_PlaneStrain (const FloatArrayF< 5 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_3d. | |
| virtual FloatArray | giveFirstPKStressVector_StressControl (const FloatArray &reducedvF, const IntArray &FControl, GaussPoint *gp, TimeStep *tStep) const |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| virtual FloatArrayF< 4 > | giveFirstPKStressVector_PlaneStress (const FloatArrayF< 4 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_StressControl. | |
| virtual FloatArrayF< 1 > | giveFirstPKStressVector_1d (const FloatArrayF< 1 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_StressControl. | |
| virtual void | giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual void | giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| double | giveReferenceTemperature () |
| virtual FloatArray | computeStressIndependentStrainVector (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const |
| FloatArrayF< 6 > | computeStressIndependentStrainVector_3d (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const |
| virtual FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | giveStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 9, 9 > | give3dMaterialStiffnessMatrix_dPdF (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| void | giveStressDependentPartOfStrainVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) const |
| int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| virtual FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 4, 4 > | givePlaneStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 4, 4 > | givePlaneStrainStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | givePlaneStrainStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 1, 1 > | give1dStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 2, 2 > | give2dBeamLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | givePlateLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 3, 3 > | giveFiberStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 3, 3 > | give2dPlateSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 6, 6 > | give3dBeamSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| Public Member Functions inherited from oofem::Material | |
| Material (int n, Domain *d) | |
| virtual | ~Material ()=default |
| Destructor. | |
| virtual double | giveCharacteristicValue (MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const |
| Returns characteristic value of the receiver. | |
| virtual 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 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) |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
| Public Member Functions inherited from oofem::RandomMaterialExtensionInterface | |
| RandomMaterialExtensionInterface () | |
| Constructor. | |
| virtual | ~RandomMaterialExtensionInterface () |
| Destructor. | |
| void | initializeFrom (InputRecord &ir) |
| void | giveInputRecord (DynamicInputRecord &ir) |
| bool | give (int key, GaussPoint *gp, double &value) const |
| Public Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
Protected Types | |
| enum | LatticePlasticityDamage_ReturnResult { RR_Unknown , RR_NotConverged , RR_Converged } |
Protected Member Functions | |
| virtual double | giveTensileStrength (GaussPoint *gp, TimeStep *tStep) const |
| virtual double | giveCompressiveStrength (GaussPoint *gp, TimeStep *tStep) const |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *atTime) override |
| int | giveIPValueSize (InternalStateType type, GaussPoint *gp) |
| int | giveIntVarCompFullIndx (IntArray &answer, InternalStateType type, MaterialMode mmode) |
| InternalStateValueType | giveIPValueType (InternalStateType type) |
| Protected Member Functions inherited from oofem::RandomMaterialExtensionInterface | |
| void | _generateStatusVariables (GaussPoint *) const |
Protected Attributes | |
| double | initialYieldStress = 0. |
| double | ft = 0. |
| tensile strength | |
| double | fc = 0. |
| compressive strength | |
| double | frictionAngleOne = 0. |
| frictional angle of the yield surface | |
| double | frictionAngleTwo = 0. |
| frictional angle of the yield surface | |
| double | flowAngleOne = 0. |
| frictional angle of the plastic potential | |
| double | flowAngleTwo = 0. |
| frictional angle of the plastic potential | |
| double | wf = 0. |
| determines the softening -> corresponds to crack opening (not strain) when tension stress vanishes | |
| int | softeningType = 0 |
| softening type determines the type of softening. 0 is exponential and 1 is bilinear. | |
| double | ftOneRatio = 0. |
| ratio of tensile stress value for bilinear stress-crack opening curve | |
| double | wfOne = 0. |
| crack opening value for bilinear stress-crack opening curve | |
| double | aHard = 0. |
| hardening parameter | |
| double | yieldTol = 0. |
| yield tolerance | |
| int | newtonIter = 0 |
| maximum number of iterations for stress return | |
| int | numberOfSubIncrements = 0 |
| int | damageFlag = 0 |
| damageFlag | |
| Protected Attributes inherited from oofem::LatticeLinearElastic | |
| double | eNormalMean = 0. |
| Normal modulus. | |
| double | alphaOne = 0. |
| Ratio of shear and normal modulus. | |
| double | alphaTwo = 0. |
| Ratio of torsion and normal modulus. | |
| double | coefficientOfVariation = 0. |
| coefficient variation of the Gaussian distribution | |
| double | localRandomType = 0. |
| flag which chooses between no distribution (0) and Gaussian distribution (1) | |
| double | cAlpha = 0. |
| parameter which allows to prescribed thermal displacement | |
| Protected Attributes inherited from oofem::LatticeStructuralMaterial | |
| double | referenceTemperature = 0. |
| Reference temperature (temperature, when material has been built into structure). | |
| 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. | |
| Protected Attributes inherited from oofem::RandomMaterialExtensionInterface | |
| IntArray | randVariables |
| Array of randomized variables (identified by a key). | |
| IntArray | randomVariableGenerators |
| Array of generators id's for corresponding randomized variables. | |
Additional Inherited Members | |
| Static Public Member Functions inherited from oofem::StructuralMaterial | |
| static int | giveSymVI (int ind1, int ind2) |
| static int | giveVI (int ind1, int ind2) |
| static FloatMatrixF< 9, 9 > | convert_dSdE_2_dPdF_3D (const FloatMatrixF< 6, 6 > &dSdE, const FloatArrayF< 6 > &S, const FloatArrayF< 9 > &F) |
| static FloatMatrixF< 5, 5 > | convert_dSdE_2_dPdF_PlaneStrain (const FloatMatrixF< 4, 4 > &dSdE, const FloatArrayF< 4 > &S, const FloatArrayF< 5 > &F) |
| static FloatMatrixF< 4, 4 > | convert_dSdE_2_dPdF_PlaneStress (const FloatMatrixF< 3, 3 > &dSdE, const FloatArrayF< 3 > &S, const FloatArrayF< 4 > &F) |
| static FloatMatrixF< 1, 1 > | convert_dSdE_2_dPdF_1D (const FloatMatrixF< 1, 1 > &dSdE, const FloatArrayF< 1 > &S, const FloatArrayF< 1 > &F) |
| static void | computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode) |
| Common functions for convenience. | |
| static FloatArrayF< 3 > | computePrincipalValues (const FloatMatrixF< 3, 3 > &s) |
| static FloatArrayF< 3 > | computePrincipalValues (double I1, double I2, double I3) |
| static void | computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode) |
| static std::pair< FloatArrayF< 3 >, FloatMatrixF< 3, 3 > > | computePrincipalValDir (const FloatMatrixF< 3, 3 > &s) |
| static FloatArrayF< 6 > | computeDeviator (const FloatArrayF< 6 > &s) |
| static std::pair< FloatArrayF< 6 >, double > | computeDeviatoricVolumetricSplit (const FloatArrayF< 6 > &s) |
| static FloatArrayF< 6 > | computeDeviatoricVolumetricSum (const FloatArrayF< 6 > &dev, double mean) |
| static FloatArrayF< 6 > | applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double GModulus) |
| static FloatArrayF< 6 > | applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double GModulus) |
| static FloatArrayF< 6 > | applyElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu) |
| static double | computeStressNorm (const FloatArrayF< 6 > &stress) |
| static double | computeFirstInvariant (const FloatArrayF< 6 > &s) |
| static double | computeSecondStressInvariant (const FloatArrayF< 6 > &s) |
| static double | computeThirdStressInvariant (const FloatArrayF< 6 > &s) |
| static double | computeFirstCoordinate (const FloatArrayF< 6 > &s) |
| static double | computeSecondCoordinate (const FloatArrayF< 6 > &s) |
| static double | computeThirdCoordinate (const FloatArrayF< 6 > &s) |
| static int | giveVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| static int | giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode) |
| static void | giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| static int | giveSizeOfVoigtVector (MaterialMode mmode) |
| static int | giveSizeOfVoigtSymVector (MaterialMode mmode) |
| static void | giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. | |
| static void | giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced deformation gradient Voigt vector (2nd order tensor). | |
| static void | giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. | |
| static void | giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. | |
| static void | giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. | |
| static void | giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. | |
| static void | giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. | |
| static void | giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. | |
| static FloatArrayF< 6 > | transformStrainVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &strain, bool transpose=false) |
| static FloatArrayF< 6 > | transformStressVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &stress, bool transpose=false) |
| static double | computeVonMisesStress (const FloatArray ¤tStress) |
| static double | computeVonMisesStress_3D (const FloatArrayF< 6 > &stress) |
| static double | computeVonMisesStress_PlaneStress (const FloatArrayF< 3 > &stress) |
| static FloatMatrixF< 6, 6 > | giveStrainVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false) |
| static FloatMatrixF< 3, 3 > | give2DStrainVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false) |
| static FloatMatrixF< 6, 6 > | giveStressVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false) |
| static FloatMatrixF< 3, 3 > | givePlaneStressVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false) |
| static void | sortPrincDirAndValCloseTo (FloatArray &pVal, FloatMatrix &pDir, const FloatMatrix &toPDir) |
| Static Public Attributes inherited from oofem::StructuralMaterial | |
| static std::array< std::array< int, 3 >, 3 > | vIindex |
| Voigt index map. | |
| static std::array< std::array< int, 3 >, 3 > | svIndex |
| Symmetric Voigt index map. | |
This class implements a local random plasticity damage model for concrete for lattice elements.
Definition at line 134 of file latticeplasticitydamage.h.
| Enumerator | |
|---|---|
| RR_Unknown | |
| RR_NotConverged | |
| RR_Converged | |
Definition at line 138 of file latticeplasticitydamage.h.
| oofem::LatticePlasticityDamage::LatticePlasticityDamage | ( | int | n, |
| Domain * | d ) |
Constructor.
Definition at line 54 of file latticeplasticitydamage.C.
References oofem::LatticeLinearElastic::LatticeLinearElastic().
Referenced by oofem::LatticePlasticityDamageViscoelastic::LatticePlasticityDamageViscoelastic().
|
virtual |
Definition at line 126 of file latticeplasticitydamage.C.
References oofem::LatticeLinearElastic::eNormalMean, frictionAngleOne, frictionAngleTwo, ftOneRatio, giveCompressiveStrength(), oofem::GaussPoint::giveElement(), oofem::LatticeStructuralElement::giveLength(), giveTensileStrength(), OOFEM_ERROR, softeningType, wf, and wfOne.
Referenced by performDamageEvaluation().
| double oofem::LatticePlasticityDamage::computeDDHardeningDDKappa | ( | const double | kappa, |
| GaussPoint * | gp ) const |
Definition at line 235 of file latticeplasticitydamage.C.
References aHard.
| double oofem::LatticePlasticityDamage::computeDHardeningDKappa | ( | const double | kappa, |
| GaussPoint * | gp ) const |
Definition at line 229 of file latticeplasticitydamage.C.
References aHard.
Referenced by computeDMMatrix(), and computeFVector().
| FloatMatrixF< 3, 3 > oofem::LatticePlasticityDamage::computeDMMatrix | ( | const FloatArrayF< 3 > & | sigma, |
| const double | deltaLambda, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 301 of file latticeplasticitydamage.C.
References oofem::FloatArrayF< N >::at(), oofem::FloatMatrixF< N, M >::at(), computeDHardeningDKappa(), computeHardening(), flowAngleOne, flowAngleTwo, ft, giveCompressiveStrength(), and giveTensileStrength().
Referenced by computeJacobian().
| double oofem::LatticePlasticityDamage::computeDuctilityMeasure | ( | FloatArray & | stress, |
| double | ductilityParameter ) const |
| const double oofem::LatticePlasticityDamage::computeEquivalentStress | ( | const FloatArray & | tempSigma | ) | const |
| FloatArrayF< 3 > oofem::LatticePlasticityDamage::computeFVector | ( | const FloatArrayF< 3 > & | sigma, |
| const double | deltaLambda, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 241 of file latticeplasticitydamage.C.
References oofem::FloatArrayF< N >::at(), computeDHardeningDKappa(), computeHardening(), frictionAngleOne, frictionAngleTwo, giveCompressiveStrength(), giveTensileStrength(), and oofem::norm().
Referenced by computeJacobian().
| double oofem::LatticePlasticityDamage::computeHardening | ( | const double | kappa, |
| GaussPoint * | gp ) const |
Definition at line 223 of file latticeplasticitydamage.C.
References aHard.
Referenced by computeDMMatrix(), computeFVector(), computeMVector(), computeYieldValue(), and performDamageEvaluation().
| FloatMatrixF< 4, 4 > oofem::LatticePlasticityDamage::computeJacobian | ( | const FloatArrayF< 3 > & | sigma, |
| const double | tempKappa, | ||
| const double | deltaLambda, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 570 of file latticeplasticitydamage.C.
References oofem::LatticeLinearElastic::alphaOne, oofem::FloatMatrixF< N, M >::at(), computeDMMatrix(), computeFVector(), computeMVector(), and oofem::LatticeLinearElastic::eNormalMean.
Referenced by performRegularReturn().
| FloatArrayF< 3 > oofem::LatticePlasticityDamage::computeMVector | ( | const FloatArrayF< 3 > & | sigma, |
| const double | deltaLambda, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 273 of file latticeplasticitydamage.C.
References oofem::FloatArrayF< N >::at(), computeHardening(), flowAngleOne, flowAngleTwo, frictionAngleTwo, giveCompressiveStrength(), giveTensileStrength(), and oofem::norm().
Referenced by computeJacobian(), and performRegularReturn().
| double oofem::LatticePlasticityDamage::computeYieldStress | ( | double | kappaP, |
| GaussPoint * | gp ) |
| double oofem::LatticePlasticityDamage::computeYieldValue | ( | const FloatArrayF< 3 > & | sigma, |
| const double | tempKappa, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 195 of file latticeplasticitydamage.C.
References oofem::FloatArrayF< N >::at(), computeHardening(), frictionAngleOne, frictionAngleTwo, giveCompressiveStrength(), giveTensileStrength(), and oofem::norm().
Referenced by performPlasticityReturn(), and performRegularReturn().
|
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::LatticeLinearElastic.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 188 of file latticeplasticitydamage.C.
References oofem::FEMComponent::domain.
|
overridevirtual |
Returns the value of material property 'aProperty'. Property must be identified by unique int id. Integration point also passed to allow for materials with spatially varying properties
| aProperty | ID of property requested. |
| gp | Integration point, |
Reimplemented from oofem::LatticeLinearElastic.
Definition at line 606 of file latticeplasticitydamage.C.
References fc_strength, ft_strength, oofem::LatticeLinearElastic::give(), oofem::RandomMaterialExtensionInterface::give(), and oofem::LatticeLinearElastic::giveStatus().
Referenced by giveCompressiveStrength(), and giveTensileStrength().
|
overridevirtual |
Reimplemented from oofem::LatticeLinearElastic.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 847 of file latticeplasticitydamage.C.
References oofem::LatticeLinearElastic::give3dLatticeStiffnessMatrix(), oofem::LatticeLinearElastic::giveStatus(), oofem::min(), and OOFEM_ERROR.
|
inlineoverridevirtual |
Reimplemented from oofem::LatticeLinearElastic.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 195 of file latticeplasticitydamage.h.
|
inlineprotectedvirtual |
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 187 of file latticeplasticitydamage.h.
References fc_strength, and give().
Referenced by computeDamageParam(), computeDMMatrix(), computeFVector(), computeMVector(), computeYieldValue(), oofem::LatticePlasticityDamageViscoelastic::giveCompressiveStrength(), performDamageEvaluation(), performPlasticityReturn(), and performRegularReturn().
|
inlineoverridevirtual |
Reimplemented from oofem::LatticeLinearElastic.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 194 of file latticeplasticitydamage.h.
References _IFT_LatticePlasticityDamage_Name.
|
protected |
|
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.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 863 of file latticeplasticitydamage.C.
References oofem::FloatArray::at(), oofem::GaussPoint::giveElement(), oofem::Material::giveIPValue(), oofem::LatticeStructuralElement::giveLength(), oofem::LatticeLinearElastic::giveStatus(), giveTensileStrength(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by oofem::LatticePlasticityDamageViscoelastic::giveIPValue().
|
protected |
|
protected |
|
overridevirtual |
Reimplemented from oofem::LatticeLinearElastic.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 624 of file latticeplasticitydamage.C.
References oofem::StructuralMaterial::computeStressIndependentStrainVector(), damageFlag, oofem::LatticeLinearElastic::giveStatus(), oofem::LatticePlasticityDamageStatus::initTempStatus(), performDamageEvaluation(), and performPlasticityReturn().
|
virtual |
Definition at line 347 of file latticeplasticitydamage.C.
References oofem::LatticeMaterialStatus::giveReducedLatticeStrain(), and oofem::LatticeLinearElastic::giveStatus().
Referenced by performPlasticityReturn().
|
inlineprotectedvirtual |
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 185 of file latticeplasticitydamage.h.
References ft_strength, and give().
Referenced by computeDamageParam(), computeDMMatrix(), computeFVector(), computeMVector(), computeYieldValue(), giveIPValue(), oofem::LatticePlasticityDamageViscoelastic::giveTensileStrength(), and performDamageEvaluation().
|
overridevirtual |
Tests if material supports material mode.
| mode | Required material mode. |
Reimplemented from oofem::LatticeLinearElastic.
Definition at line 59 of file latticeplasticitydamage.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::LatticeLinearElastic.
Reimplemented in oofem::LatticePlasticityDamageViscoelastic.
Definition at line 65 of file latticeplasticitydamage.C.
References _IFT_LatticePlasticityDamage_ahard, _IFT_LatticePlasticityDamage_angle1, _IFT_LatticePlasticityDamage_angle2, _IFT_LatticePlasticityDamage_damage, _IFT_LatticePlasticityDamage_fc, _IFT_LatticePlasticityDamage_flow, _IFT_LatticePlasticityDamage_ft, _IFT_LatticePlasticityDamage_ft1, _IFT_LatticePlasticityDamage_iter, _IFT_LatticePlasticityDamage_stype, _IFT_LatticePlasticityDamage_sub, _IFT_LatticePlasticityDamage_tol, _IFT_LatticePlasticityDamage_wf, _IFT_LatticePlasticityDamage_wf1, aHard, damageFlag, fc, flowAngleOne, flowAngleTwo, frictionAngleOne, frictionAngleTwo, ft, ftOneRatio, oofem::LatticeLinearElastic::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, newtonIter, numberOfSubIncrements, OOFEM_ERROR, OOFEM_WARNING, softeningType, wf, wfOne, and yieldTol.
Referenced by oofem::LatticePlasticityDamageViscoelastic::initializeFrom().
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::LatticeLinearElastic.
Definition at line 200 of file latticeplasticitydamage.h.
| void oofem::LatticePlasticityDamage::performDamageEvaluation | ( | GaussPoint * | gp, |
| FloatArrayF< 6 > & | reducedStrain, | ||
| TimeStep * | tStep ) const |
Definition at line 654 of file latticeplasticitydamage.C.
References computeDamageParam(), computeHardening(), oofem::LatticeLinearElastic::eNormalMean, frictionAngleOne, frictionAngleTwo, giveCompressiveStrength(), oofem::GaussPoint::giveElement(), oofem::LatticeStructuralElement::giveLength(), oofem::LatticeLinearElastic::giveStatus(), giveTensileStrength(), and oofem::norm().
Referenced by giveLatticeStress3d(), and oofem::LatticePlasticityDamageViscoelastic::giveLatticeStress3d().
| FloatArrayF< 6 > oofem::LatticePlasticityDamage::performPlasticityReturn | ( | GaussPoint * | gp, |
| const FloatArrayF< 6 > & | reducedStrain, | ||
| TimeStep * | tStep ) const |
Definition at line 355 of file latticeplasticitydamage.C.
References oofem::LatticeLinearElastic::alphaOne, oofem::LatticeLinearElastic::alphaTwo, oofem::assemble(), oofem::FloatArrayF< N >::at(), computeYieldValue(), oofem::LatticeLinearElastic::eNormalMean, giveCompressiveStrength(), oofem::GaussPoint::giveElement(), oofem::Element::giveGlobalNumber(), oofem::LatticePlasticityDamageStatus::giveKappaP(), giveReducedStrain(), oofem::LatticeLinearElastic::giveStatus(), oofem::mult(), numberOfSubIncrements, OOFEM_ERROR, OOFEM_LOG_INFO, performRegularReturn(), RR_Converged, RR_NotConverged, RR_Unknown, and yieldTol.
Referenced by giveLatticeStress3d(), and oofem::LatticePlasticityDamageViscoelastic::giveLatticeStress3d().
| double oofem::LatticePlasticityDamage::performRegularReturn | ( | FloatArrayF< 3 > & | stress, |
| LatticePlasticityDamage_ReturnResult & | returnResult, | ||
| double | yieldValue, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Definition at line 460 of file latticeplasticitydamage.C.
References oofem::LatticeLinearElastic::alphaOne, oofem::FloatArrayF< N >::at(), computeJacobian(), computeMVector(), computeYieldValue(), oofem::LatticeLinearElastic::eNormalMean, giveCompressiveStrength(), oofem::LatticeLinearElastic::giveStatus(), oofem::max(), newtonIter, oofem::norm(), RR_Converged, RR_NotConverged, oofem::solve_check(), and yieldTol.
Referenced by performPlasticityReturn().
|
protected |
hardening parameter
Definition at line 173 of file latticeplasticitydamage.h.
Referenced by computeDDHardeningDDKappa(), computeDHardeningDKappa(), computeHardening(), and initializeFrom().
|
protected |
damageFlag
Definition at line 183 of file latticeplasticitydamage.h.
Referenced by giveLatticeStress3d(), and initializeFrom().
|
protected |
compressive strength
Definition at line 149 of file latticeplasticitydamage.h.
Referenced by initializeFrom().
|
protected |
frictional angle of the plastic potential
Definition at line 155 of file latticeplasticitydamage.h.
Referenced by computeDMMatrix(), computeMVector(), and initializeFrom().
|
protected |
frictional angle of the plastic potential
Definition at line 158 of file latticeplasticitydamage.h.
Referenced by computeDMMatrix(), computeMVector(), and initializeFrom().
|
protected |
frictional angle of the yield surface
Definition at line 151 of file latticeplasticitydamage.h.
Referenced by computeDamageParam(), computeFVector(), computeYieldValue(), initializeFrom(), and performDamageEvaluation().
|
protected |
frictional angle of the yield surface
Definition at line 153 of file latticeplasticitydamage.h.
Referenced by computeDamageParam(), computeFVector(), computeMVector(), computeYieldValue(), initializeFrom(), and performDamageEvaluation().
|
protected |
tensile strength
Definition at line 147 of file latticeplasticitydamage.h.
Referenced by computeDMMatrix(), and initializeFrom().
|
protected |
ratio of tensile stress value for bilinear stress-crack opening curve
Definition at line 167 of file latticeplasticitydamage.h.
Referenced by computeDamageParam(), and initializeFrom().
|
protected |
Definition at line 144 of file latticeplasticitydamage.h.
|
protected |
maximum number of iterations for stress return
Definition at line 179 of file latticeplasticitydamage.h.
Referenced by initializeFrom(), and performRegularReturn().
|
protected |
Definition at line 180 of file latticeplasticitydamage.h.
Referenced by initializeFrom(), and performPlasticityReturn().
|
protected |
softening type determines the type of softening. 0 is exponential and 1 is bilinear.
Definition at line 164 of file latticeplasticitydamage.h.
Referenced by computeDamageParam(), and initializeFrom().
|
protected |
determines the softening -> corresponds to crack opening (not strain) when tension stress vanishes
Definition at line 161 of file latticeplasticitydamage.h.
Referenced by computeDamageParam(), and initializeFrom().
|
protected |
crack opening value for bilinear stress-crack opening curve
Definition at line 170 of file latticeplasticitydamage.h.
Referenced by computeDamageParam(), and initializeFrom().
|
protected |
yield tolerance
Definition at line 176 of file latticeplasticitydamage.h.
Referenced by initializeFrom(), performPlasticityReturn(), and performRegularReturn().