OOFEM 3.0
Loading...
Searching...
No Matches
oofem::IsotropicDamageMaterial Class Referenceabstract

#include <isodamagemodel.h>

Inheritance diagram for oofem::IsotropicDamageMaterial:
Collaboration diagram for oofem::IsotropicDamageMaterial:

Public Member Functions

 IsotropicDamageMaterial (int n, Domain *d)
 Constructor.
virtual ~IsotropicDamageMaterial ()
 Destructor.
bool hasMaterialModeCapability (MaterialMode mode) const override
const char * giveClassName () const override
LinearElasticMaterialgiveLinearElasticMaterial () 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
virtual double evaluatePermanentStrain (double kappa, double omega) const
double give (int aProperty, GaussPoint *gp) const override
virtual double computeEquivalentStrain (const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const =0
virtual void computeEta (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep) const
virtual double computeDamageParam (double kappa, const FloatArray &strain, GaussPoint *gp) const =0
void initializeFrom (InputRecord &ir) override
void giveInputRecord (DynamicInputRecord &input) override
std::unique_ptr< MaterialStatusCreateStatus (GaussPoint *gp) const override
FloatMatrixF< 1, 1 > give1dStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) 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_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 bool isCharacteristicMtrxSymmetric (MatResponseMode rMode) const
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 MaterialStatusgiveStatus (GaussPoint *gp) const
virtual int packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
virtual int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
virtual int estimatePackSize (DataStream &buff, GaussPoint *ip)
virtual double predictRelativeComputationalCost (GaussPoint *gp)
virtual double predictRelativeRedistributionCost (GaussPoint *gp)
virtual void initTempStatus (GaussPoint *gp) const
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
virtual ~FEMComponent ()=default
 Virtual destructor.
virtual const char * giveInputRecordName () const =0
DomaingiveDomain () 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 InterfacegiveInterface (InterfaceType t)
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros).

Protected Types

enum  loaUnloCriterium { idm_strainLevelCR , idm_damageLevelCR }

Protected Member Functions

virtual void initDamaged (double kappa, FloatArray &totalStrainVector, GaussPoint *gp) const
virtual double damageFunctionPrime (double kappa, GaussPoint *gp) const
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 Attributes

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).
LinearElasticMateriallinearElasticMaterial = 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.
Domaindomain
 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 &currentStress)
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.

Detailed Description

Base class representing general isotropic damage model. It is based on isotropic damage concept, assuming that damage evolution law is postulated in explicit form, relation damage parameter (omega) to scalar measure of the largest strain level ever reached in material (kappa).

Definition at line 160 of file isodamagemodel.h.

Member Enumeration Documentation

◆ loaUnloCriterium

Variable controlling type of loading/unloading law, default set to idm_strainLevel defines the two two possibilities:

  • idm_strainLevelCR the unloading takes place, when strain level is smaller than the largest level ever reached;
  • idm_damageLevelCR the unloading takes place, when damage level is smaller than the largest damage ever reached;
Enumerator
idm_strainLevelCR 
idm_damageLevelCR 

Definition at line 180 of file isodamagemodel.h.

Constructor & Destructor Documentation

◆ IsotropicDamageMaterial()

◆ ~IsotropicDamageMaterial()

oofem::IsotropicDamageMaterial::~IsotropicDamageMaterial ( )
virtual

Destructor.

Definition at line 50 of file isodamagemodel.C.

References linearElasticMaterial.

Member Function Documentation

◆ computeDamageParam()

virtual double oofem::IsotropicDamageMaterial::computeDamageParam ( double kappa,
const FloatArray & strain,
GaussPoint * gp ) const
pure virtual

Computes the value of damage parameter omega, based on given value of equivalent strain.

Parameters
[out]omegaContains result.
kappaEquivalent strain measure.
strainTotal strain in full form.
gpIntegration point.

Implemented in oofem::IDNLMaterial, oofem::IsotropicDamageMaterial1, and oofem::MazarsMaterial.

Referenced by giveRealStressVector().

