|
OOFEM 3.0
|
#include <isodamagemodel.h>
Public Member Functions | |
| IsotropicDamageMaterialStatus (GaussPoint *g) | |
| Constructor. | |
| void | printOutputAt (FILE *file, TimeStep *tStep) const override |
| Print receiver's output to given stream. | |
| double | giveKappa () const |
| Returns the last equilibrated scalar measure of the largest strain level. | |
| double | giveTempKappa () const |
| Returns the temp. scalar measure of the largest strain level. | |
| void | setTempKappa (double newKappa) |
| Sets the temp scalar measure of the largest strain level to given value. | |
| double | giveDamage () const |
| Returns the last equilibrated damage level. | |
| double | giveTempDamage () const |
| Returns the temp. damage level. | |
| void | setTempDamage (double newDamage) |
| Sets the temp damage level to given value. | |
| double | giveLe () const |
| Returns characteristic length stored in receiver. | |
| void | setLe (double ls) |
| Sets characteristic length to given value. | |
| double | giveCrackAngle () const |
| Returns crack angle stored in receiver. | |
| void | setCrackAngle (double ca) |
| Sets crack angle to given value. | |
| FloatArrayF< 3 > | giveCrackVector () const |
| Returns crack vector stored in receiver. This is useful for plotting cracks as a vector field (paraview etc.). | |
| void | setCrackVector (const FloatArrayF< 3 > &cv) |
| Sets crack vector to given value. This is useful for plotting cracks as a vector field (paraview etc.). | |
| double | giveStressWork () const |
| Returns the density of total work of stress on strain increments. | |
| double | giveTempStressWork () const |
| Returns the temp density of total work of stress on strain increments. | |
| void | setTempStressWork (double w) |
| Sets the density of total work of stress on strain increments to given value. | |
| double | giveDissWork () const |
| Returns the density of dissipated work. | |
| double | giveTempDissWork () const |
| Returns the density of temp dissipated work. | |
| void | setTempDissWork (double w) |
| Sets the density of dissipated work to given value. | |
| void | computeWork (GaussPoint *gp) |
| Computes the increment of total stress work and of dissipated work. | |
| const char * | giveClassName () const override |
| void | initTempStatus () override |
| void | updateYourself (TimeStep *tStep) override |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| Public Member Functions inherited from oofem::StructuralMaterialStatus | |
| StructuralMaterialStatus (GaussPoint *g) | |
| Constructor. Creates new StructuralMaterialStatus with IntegrationPoint g. | |
| void | printOutputAt (FILE *file, TimeStep *tStep) const override |
| Print receiver's output to given stream. | |
| void | initTempStatus () override |
| void | updateYourself (TimeStep *tStep) override |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| const FloatArray & | giveStrainVector () const |
| Returns the const pointer to receiver's strain vector. | |
| const FloatArray & | giveStressVector () const |
| Returns the const pointer to receiver's stress vector. | |
| const FloatArray & | givePVector () const |
| Returns the const pointer to receiver's first Piola-Kirchhoff stress vector. | |
| const FloatArray & | giveCVector () const |
| Returns the const pointer to receiver's Cauchy stress vector. | |
| const FloatArray & | giveFVector () const |
| Returns the const pointer to receiver's deformation gradient vector. | |
| const FloatArray & | giveTempStrainVector () const |
| Returns the const pointer to receiver's temporary strain vector. | |
| const FloatArray & | giveTempStressVector () const |
| Returns the const pointer to receiver's temporary stress vector. | |
| const FloatArray & | giveTempPVector () const |
| Returns the const pointer to receiver's temporary first Piola-Kirchhoff stress vector. | |
| const FloatArray & | giveTempCVector () const |
| Returns the const pointer to receiver's temporary Cauchy stress vector. | |
| const FloatArray & | giveTempFVector () const |
| Returns the const pointer to receiver's temporary deformation gradient vector. | |
| void | letStrainVectorBe (const FloatArray &v) |
| Assigns strain vector to given vector v. | |
| void | letStressVectorBe (const FloatArray &v) |
| Assigns stressVector to given vector v. | |
| void | letPVectorBe (const FloatArray &v) |
| Assigns PVector to given vector v. | |
| void | letCVectorBe (const FloatArray &v) |
| Assigns CVector to given vector v. | |
| void | letFVectorBe (const FloatArray &v) |
| Assigns FVector to given vector v. | |
| void | letTempStressVectorBe (const FloatArray &v) |
| Assigns tempStressVector to given vector v. | |
| void | letTempStrainVectorBe (const FloatArray &v) |
| Assigns tempStrainVector to given vector v. | |
| void | letTempPVectorBe (const FloatArray &v) |
| Assigns tempPVector to given vector v. | |
| void | letTempCVectorBe (const FloatArray &v) |
| Assigns tempPVector to given vector v. | |
| void | letTempFVectorBe (const FloatArray &v) |
| Assigns tempFVector to given vector v. | |
| const char * | giveClassName () const override |
| void | copyStateVariables (const MaterialStatus &iStatus) override |
| Functions for MaterialStatusMapperInterface. | |
| void | addStateVariables (const MaterialStatus &iStatus) override |
| Public Member Functions inherited from oofem::MaterialStatus | |
| MaterialStatus (GaussPoint *g) | |
| virtual bool | giveMaterialProperty (int propID, double &value) |
| virtual void | setMaterialProperty (int propID, double value) |
| Public Member Functions inherited from oofem::IntegrationPointStatus | |
| IntegrationPointStatus (GaussPoint *g) | |
| virtual | ~IntegrationPointStatus ()=default |
| Destructor. | |
| virtual void | setStatusVariable (int varID, double value) |
| virtual Interface * | giveInterface (InterfaceType t) |
| Public Member Functions inherited from oofem::MaterialStatusMapperInterface | |
| MaterialStatusMapperInterface () | |
| virtual | ~MaterialStatusMapperInterface () |
| virtual int | MSMI_map (const GaussPoint &iGP, const Domain &iOldDom, Set &sourceSet, const TimeStep &iTStep, MaterialStatus &oStatus) |
| virtual int | MSMI_map_cz (const GaussPoint &iGP, const Domain &iOldDom, Set &sourceSet, const TimeStep &iTStep, MaterialStatus &oStatus) |
| virtual int | MSMI_update (const GaussPoint &iGP, const TimeStep &iTStep) |
| virtual int | MSMI_finish (const TimeStep &iTStep) |
Protected Attributes | |
| double | kappa = 0. |
| Scalar measure of the largest strain level ever reached in material. | |
| double | tempKappa = 0. |
| Non-equilibrated scalar measure of the largest strain level. | |
| double | damage = 0. |
| Damage level of material. | |
| double | tempDamage = 0. |
| Non-equilibrated damage level of material. | |
| double | le = 0. |
| double | crack_angle = -1000.0 |
| Angle characterizing the crack direction. | |
| FloatArrayF< 3 > | crackVector |
| Crack orientation normalized to damage magnitude. This is useful for plotting cracks as a vector field (paraview etc.). | |
| double | stressWork = 0. |
| Density of total work done by stresses on strain increments. | |
| double | tempStressWork = 0. |
| Non-equilibrated density of total work done by stresses on strain increments. | |
| double | dissWork = 0. |
| Density of dissipated work. | |
| double | tempDissWork = 0. |
| Non-equilibrated density of dissipated work. | |
| Protected Attributes inherited from oofem::StructuralMaterialStatus | |
| FloatArray | strainVector |
| Equilibrated strain vector in reduced form. | |
| FloatArray | stressVector |
| Equilibrated stress vector in reduced form. | |
| FloatArray | tempStressVector |
| Temporary stress vector in reduced form (increments are used mainly in nonlinear analysis). | |
| FloatArray | tempStrainVector |
| Temporary strain vector in reduced form (to find balanced state). | |
| FloatArray | PVector |
| Equilibrated first Piola-Kirchhoff stress vector. | |
| FloatArray | tempPVector |
| Temporary first Piola-Kirchhoff stress vector (to find balanced state). | |
| FloatArray | CVector |
| Equilibrated Cauchy stress vector. | |
| FloatArray | tempCVector |
| Temporary Cauchy stress vector (to find balanced state). | |
| FloatArray | FVector |
| Equilibrated deformation gradient in reduced form. | |
| FloatArray | tempFVector |
| Temporary deformation gradient in reduced form (to find balanced state). | |
| Protected Attributes inherited from oofem::IntegrationPointStatus | |
| GaussPoint * | gp |
| Associated integration point. | |
| Protected Attributes inherited from oofem::MaterialStatusMapperInterface | |
| std::unique_ptr< MaterialMappingAlgorithm > | mpMaterialMapper |
This class implements associated Material Status to IsotropicDamageMaterial. Stores a scalar damage and hardening variable (and possible extra information).
Definition at line 62 of file isodamagemodel.h.
| oofem::IsotropicDamageMaterialStatus::IsotropicDamageMaterialStatus | ( | GaussPoint * | g | ) |
Constructor.
Definition at line 418 of file isodamagemodel.C.
References oofem::StructuralMaterialStatus::StructuralMaterialStatus().
Referenced by oofem::IsotropicDamageMaterial1Status::IsotropicDamageMaterial1Status().
| void oofem::IsotropicDamageMaterialStatus::computeWork | ( | GaussPoint * | gp | ) |
Computes the increment of total stress work and of dissipated work.
Definition at line 523 of file isodamagemodel.C.
References oofem::FloatArray::beDifferenceOf(), oofem::IntegrationPointStatus::gp, oofem::StructuralMaterialStatus::strainVector, oofem::StructuralMaterialStatus::stressVector, stressWork, tempDissWork, oofem::StructuralMaterialStatus::tempStrainVector, oofem::StructuralMaterialStatus::tempStressVector, and tempStressWork.
Referenced by oofem::IsotropicDamageMaterial::giveRealStressVector().
|
inlineoverridevirtual |
Implements oofem::IntegrationPointStatus.
Reimplemented in oofem::IsotropicGradientDamageMaterialStatus, oofem::MazarsMaterialStatus, and oofem::MazarsNLMaterialStatus.
Definition at line 144 of file isodamagemodel.h.
|
inline |
Returns crack angle stored in receiver.
Definition at line 119 of file isodamagemodel.h.
References crack_angle.
Referenced by oofem::IsotropicDamageMaterial::giveIPValue().
|
inline |
Returns crack vector stored in receiver. This is useful for plotting cracks as a vector field (paraview etc.).
Definition at line 123 of file isodamagemodel.h.
References crackVector, and damage.
Referenced by oofem::IsotropicDamageMaterial::giveIPValue().
|
inline |
Returns the last equilibrated damage level.
Definition at line 108 of file isodamagemodel.h.
References damage.
Referenced by oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_du(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_ud(), oofem::IsotropicDamageMaterial::giveIPValue(), oofem::IsotropicDamageMaterial::giveRealStressVector(), and oofem::MazarsMaterial::initDamaged().
|
inline |
Returns the density of dissipated work.
Definition at line 135 of file isodamagemodel.h.
References dissWork.
Referenced by oofem::IsotropicDamageMaterial::giveIPValue().
|
inline |
Returns the last equilibrated scalar measure of the largest strain level.
Definition at line 102 of file isodamagemodel.h.
References kappa.
Referenced by oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_dd_BN(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_dd_NN(), oofem::IsotropicDamageMaterial::giveIPValue(), and oofem::IsotropicDamageMaterial::giveRealStressVector().
|
inline |
Returns characteristic length stored in receiver.
Definition at line 115 of file isodamagemodel.h.
References le.
Referenced by oofem::IsotropicDamageMaterial1::computeDamageParamForCohesiveCrack(), oofem::MazarsMaterial::computeGt(), oofem::IsotropicDamageMaterial1::damageFunctionPrime(), and oofem::IsotropicDamageMaterial::giveIPValue().
|
inline |
Returns the density of total work of stress on strain increments.
Definition at line 129 of file isodamagemodel.h.
References stressWork.
Referenced by oofem::IsotropicDamageMaterial::giveIPValue().
|
inline |
Returns the temp. damage level.
Definition at line 110 of file isodamagemodel.h.
References tempDamage.
Referenced by oofem::IsotropicGradientDamageMaterial::computeEikonalInternalLength_a(), oofem::IsotropicGradientDamageMaterial::computeEikonalInternalLength_aPrime(), oofem::IsotropicGradientDamageMaterial::computeEikonalInternalLength_b(), oofem::IsotropicGradientDamageMaterial::computeEikonalInternalLength_bPrime(), oofem::IsotropicGradientDamageMaterial::computeInternalLength(), oofem::IsotropicDamageMaterial1::damageFunctionPrime(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_dd_BN(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_du(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_ud(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_uu(), oofem::IsotropicDamageMaterial::giveIPValue(), and oofem::IDNLMaterial::giveNonlocalMetricModifierAt().
|
inline |
Returns the density of temp dissipated work.
Definition at line 137 of file isodamagemodel.h.
References tempDissWork.
|
inline |
Returns the temp. scalar measure of the largest strain level.
Definition at line 104 of file isodamagemodel.h.
References tempKappa.
Referenced by oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_dd_BN(), oofem::IsotropicGradientDamageMaterial::giveGradientDamageStiffnessMatrix_dd_NN(), and oofem::IDNLMaterial::NonlocalMaterialStiffnessInterface_showSparseMtrxStructure().
|
inline |
Returns the temp density of total work of stress on strain increments.
Definition at line 131 of file isodamagemodel.h.
References tempStressWork.
|
overridevirtual |
Initializes the temporary internal variables, describing the current state according to previously reached equilibrium internal variables.
Reimplemented from oofem::MaterialStatus.
Reimplemented in oofem::IsotropicGradientDamageMaterialStatus, and oofem::MazarsNLMaterialStatus.
Definition at line 445 of file isodamagemodel.C.
References damage, dissWork, kappa, stressWork, tempDamage, tempDissWork, tempKappa, and tempStressWork.
|
overridevirtual |
Print receiver's output to given stream.
Reimplemented from oofem::IntegrationPointStatus.
Reimplemented in oofem::MazarsNLMaterialStatus.
Definition at line 424 of file isodamagemodel.C.
References crackVector, damage, dissWork, kappa, and stressWork.
|
overridevirtual |
Restores the receiver state previously written in stream.
| stream | Input stream. |
| mode | Determines amount of info available in stream (state, definition, ...). |
| throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::IntegrationPointStatus.
Reimplemented in oofem::MazarsMaterialStatus, and oofem::MazarsNLMaterialStatus.
Definition at line 497 of file isodamagemodel.C.
References oofem::CIO_IOERR, damage, dissWork, kappa, oofem::DataStream::read(), stressWork, and THROW_CIOERR.
|
overridevirtual |
Stores receiver state to output stream.
| stream | Output stream. |
| mode | Determines amount of info required in stream (state, definition, ...). |
| throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::IntegrationPointStatus.
Reimplemented in oofem::MazarsMaterialStatus, and oofem::MazarsNLMaterialStatus.
Definition at line 472 of file isodamagemodel.C.
References oofem::CIO_IOERR, damage, dissWork, kappa, stressWork, THROW_CIOERR, and oofem::DataStream::write().
|
inline |
Sets crack angle to given value.
Definition at line 121 of file isodamagemodel.h.
References crack_angle.
|
inline |
Sets crack vector to given value. This is useful for plotting cracks as a vector field (paraview etc.).
Definition at line 125 of file isodamagemodel.h.
References crackVector.
|
inline |
Sets characteristic length to given value.
Definition at line 117 of file isodamagemodel.h.
References le.
Referenced by oofem::MazarsMaterial::initDamaged(), and oofem::MazarsNLMaterial::initDamaged().
|
inline |
Sets the temp damage level to given value.
Definition at line 112 of file isodamagemodel.h.
References tempDamage.
Referenced by oofem::IsotropicDamageMaterial::giveRealStressVector(), and oofem::IsotropicDamageMaterial1::MMI_map().
|
inline |
Sets the density of dissipated work to given value.
Definition at line 139 of file isodamagemodel.h.
References tempDissWork.
|
inline |
Sets the temp scalar measure of the largest strain level to given value.
Definition at line 106 of file isodamagemodel.h.
References tempKappa.
Referenced by oofem::IsotropicDamageMaterial::giveRealStressVector(), and oofem::IsotropicDamageMaterial1::MMI_map().
|
inline |
Sets the density of total work of stress on strain increments to given value.
Definition at line 133 of file isodamagemodel.h.
References tempStressWork.
|
overridevirtual |
Update equilibrium history variables according to temp-variables. Invoked, after new equilibrium state has been reached.
Reimplemented from oofem::IntegrationPointStatus.
Reimplemented in oofem::IsotropicGradientDamageMaterialStatus, and oofem::MazarsNLMaterialStatus.
Definition at line 459 of file isodamagemodel.C.
References damage, dissWork, kappa, stressWork, tempDamage, tempDissWork, tempKappa, and tempStressWork.
Referenced by oofem::IsotropicDamageMaterial1::MMI_map().
|
protected |
Angle characterizing the crack direction.
Definition at line 80 of file isodamagemodel.h.
Referenced by giveCrackAngle(), and setCrackAngle().
|
protected |
Crack orientation normalized to damage magnitude. This is useful for plotting cracks as a vector field (paraview etc.).
Definition at line 82 of file isodamagemodel.h.
Referenced by giveCrackVector(), printOutputAt(), and setCrackVector().
|
protected |
Damage level of material.
Definition at line 70 of file isodamagemodel.h.
Referenced by giveCrackVector(), giveDamage(), initTempStatus(), oofem::IsotropicGradientDamageMaterialStatus::initTempStatus(), oofem::IDNLMaterialStatus::printOutputAt(), printOutputAt(), oofem::MazarsNLMaterialStatus::printOutputAt(), restoreContext(), saveContext(), and updateYourself().
|
protected |
Density of dissipated work.
Definition at line 90 of file isodamagemodel.h.
Referenced by giveDissWork(), initTempStatus(), oofem::IDNLMaterialStatus::printOutputAt(), printOutputAt(), restoreContext(), saveContext(), and updateYourself().
|
protected |
Scalar measure of the largest strain level ever reached in material.
Definition at line 66 of file isodamagemodel.h.
Referenced by giveKappa(), initTempStatus(), oofem::IDNLMaterialStatus::printOutputAt(), printOutputAt(), oofem::MazarsNLMaterialStatus::printOutputAt(), restoreContext(), saveContext(), and updateYourself().
|
protected |
Characteristic element length, computed when damage initialized from direction of maximum positive principal strain. Fixed during further loading.
Definition at line 78 of file isodamagemodel.h.
|
protected |
Density of total work done by stresses on strain increments.
Definition at line 86 of file isodamagemodel.h.
Referenced by computeWork(), giveStressWork(), initTempStatus(), oofem::IDNLMaterialStatus::printOutputAt(), printOutputAt(), restoreContext(), saveContext(), and updateYourself().
|
protected |
Non-equilibrated damage level of material.
Definition at line 72 of file isodamagemodel.h.
Referenced by giveTempDamage(), initTempStatus(), oofem::IsotropicGradientDamageMaterialStatus::initTempStatus(), setTempDamage(), and updateYourself().
|
protected |
Non-equilibrated density of dissipated work.
Definition at line 92 of file isodamagemodel.h.
Referenced by computeWork(), giveTempDissWork(), initTempStatus(), setTempDissWork(), and updateYourself().
|
protected |
Non-equilibrated scalar measure of the largest strain level.
Definition at line 68 of file isodamagemodel.h.
Referenced by giveTempKappa(), initTempStatus(), setTempKappa(), and updateYourself().
|
protected |
Non-equilibrated density of total work done by stresses on strain increments.
Definition at line 88 of file isodamagemodel.h.
Referenced by computeWork(), giveTempStressWork(), initTempStatus(), setTempStressWork(), and updateYourself().