|
OOFEM 3.0
|
#include <idmgrad.h>
Public Member Functions | |
| IsotropicGradientDamageMaterial (int n, Domain *d) | |
| Constructor. | |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| void | initializeFrom (InputRecord &ir) override |
| Interface * | giveInterface (InterfaceType t) override |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| void | giveGradientDamageStiffnessMatrix_uu (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override |
| Left upper block. | |
| void | giveGradientDamageStiffnessMatrix_du (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override |
| Left lower block. | |
| void | giveGradientDamageStiffnessMatrix_ud (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override |
| Right upper block. | |
| void | giveGradientDamageStiffnessMatrix_dd_NN (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override |
| Right lower block. | |
| void | giveGradientDamageStiffnessMatrix_dd_BB (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override |
| void | giveGradientDamageStiffnessMatrix_dd_BN (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override |
| void | giveRealStressVectorGradientDamage (FloatArray &answer1, double &answer2, GaussPoint *gp, const FloatArray &totalStrain, double nonlocalCumulatedStrain, TimeStep *tStep) override |
| gradient - based giveRealStressVector | |
| void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const override |
| void | computeLocalDamageDrivingVariable (double &answer, GaussPoint *gp, TimeStep *tStep) override |
| void | giveNonlocalInternalForces_N_factor (double &answer, double nlddv, GaussPoint *gp, TimeStep *tStep) override |
| void | giveNonlocalInternalForces_B_factor (FloatArray &answer, const FloatArray &nlddv, GaussPoint *gp, TimeStep *tStep) override |
| Public Member Functions inherited from oofem::IsotropicDamageMaterial1 | |
| IsotropicDamageMaterial1 (int n, Domain *d) | |
| Constructor. | |
| virtual | ~IsotropicDamageMaterial1 () |
| Destructor. | |
| void | giveInputRecord (DynamicInputRecord &input) override |
| bool | isCrackBandApproachUsed () const |
| double | computeEquivalentStrain (const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const override |
| void | computeEta (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const override |
| double | computeDamageParam (double kappa, const FloatArray &strain, GaussPoint *gp) const override |
| double | computeDamageParamForCohesiveCrack (double kappa, GaussPoint *gp) const |
| double | damageFunction (double kappa, GaussPoint *gp) const |
| double | damageFunctionPrime (double kappa, GaussPoint *gp) const override |
| double | complianceFunction (double kappa, GaussPoint *gp) const |
| double | evaluatePermanentStrain (double kappa, double omega) const override |
| int | MMI_map (GaussPoint *gp, Domain *oldd, TimeStep *tStep) override |
| int | MMI_update (GaussPoint *gp, TimeStep *tStep, FloatArray *estrain=nullptr) override |
| int | MMI_finish (TimeStep *tStep) override |
| MaterialStatus * | giveStatus (GaussPoint *gp) const override |
| double | give (int aProperty, GaussPoint *gp) const override |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| Public Member Functions inherited from oofem::IsotropicDamageMaterial | |
| IsotropicDamageMaterial (int n, Domain *d) | |
| Constructor. | |
| virtual | ~IsotropicDamageMaterial () |
| Destructor. | |
| LinearElasticMaterial * | giveLinearElasticMaterial () const |
| Returns reference to undamaged (bulk) material. | |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override |
| void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, 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. | |
| FloatArray | giveRealStressVector_StressControl (const FloatArray &strain, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const override |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| 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. | |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) 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 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< 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 () |
| Public Member Functions inherited from oofem::MaterialModelMapperInterface | |
| MaterialModelMapperInterface () | |
| Constructor. | |
| virtual | ~MaterialModelMapperInterface () |
| Destructor. | |
Protected Types | |
| enum | GradientDamageFormulationType { GDFT_Standard = 0 , GDFT_DecreasingInteractions = 1 , GDFT_Eikonal = 2 } |
| Protected Types inherited from oofem::IsotropicDamageMaterial1 | |
| enum | EquivStrainType { EST_Mazars =0 , EST_Rankine_Smooth =1 , EST_ElasticEnergy =2 , EST_Mises =3 , EST_Rankine_Standard =4 , EST_ElasticEnergyPositiveStress =5 , EST_ElasticEnergyPositiveStrain =6 , EST_Griffith =7 , EST_Unknown = 100 } |
| enum | SofteningType { ST_Unknown , ST_Exponential , ST_Linear , ST_Mazars , ST_Smooth , ST_SmoothExtended , ST_Exponential_Cohesive_Crack , ST_Linear_Cohesive_Crack , ST_BiLinear_Cohesive_Crack , ST_Disable_Damage , ST_PowerExponential , ST_DoubleExponential , ST_Hordijk_Cohesive_Crack , ST_ModPowerExponential , ST_Trilinear_Cohesive_Crack } |
| Protected Types inherited from oofem::IsotropicDamageMaterial | |
| enum | loaUnloCriterium { idm_strainLevelCR , idm_damageLevelCR } |
Protected Member Functions | |
| double | computeInternalLength (GaussPoint *gp) |
| int | giveDimension (GaussPoint *gp) |
| double | computeEikonalInternalLength_a (GaussPoint *gp) |
| double | computeEikonalInternalLength_b (GaussPoint *gp) |
| double | computeEikonalInternalLength_aPrime (GaussPoint *gp) |
| double | computeEikonalInternalLength_bPrime (GaussPoint *gp) |
| Protected Member Functions inherited from oofem::IsotropicDamageMaterial1 | |
| void | initDamaged (double kappa, FloatArray &totalStrainVector, GaussPoint *gp) const override |
| Protected Member Functions inherited from oofem::IsotropicDamageMaterial | |
| FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 4, 4 > | givePlaneStrainStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| Protected Member Functions inherited from oofem::RandomMaterialExtensionInterface | |
| void | _generateStatusVariables (GaussPoint *) const |
Protected Attributes | |
| GradientDamageFormulationType | gradientDamageFormulationType = GDFT_Standard |
| double | di_rho = 0. |
| double | di_eta = 0. |
| Protected Attributes inherited from oofem::IsotropicDamageMaterial1 | |
| double | e0 = 0. |
| Equivalent strain at stress peak (or a similar parameter). | |
| double | ef = 0. |
| Determines ductility -> corresponds to fracturing strain. | |
| double | wf = 0. |
| Determines ductility -> corresponds to crack opening in the cohesive crack model. | |
| double | gf = 0. |
| double | gft = 0. |
| Determines the softening for the bilinear law -> corresponds to the total fracture energy. | |
| double | ek = 0. |
| Determines the softening for the bilinear law -> corresponds to the strain at the knee point. | |
| double | wk = 0. |
| Determines the softening for the bilinear law -> corresponds to the crack opening at the knee point. | |
| double | sk = 0. |
| Determines the softening for the bilinear law -> corresponds to the stress at the knee point. | |
| double | c1 = 3. |
| Parameters used in Hordijk's softening law. | |
| double | c2 = 6.93 |
| double | w_k = 0. |
| Parameters used in Trilinear_Cohesive_Crack softening law. | |
| double | w_r = 0. |
| double | w_f = 0. |
| double | f_k = 0. |
| double | f_r = 0. |
| EquivStrainType | equivStrainType = EST_Unknown |
| Parameter specifying the definition of equivalent strain. | |
| double | k = 0. |
| Parameter used in Mises definition of equivalent strain. | |
| double | griff_n = 8. |
| Parameter used in Griffith's criterion. | |
| int | damageLaw = 0 |
| Temporary parameter reading type of softening law, used in other isotropic damage material models. | |
| SofteningType | softType = ST_Unknown |
| Parameter specifying the type of softening (damage law). | |
| double | At = 0. |
| Parameters used in Mazars damage law. | |
| double | Bt = 0. |
| double | md = 1. |
| Parameter used in "smooth damage law". | |
| double | e1 = 0. |
| Parameters used if softType = 7 (extended smooth damage law). | |
| double | e2 = 0. |
| double | s1 = 0. |
| double | nd = 0. |
| int | checkSnapBack = 0 |
| Check possible snap back flag. | |
| double | ep = 0. |
| auxiliary input variablesfor softType == ST_SmoothExtended | |
| double | ft = 0. |
| double | ps_alpha = 0. |
| Parameters used by the model with permanent strain. | |
| double | ps_H = 0. |
| ElementCharSizeMethod | ecsMethod = ECSM_Unknown |
| Method used for evaluation of characteristic element size. | |
| Set * | sourceElemSet = nullptr |
| Cached source element set used to map internal variables (adaptivity), created on demand. | |
| Protected Attributes inherited from oofem::IsotropicDamageMaterial | |
| double | tempDillatCoeff = 0. |
| Coefficient of thermal dilatation. | |
| double | maxOmega = 0.999999 |
| Maximum limit on omega. The purpose is elimination of a too compliant material which may cause convergence problems. Set to something like 0.99 if needed. | |
| int | permStrain = 0 |
| Indicator of the type of permanent strain formulation (0 = standard damage with no permanent strain). | |
| LinearElasticMaterial * | linearElasticMaterial = nullptr |
| Reference to bulk (undamaged) material. | |
| enum oofem::IsotropicDamageMaterial::loaUnloCriterium | llcriteria = idm_strainLevelCR |
| 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::IsotropicDamageMaterial1 | |
| static void | computeStrainInvariants (const FloatArray &strainVector, double &I1e, double &J2e) |
| 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. | |
| Static Protected Attributes inherited from oofem::IsotropicDamageMaterial1 | |
| static MMAContainingElementProjection | mapper |
| Mapper used to map internal variables in adaptivity. | |
| Private Member Functions inherited from oofem::GradientDamageMaterialExtensionInterface | |
| GradientDamageMaterialExtensionInterface (Domain *d) | |
| virtual | ~GradientDamageMaterialExtensionInterface () |
| Destructor. | |
| virtual void | giveFirstPKStressVectorGradientDamage (FloatArray &answer1, double &answer2, GaussPoint *gp, const FloatArray &totalStrain, double nonlocalDamageDrivningVariable, TimeStep *tStep) |
| virtual void | giveCauchyStressVectorGradientDamage (FloatArray &answer1, double &answer2, GaussPoint *gp, const FloatArray &totalStrain, double nonlocalDamageDrivningVariable, TimeStep *tStep) |
| virtual void | computeInternalForcesRegularizationTerm (double &answer, GaussPoint *gp, TimeStep *tStep) |
| virtual void | computeStiffnessRegularizationTerm (double &answer, GaussPoint *gp, TimeStep *tStep) |
| Private Member Functions inherited from oofem::Interface | |
| Interface () | |
| Constructor. | |
| virtual | ~Interface () |
| Private Attributes inherited from oofem::GradientDamageMaterialExtensionInterface | |
| Domain * | dom = nullptr |
| double | internalLength = 0. |
Gradient-enhanced Isotropic Damage model for concrete in tension,
Type characterizing the dependence of the internal lenght on variable of the state Note that the assigned numbers to enum values have to correspond to values used in initializeFrom to resolve internalLenghtDependence. If not, the consistency between initializeFrom and giveInputRecord methods is lost.
| Enumerator | |
|---|---|
| GDFT_Standard | |
| GDFT_DecreasingInteractions | |
| GDFT_Eikonal | |
| oofem::IsotropicGradientDamageMaterial::IsotropicGradientDamageMaterial | ( | int | n, |
| Domain * | d ) |
Constructor.
Definition at line 53 of file idmgrad.C.
References oofem::GradientDamageMaterialExtensionInterface::GradientDamageMaterialExtensionInterface(), and oofem::IsotropicDamageMaterial1::IsotropicDamageMaterial1().
|
protected |
Definition at line 265 of file idmgrad.C.
References oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), and oofem::GradientDamageMaterialExtensionInterface::internalLength.
Referenced by giveGradientDamageStiffnessMatrix_dd_NN(), giveGradientDamageStiffnessMatrix_du(), and giveNonlocalInternalForces_N_factor().
|
protected |
Definition at line 284 of file idmgrad.C.
References oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), and oofem::GradientDamageMaterialExtensionInterface::internalLength.
Referenced by giveGradientDamageStiffnessMatrix_dd_NN().
|
protected |
Definition at line 274 of file idmgrad.C.
References oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), and oofem::GradientDamageMaterialExtensionInterface::internalLength.
Referenced by giveGradientDamageStiffnessMatrix_dd_BB(), and giveNonlocalInternalForces_B_factor().
|
protected |
Definition at line 293 of file idmgrad.C.
References oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), and oofem::GradientDamageMaterialExtensionInterface::internalLength.
Referenced by giveGradientDamageStiffnessMatrix_dd_BN().
|
protected |
Definition at line 413 of file idmgrad.C.
References di_eta, di_rho, GDFT_DecreasingInteractions, GDFT_Standard, oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), gradientDamageFormulationType, oofem::GradientDamageMaterialExtensionInterface::internalLength, and OOFEM_WARNING.
Referenced by giveGradientDamageStiffnessMatrix_dd_BB(), giveGradientDamageStiffnessMatrix_dd_BN(), and giveNonlocalInternalForces_B_factor().
|
inlineoverridevirtual |
Implements oofem::GradientDamageMaterialExtensionInterface.
|
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::IsotropicDamageMaterial1.
|
inlineoverridevirtual |
Reimplemented from oofem::IsotropicDamageMaterial1.
|
protected |
Definition at line 369 of file idmgrad.C.
References oofem::GaussPoint::giveMaterialMode().
Referenced by giveGradientDamageStiffnessMatrix_dd_BB().
|
overridevirtual |
Implements oofem::GradientDamageMaterialExtensionInterface.
Definition at line 304 of file idmgrad.C.
References oofem::FloatMatrix::beUnitMatrix(), computeEikonalInternalLength_b(), computeInternalLength(), GDFT_DecreasingInteractions, GDFT_Eikonal, GDFT_Standard, giveDimension(), gradientDamageFormulationType, oofem::GradientDamageMaterialExtensionInterface::internalLength, OOFEM_WARNING, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::times().
|
overridevirtual |
Reimplemented from oofem::GradientDamageMaterialExtensionInterface.
Definition at line 325 of file idmgrad.C.
References oofem::FloatMatrix::clear(), computeEikonalInternalLength_bPrime(), computeInternalLength(), oofem::IsotropicDamageMaterial1::damageFunctionPrime(), di_eta, di_rho, oofem::FloatMatrix::fromArray(), GDFT_DecreasingInteractions, GDFT_Eikonal, GDFT_Standard, oofem::IsotropicDamageMaterialStatus::giveKappa(), oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), oofem::IsotropicDamageMaterialStatus::giveTempKappa(), oofem::GradientDamageMaterialStatusExtensionInterface::giveTempNonlocalDamageDrivingVariableGrad(), gradientDamageFormulationType, oofem::GradientDamageMaterialExtensionInterface::internalLength, OOFEM_WARNING, and oofem::FloatMatrix::times().
|
overridevirtual |
Right lower block.
Reimplemented from oofem::GradientDamageMaterialExtensionInterface.
Definition at line 236 of file idmgrad.C.
References oofem::FloatMatrix::at(), oofem::FloatMatrix::clear(), computeEikonalInternalLength_a(), computeEikonalInternalLength_aPrime(), oofem::IsotropicDamageMaterial1::damageFunctionPrime(), GDFT_Eikonal, oofem::IsotropicDamageMaterialStatus::giveKappa(), oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempKappa(), oofem::GradientDamageMaterialStatusExtensionInterface::giveTempLocalDamageDrivingVariable(), oofem::GradientDamageMaterialStatusExtensionInterface::giveTempNonlocalDamageDrivingVariable(), gradientDamageFormulationType, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
overridevirtual |
Left lower block.
Implements oofem::GradientDamageMaterialExtensionInterface.
Definition at line 188 of file idmgrad.C.
References computeEikonalInternalLength_a(), oofem::IsotropicDamageMaterial1::computeEta(), oofem::FloatMatrix::fromArray(), GDFT_Eikonal, oofem::IsotropicDamageMaterialStatus::giveDamage(), oofem::GaussPoint::giveMaterialMode(), oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), oofem::StructuralMaterialStatus::giveTempStrainVector(), gradientDamageFormulationType, and oofem::FloatMatrix::times().
|
overridevirtual |
Right upper block.
Implements oofem::GradientDamageMaterialExtensionInterface.
Definition at line 145 of file idmgrad.C.
References oofem::IsotropicDamageMaterial1::damageFunctionPrime(), oofem::FloatMatrix::fromArray(), oofem::IsotropicDamageMaterialStatus::giveDamage(), oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), oofem::GradientDamageMaterialStatusExtensionInterface::giveTempNonlocalDamageDrivingVariable(), oofem::StructuralMaterialStatus::giveTempStressVector(), oofem::FloatArray::times(), and oofem::FloatMatrix::times().
|
overridevirtual |
Left upper block.
Implements oofem::GradientDamageMaterialExtensionInterface.
Definition at line 107 of file idmgrad.C.
References oofem::IsotropicDamageMaterial::giveLinearElasticMaterial(), oofem::IsotropicDamageMaterial1::giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), oofem::IsotropicDamageMaterial::maxOmega, oofem::min(), and oofem::FloatMatrix::times().
|
inlineoverridevirtual |
Reimplemented from oofem::IsotropicDamageMaterial1.
Definition at line 79 of file idmgrad.h.
References _IFT_IsotropicGradientDamageMaterial_Name.
|
inlineoverridevirtual |
Interface requesting service.
Reimplemented from oofem::IsotropicDamageMaterial1.
Definition at line 82 of file idmgrad.h.
References oofem::GradientDamageMaterialExtensionInterface::GradientDamageMaterialExtensionInterface(), and oofem::GradientDamageMaterialExtensionInterfaceType.
|
overridevirtual |
Implements oofem::GradientDamageMaterialExtensionInterface.
Definition at line 399 of file idmgrad.C.
References computeEikonalInternalLength_b(), computeInternalLength(), GDFT_Eikonal, gradientDamageFormulationType, and oofem::FloatArray::times().
|
overridevirtual |
Implements oofem::GradientDamageMaterialExtensionInterface.
Definition at line 385 of file idmgrad.C.
References computeEikonalInternalLength_a(), GDFT_Eikonal, oofem::IsotropicDamageMaterial1::giveStatus(), oofem::GradientDamageMaterialStatusExtensionInterface::giveTempLocalDamageDrivingVariable(), and gradientDamageFormulationType.
|
overridevirtual |
gradient - based giveRealStressVector
Reimplemented from oofem::GradientDamageMaterialExtensionInterface.
Definition at line 430 of file idmgrad.C.
References oofem::FloatArray::beProductOf(), oofem::IsotropicDamageMaterial1::computeDamageParam(), oofem::IsotropicDamageMaterial1::computeEquivalentStrain(), oofem::IsotropicDamageMaterial::giveLinearElasticMaterial(), oofem::IsotropicDamageMaterial1::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::IsotropicDamageMaterial::idm_damageLevelCR, oofem::IsotropicDamageMaterial::idm_strainLevelCR, oofem::IsotropicDamageMaterial1::initDamaged(), oofem::Material::initTempStatus(), oofem::IsotropicDamageMaterial::llcriteria, OOFEM_ERROR, and oofem::FloatMatrix::times().
|
overridevirtual |
Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history. The algorithm should use temporary or equilibrium history variables stored in integration point status to compute and return required result.
| answer | Contains result. |
| 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 96 of file idmgrad.C.
References OOFEM_ERROR.
|
overridevirtual |
Tests if material supports material mode.
| mode | Required material mode. |
Reimplemented from oofem::IsotropicDamageMaterial.
|
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::IsotropicDamageMaterial1.
Definition at line 60 of file idmgrad.C.
References _IFT_IsotropicGradientDamageMaterial_di_eta, _IFT_IsotropicGradientDamageMaterial_di_rho, _IFT_IsotropicGradientDamageMaterial_formulationType, di_eta, di_rho, GDFT_DecreasingInteractions, GDFT_Eikonal, GDFT_Standard, gradientDamageFormulationType, IR_GIVE_OPTIONAL_FIELD, and oofem::IsotropicDamageMaterial1::mapper.
|
protected |
Definition at line 69 of file idmgrad.h.
Referenced by computeInternalLength(), giveGradientDamageStiffnessMatrix_dd_BN(), and initializeFrom().
|
protected |
Definition at line 68 of file idmgrad.h.
Referenced by computeInternalLength(), giveGradientDamageStiffnessMatrix_dd_BN(), and initializeFrom().
|
protected |
Definition at line 66 of file idmgrad.h.
Referenced by computeInternalLength(), giveGradientDamageStiffnessMatrix_dd_BB(), giveGradientDamageStiffnessMatrix_dd_BN(), giveGradientDamageStiffnessMatrix_dd_NN(), giveGradientDamageStiffnessMatrix_du(), giveNonlocalInternalForces_B_factor(), giveNonlocalInternalForces_N_factor(), and initializeFrom().