◆ computeEquivalentStrain()

virtual double oofem::IsotropicDamageMaterial::computeEquivalentStrain ( const FloatArray & strain,
GaussPoint * gp,
TimeStep * tStep ) const
pure virtual

Computes the equivalent strain measure from given strain vector (full form).

Parameters
[out]kappaReturn parameter, containing the corresponding equivalent strain.
strainTotal strain vector in full form.
gpIntegration point.
tStepTime step.

Implemented in oofem::IDNLMaterial, oofem::IsotropicDamageMaterial1, oofem::MazarsMaterial, and oofem::MazarsNLMaterial.

Referenced by giveRealStressVector().

◆ computeEta()

virtual void oofem::IsotropicDamageMaterial::computeEta ( FloatArray & answer,
const FloatArray & strain,
GaussPoint * gp,
TimeStep * tStep ) const
inlinevirtual

Computes derivative of the equivalent strain with regards to strain

Parameters
[out]answerContains the resulting derivative.
strainStrain vector.
gpIntegration point.
tStepTime step.

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 258 of file isodamagemodel.h.

References OOFEM_ERROR.

Referenced by give1dStressStiffMtrx(), and givePlaneStressStiffMtrx().

◆ CreateStatus()

std::unique_ptr< MaterialStatus > oofem::IsotropicDamageMaterial::CreateStatus ( GaussPoint * gp) const
inlineoverridevirtual

Creates new copy of associated status and inserts it into given integration point.

Parameters
gpIntegration point where newly created status will be stored.
Returns
Reference to new status.

Reimplemented from oofem::Material.

Reimplemented in oofem::IsotropicDamageMaterial1, oofem::IsotropicGradientDamageMaterial, oofem::MazarsMaterial, and oofem::MazarsNLMaterial.

Definition at line 271 of file isodamagemodel.h.

◆ damageFunctionPrime()

virtual double oofem::IsotropicDamageMaterial::damageFunctionPrime ( double kappa,
GaussPoint * gp ) const
inlineprotectedvirtual

Returns the value of derivative of damage function wrt damage-driving variable kappa corresponding to a given value of the kappa, depending on the type of selected damage law.

Parameters
kappaEquivalent strain measure.
gpIntegration point.

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 296 of file isodamagemodel.h.

References OOFEM_ERROR.

Referenced by give1dStressStiffMtrx(), and givePlaneStressStiffMtrx().

◆ evaluatePermanentStrain()

virtual double oofem::IsotropicDamageMaterial::evaluatePermanentStrain ( double kappa,
double omega ) const
inlinevirtual

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 233 of file isodamagemodel.h.

Referenced by giveIPValue(), and giveRealStressVector().

◆ give()

double oofem::IsotropicDamageMaterial::give ( int aProperty,
GaussPoint * gp ) const
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

Parameters
aPropertyID of property requested.
gpIntegration point,
Returns
Property value.

Reimplemented from oofem::Material.

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 367 of file isodamagemodel.C.

References linearElasticMaterial.

◆ give1dStressStiffMtrx()

FloatMatrixF< 1, 1 > oofem::IsotropicDamageMaterial::give1dStressStiffMtrx ( MatResponseMode mmode,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

Method for computing 1d stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).
Returns
Stiffness matrix.

Reimplemented from oofem::StructuralMaterial.

Definition at line 228 of file isodamagemodel.C.

References oofem::FloatMatrix::beDyadicProductOf(), computeEta(), damageFunctionPrime(), oofem::Material::giveStatus(), linearElasticMaterial, maxOmega, oofem::min(), oofem::FloatArray::times(), and oofem::FloatMatrix::times().

◆ give3dMaterialStiffnessMatrix()

FloatMatrixF< 6, 6 > oofem::IsotropicDamageMaterial::give3dMaterialStiffnessMatrix ( MatResponseMode mode,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.

Parameters
answerComputed results.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented from oofem::StructuralMaterial.

Definition at line 63 of file isodamagemodel.C.

References oofem::Material::giveStatus(), linearElasticMaterial, maxOmega, and oofem::min().

◆ giveClassName()

const char * oofem::IsotropicDamageMaterial::giveClassName ( ) const
inlineoverridevirtual
Returns
Class name of the receiver.

Implements oofem::FEMComponent.

Reimplemented in oofem::IsotropicDamageMaterial1, oofem::IsotropicGradientDamageMaterial, oofem::MazarsMaterial, and oofem::MazarsNLMaterial.

Definition at line 189 of file isodamagemodel.h.

◆ giveInputRecord()

void oofem::IsotropicDamageMaterial::giveInputRecord ( DynamicInputRecord & input)
overridevirtual

Setups the input record string of receiver.

Parameters
inputDynamic input record to be filled by receiver.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 390 of file isodamagemodel.C.

References _IFT_IsotropicDamageMaterial_maxOmega, _IFT_IsotropicDamageMaterial_talpha, maxOmega, oofem::DynamicInputRecord::setField(), and tempDillatCoeff.

◆ giveIPValue()

int oofem::IsotropicDamageMaterial::giveIPValue ( FloatArray & answer,
GaussPoint * gp,
InternalStateType type,
TimeStep * tStep )
overridevirtual

◆ giveLinearElasticMaterial()

◆ givePlaneStrainStiffMtrx()

FloatMatrixF< 4, 4 > oofem::IsotropicDamageMaterial::givePlaneStrainStiffMtrx ( MatResponseMode mmode,
GaussPoint * gp,
TimeStep * tStep ) const
overrideprotectedvirtual

Method for computing plane strain stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane strain stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix. Note: as already described, if zero strain component is imposed (Plane strain, ..) this condition must be taken into account in geometrical relations, and corresponding component has to be included in reduced vector. (So plane strain conditions are \( \epsilon_z = \gamma_{xz} = \gamma_{yz} = 0 \), but relations for \( \epsilon_z\) and \(\sigma_z\) are included).

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented from oofem::StructuralMaterial.

Definition at line 211 of file isodamagemodel.C.

References oofem::Material::giveStatus(), linearElasticMaterial, maxOmega, and oofem::min().

◆ givePlaneStressStiffMtrx()

FloatMatrixF< 3, 3 > oofem::IsotropicDamageMaterial::givePlaneStressStiffMtrx ( MatResponseMode mmode,
GaussPoint * gp,
TimeStep * tStep ) const
overrideprotectedvirtual

Method for computing plane stress stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane stress stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.

Parameters
answerStiffness matrix.
mmodeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented from oofem::StructuralMaterial.

Definition at line 170 of file isodamagemodel.C.

References oofem::FloatMatrix::beDyadicProductOf(), computeEta(), damageFunctionPrime(), oofem::Material::giveStatus(), linearElasticMaterial, maxOmega, oofem::min(), oofem::FloatArray::times(), and oofem::FloatMatrix::times().

◆ giveRealStressVector()

void oofem::IsotropicDamageMaterial::giveRealStressVector ( FloatArray & answer,
GaussPoint * gp,
const FloatArray & reducedStrain,
TimeStep * tStep ) const
overridevirtual

Computes the real stress vector for given total strain and integration point. The total strain is defined as strain computed directly from displacement field at given time. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process.

Parameters
answerStress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress.
gpIntegration point.
reducedStrainStrain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain.
tStepCurrent time step (most models are able to respond only when tStep is current time step).
Todo
Move this to StructuralCrossSection ?

Reimplemented from oofem::StructuralMaterial.

Definition at line 86 of file isodamagemodel.C.

References oofem::FloatArray::at(), computeDamageParam(), computeEquivalentStrain(), oofem::IsotropicDamageMaterialStatus::computeWork(), evaluatePermanentStrain(), oofem::IsotropicDamageMaterialStatus::giveDamage(), oofem::IsotropicDamageMaterialStatus::giveKappa(), giveLinearElasticMaterial(), oofem::FloatArray::giveSize(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), idm_damageLevelCR, idm_strainLevelCR, initDamaged(), oofem::Material::initTempStatus(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), llcriteria, OOFEM_ERROR, permStrain, oofem::IsotropicDamageMaterialStatus::setTempDamage(), oofem::IsotropicDamageMaterialStatus::setTempKappa(), and oofem::FloatMatrix::times().

Referenced by giveRealStressVector_1d(), giveRealStressVector_3d(), giveRealStressVector_PlaneStrain(), giveRealStressVector_PlaneStress(), giveRealStressVector_StressControl(), and oofem::IsotropicDamageMaterial1::MMI_update().

◆ giveRealStressVector_1d()

FloatArrayF< 1 > oofem::IsotropicDamageMaterial::giveRealStressVector_1d ( const FloatArrayF< 1 > & reducedE,
GaussPoint * gp,
TimeStep * tStep ) const
inlineoverridevirtual

Default implementation relies on giveRealStressVector_StressControl.

Reimplemented from oofem::StructuralMaterial.

Definition at line 223 of file isodamagemodel.h.

References giveRealStressVector(), and IsotropicDamageMaterial().

◆ giveRealStressVector_3d()

FloatArrayF< 6 > oofem::IsotropicDamageMaterial::giveRealStressVector_3d ( const FloatArrayF< 6 > & strain,
GaussPoint * gp,
TimeStep * tStep ) const
inlineoverridevirtual

Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.

Reimplemented from oofem::StructuralMaterial.

Definition at line 199 of file isodamagemodel.h.

References giveRealStressVector(), and IsotropicDamageMaterial().

◆ giveRealStressVector_PlaneStrain()

FloatArrayF< 4 > oofem::IsotropicDamageMaterial::giveRealStressVector_PlaneStrain ( const FloatArrayF< 4 > & strain,
GaussPoint * gp,
TimeStep * tStep ) const
inlineoverridevirtual

Default implementation relies on giveRealStressVector_3d.

Reimplemented from oofem::StructuralMaterial.

Definition at line 205 of file isodamagemodel.h.

References giveRealStressVector(), and IsotropicDamageMaterial().

◆ giveRealStressVector_PlaneStress()

FloatArrayF< 3 > oofem::IsotropicDamageMaterial::giveRealStressVector_PlaneStress ( const FloatArrayF< 3 > & reducedE,
GaussPoint * gp,
TimeStep * tStep ) const
inlineoverridevirtual

Default implementation relies on giveRealStressVector_StressControl.

Reimplemented from oofem::StructuralMaterial.

Definition at line 217 of file isodamagemodel.h.

References giveRealStressVector(), and IsotropicDamageMaterial().

◆ giveRealStressVector_StressControl()

FloatArray oofem::IsotropicDamageMaterial::giveRealStressVector_StressControl ( const FloatArray & reducedE,
const IntArray & strainControl,
GaussPoint * gp,
TimeStep * tStep ) const
inlineoverridevirtual

Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero·

Todo
We need a relative tolerance here! A relative tolerance like this could work, but if a really small increment is performed it won't work (it will be limited by machine precision)

Reimplemented from oofem::StructuralMaterial.

Definition at line 211 of file isodamagemodel.h.

References giveRealStressVector(), and IsotropicDamageMaterial().

◆ giveThermalDilatationVector()

FloatArrayF< 6 > oofem::IsotropicDamageMaterial::giveThermalDilatationVector ( GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis.

Parameters
answerVector of thermal dilatation coefficients.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented from oofem::StructuralMaterial.

Definition at line 351 of file isodamagemodel.C.

References tempDillatCoeff.

◆ hasMaterialModeCapability()

bool oofem::IsotropicDamageMaterial::hasMaterialModeCapability ( MaterialMode mode) const
overridevirtual

Tests if material supports material mode.

Parameters
modeRequired material mode.
Returns
Nonzero if supported, zero otherwise.

Reimplemented from oofem::Material.

Reimplemented in oofem::IsotropicGradientDamageMaterial.

Definition at line 56 of file isodamagemodel.C.

◆ initDamaged()

virtual void oofem::IsotropicDamageMaterial::initDamaged ( double kappa,
FloatArray & totalStrainVector,
GaussPoint * gp ) const
inlineprotectedvirtual

Abstract service allowing to perform some initialization, when damage first appear.

Parameters
kappaScalar measure of strain level.
totalStrainVectorCurrent total strain vector.
gpIntegration point.

Reimplemented in oofem::IDNLMaterial, oofem::IsotropicDamageMaterial1, oofem::MazarsMaterial, and oofem::MazarsNLMaterial.

Definition at line 286 of file isodamagemodel.h.

Referenced by giveRealStressVector().

◆ initializeFrom()

void oofem::IsotropicDamageMaterial::initializeFrom ( InputRecord & ir)
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.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
priorityPriority of the input record. This is used to determine the order of initialization

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::IsotropicDamageMaterial1, oofem::IsotropicGradientDamageMaterial, oofem::MazarsMaterial, and oofem::MazarsNLMaterial.

Definition at line 373 of file isodamagemodel.C.

References _IFT_IsotropicDamageMaterial_maxOmega, _IFT_IsotropicDamageMaterial_permstrain, _IFT_IsotropicDamageMaterial_talpha, IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::max(), maxOmega, oofem::min(), permStrain, and tempDillatCoeff.

◆ restoreContext()

void oofem::IsotropicDamageMaterial::restoreContext ( DataStream & stream,
ContextMode mode )
overridevirtual

Restores the receiver state previously written in stream.

See also
saveContext
Parameters
streamInput stream.
modeDetermines amount of info available in stream (state, definition, ...).
Exceptions
throwsan ContextIOERR exception if error encountered.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 409 of file isodamagemodel.C.

References CM_Definition, linearElasticMaterial, and oofem::FEMComponent::restoreContext().

◆ saveContext()

void oofem::IsotropicDamageMaterial::saveContext ( DataStream & stream,
ContextMode mode )
overridevirtual

Stores receiver state to output stream.

Parameters
streamOutput stream.
modeDetermines amount of info required in stream (state, definition, ...).
Exceptions
throwsan ContextIOERR exception if error encountered.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::IsotropicDamageMaterial1.

Definition at line 399 of file isodamagemodel.C.

References CM_Definition, linearElasticMaterial, and oofem::FEMComponent::saveContext().

Member Data Documentation

◆ linearElasticMaterial

◆ llcriteria

◆ maxOmega

double oofem::IsotropicDamageMaterial::maxOmega = 0.999999
protected

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.

Definition at line 167 of file isodamagemodel.h.

Referenced by oofem::IsotropicDamageMaterial1::computeDamageParamForCohesiveCrack(), give1dStressStiffMtrx(), give3dMaterialStiffnessMatrix(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_uu(), giveInputRecord(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), and initializeFrom().

◆ permStrain

int oofem::IsotropicDamageMaterial::permStrain = 0
protected

Indicator of the type of permanent strain formulation (0 = standard damage with no permanent strain).

Definition at line 170 of file isodamagemodel.h.

Referenced by oofem::IsotropicDamageMaterial1::damageFunction(), oofem::IsotropicDamageMaterial1::evaluatePermanentStrain(), giveIPValue(), giveRealStressVector(), oofem::IsotropicDamageMaterial1::initializeFrom(), and initializeFrom().

◆ tempDillatCoeff

double oofem::IsotropicDamageMaterial::tempDillatCoeff = 0.
protected

Coefficient of thermal dilatation.

Definition at line 164 of file isodamagemodel.h.

Referenced by giveInputRecord(), giveThermalDilatationVector(), and initializeFrom().


The documentation for this class was generated from the following files:

This page is part of the OOFEM-3.0 documentation. Copyright Copyright (C) 1994-2025 Borek Patzak Bořek Patzák
Project e-mail: oofem@fsv.cvut.cz
Generated at for OOFEM by doxygen 1.15.0 written by Dimitri van Heesch, © 1997-2011