OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::IsotropicDamageMaterial1 Class Reference

This class implements a simple local isotropic damage model for concrete in tension. More...

#include <idm1.h>

+ Inheritance diagram for oofem::IsotropicDamageMaterial1:
+ Collaboration diagram for oofem::IsotropicDamageMaterial1:

Public Member Functions

 IsotropicDamageMaterial1 (int n, Domain *d)
 Constructor. More...
 
virtual ~IsotropicDamageMaterial1 ()
 Destructor. More...
 
virtual const char * giveClassName () const
 
virtual const char * giveInputRecordName () const
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual void giveInputRecord (DynamicInputRecord &input)
 Setups the input record string of receiver. More...
 
bool isCrackBandApproachUsed ()
 
virtual void computeEquivalentStrain (double &kappa, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
 Computes the equivalent strain measure from given strain vector (full form). More...
 
virtual void computeEta (FloatArray &answer, const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
 Computes derivative of the equivalent strain with regards to strain. More...
 
virtual void computeDamageParam (double &omega, double kappa, const FloatArray &strain, GaussPoint *gp)
 Computes the value of damage parameter omega, based on given value of equivalent strain. More...
 
void computeDamageParamForCohesiveCrack (double &omega, double kappa, GaussPoint *gp)
 computes the value of damage parameter omega, based on a given value of equivalent strain, using iterations to achieve objectivity, based on the crack band concept (effective element size used) More...
 
double damageFunction (double kappa, GaussPoint *gp)
 Returns the value of damage parameter corresponding to a given value of the damage-driving variable kappa, depending on the type of selected damage law, using a simple dependence (no adjustment for element size). More...
 
double damageFunctionPrime (double kappa, GaussPoint *gp)
 Returns the value of compliance parameter corresponding to a given value of the damage-driving variable kappa, depending on the type of selected damage law, using a simple dependence (no adjustment for element size). More...
 
double complianceFunction (double kappa, GaussPoint *gp)
 Returns the value of compliance parameter corresponding to a given value of the damage-driving variable kappa, depending on the type of selected damage law, using a simple dependence (no adjustment for element size). More...
 
double evaluatePermanentStrain (double kappa, double omega)
 
virtual InterfacegiveInterface (InterfaceType it)
 Interface requesting service. More...
 
virtual int MMI_map (GaussPoint *gp, Domain *oldd, TimeStep *tStep)
 Maps the required internal state variables from old mesh oldd to given ip. More...
 
virtual int MMI_update (GaussPoint *gp, TimeStep *tStep, FloatArray *estrain=NULL)
 Updates the required internal state variables from previously mapped values. More...
 
virtual int MMI_finish (TimeStep *tStep)
 Finishes the mapping for given time step. More...
 
virtual MaterialStatusCreateStatus (GaussPoint *gp) const
 Creates new copy of associated status and inserts it into given integration point. More...
 
virtual MaterialStatusgiveStatus (GaussPoint *gp) const
 Returns material status of receiver in given integration point. More...
 
virtual double give (int aProperty, GaussPoint *gp)
 Returns the value of material property 'aProperty'. More...
 
virtual bool isCharacteristicMtrxSymmetric (MatResponseMode rMode)
 Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More...
 
- Public Member Functions inherited from oofem::IsotropicDamageMaterial
 IsotropicDamageMaterial (int n, Domain *d)
 Constructor. More...
 
virtual ~IsotropicDamageMaterial ()
 Destructor. More...
 
virtual int hasNonLinearBehaviour ()
 Returns nonzero if receiver is non linear. More...
 
virtual int hasMaterialModeCapability (MaterialMode mode)
 Tests if material supports material mode. More...
 
LinearElasticMaterialgiveLinearElasticMaterial ()
 Returns reference to undamaged (bulk) material. More...
 
virtual void give3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point. More...
 
virtual void giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
 Computes the real stress vector for given total strain and integration point. More...
 
virtual void giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More...
 
virtual void giveRealStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_3d. More...
 
virtual void giveRealStressVector_StressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
 Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· More...
 
virtual void giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. More...
 
virtual void giveThermalDilatationVector (FloatArray &answer, GaussPoint *, TimeStep *)
 Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis. More...
 
MaterialStatusCreateStatus (GaussPoint *gp) const
 Creates new copy of associated status and inserts it into given integration point. More...
 
- Public Member Functions inherited from oofem::StructuralMaterial
 StructuralMaterial (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralMaterial ()
 Destructor. More...
 
virtual void giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history. More...
 
virtual void giveRealStressVector_ShellStressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
 
virtual void giveRealStressVector_Warping (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_2dBeamLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_PlateLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_Fiber (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_Lattice2d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 
virtual void giveRealStressVector_Lattice3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 
virtual void giveRealStressVector_2dPlateSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation is not provided. More...
 
virtual void giveRealStressVector_3dBeamSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 
virtual void giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Prototype for computation of Eshelby stress. More...
 
void give_dPdF_from (const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp)
 
void convert_dSdE_2_dPdF (FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode)
 
double giveReferenceTemperature ()
 Returns the reference temperature of receiver. More...
 
virtual void computeStressIndependentStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
 Computes reduced strain vector in given integration point, generated by internal processes in material, which are independent on loading in particular integration point. More...
 
virtual void computeStressIndependentStrainVector_3d (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
 
virtual void give3dMaterialStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 
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)
 Method for subtracting from reduced space strain vector its stress-independent parts (caused by temperature, shrinkage, creep and possibly by other phenomena). More...
 
void giveStressDependentPartOfStrainVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode)
 
virtual int setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type)
 Sets the value of a certain variable at a given integration point to the given value. More...
 
virtual void give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d beam layer stiffness matrix of receiver. More...
 
virtual void givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d plate layer stiffness matrix of receiver. More...
 
virtual void giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d fiber stiffness matrix of receiver. More...
 
virtual void give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d lattice stiffness matrix of receiver. More...
 
virtual void give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 3d lattice stiffness matrix of receiver. More...
 
virtual void give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of plate subsoil model. More...
 
virtual void give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of beam3d subsoil model. More...
 
virtual void giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More...
 
virtual void giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
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 givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void give1dStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
- Public Member Functions inherited from oofem::Material
 Material (int n, Domain *d)
 Constructor. More...
 
virtual ~Material ()
 Destructor. More...
 
virtual bool hasProperty (int aProperty, GaussPoint *gp)
 Returns true if 'aProperty' exists on material. More...
 
virtual void modifyProperty (int aProperty, double value, GaussPoint *gp)
 Modify 'aProperty', which already exists on material. More...
 
double giveCastingTime ()
 
virtual bool isActivated (TimeStep *tStep)
 
virtual int hasCastingTimeSupport ()
 Tests if material supports casting time. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual contextIOResultType saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Stores integration point state to output stream. More...
 
virtual contextIOResultType restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Reads integration point state to output stream. More...
 
virtual int checkConsistency ()
 Allows programmer to test some internal data, before computation begins. More...
 
virtual int initMaterial (Element *element)
 Optional function to call specific procedures when initializing a material. More...
 
virtual int packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
 Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
 Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int estimatePackSize (DataStream &buff, GaussPoint *ip)
 Estimates the necessary pack size to hold all packed data of receiver. More...
 
virtual double predictRelativeComputationalCost (GaussPoint *gp)
 Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model. More...
 
virtual double predictRelativeRedistributionCost (GaussPoint *gp)
 Returns the relative redistribution cost of the receiver. More...
 
virtual void initTempStatus (GaussPoint *gp)
 Initializes temporary variables stored in integration point status at the beginning of new time step. More...
 
- Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
 Regular constructor, creates component with given number and belonging to given domain. More...
 
virtual ~FEMComponent ()
 Virtual destructor. More...
 
DomaingiveDomain () const
 
virtual void setDomain (Domain *d)
 Sets associated Domain. More...
 
int giveNumber () const
 
void setNumber (int num)
 Sets number of receiver. More...
 
virtual void updateLocalNumbering (EntityRenumberingFunctor &f)
 Local renumbering support. More...
 
virtual contextIOResultType saveContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Stores receiver state to output stream. More...
 
virtual contextIOResultType restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Restores the receiver state previously written in stream. More...
 
virtual void printOutputAt (FILE *file, TimeStep *tStep)
 Prints output of receiver to stream, for given time step. More...
 
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros). More...
 
- Public Member Functions inherited from oofem::RandomMaterialExtensionInterface
 RandomMaterialExtensionInterface ()
 Constructor. More...
 
virtual ~RandomMaterialExtensionInterface ()
 Destructor. More...
 
IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
void giveInputRecord (DynamicInputRecord &ir)
 
bool give (int key, GaussPoint *gp, double &value)
 Returns the property in associated status of given integration point if defined. More...
 
- Public Member Functions inherited from oofem::Interface
 Interface ()
 Constructor. More...
 
virtual ~Interface ()
 
- Public Member Functions inherited from oofem::MaterialModelMapperInterface
 MaterialModelMapperInterface ()
 Constructor. More...
 
virtual ~MaterialModelMapperInterface ()
 Destructor. More...
 

Static Public Member Functions

static void computeStrainInvariants (const FloatArray &strainVector, double &I1e, double &J2e)
 Computes invariants I1 and J2 of the strain tensor from the strain components stored in a vector. More...
 
- Static Public Member Functions inherited from oofem::StructuralMaterial
static int giveSymVI (int ind1, int ind2)
 
static int giveVI (int ind1, int ind2)
 
static int giveVoigtVectorMask (IntArray &answer, MaterialMode mmode)
 Returns a mask of the vector indicies corresponding to components in a general (non-symmetric) second order tensor of some stress/strain/deformation measure that performes work. More...
 
static int giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode)
 The same as giveVoigtVectorMask but returns a mask corresponding to a symmetric second order tensor. More...
 
static void giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode)
 Gives the inverted version of giveVoigtVectorMask. More...
 
static int giveSizeOfVoigtVector (MaterialMode mmode)
 Returns the size of reduced stress/strain vector according to given mode. More...
 
static int giveSizeOfVoigtSymVector (MaterialMode mmode)
 Returns the size of symmetric part of a reduced stress/strain vector according to given mode. More...
 
static void giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
 Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. More...
 
static void giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
 Converts the reduced deformation gradient Voigt vector (2nd order tensor). More...
 
static void giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. More...
 
static void giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. More...
 
static void giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. More...
 
static void giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode)
 Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More...
 
static void giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
 Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. More...
 
static void giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
 Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. More...
 
static void transformStrainVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &strainVector, bool transpose=false)
 Transforms 3d strain vector into another coordinate system. More...
 
static void transformStressVectorTo (FloatArray &answer, const FloatMatrix &base, const FloatArray &stressVector, bool transpose=false)
 Transforms 3d stress vector into another coordinate system. More...
 
static double computeVonMisesStress (const FloatArray *currentStress)
 Computes equivalent of von Mises stress. More...
 
static void giveStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 3d strain vector transformation matrix from standard vector transformation matrix. More...
 
static void give2DStrainVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 2d strain vector transformation matrix from standard vector transformation matrix. More...
 
static void giveStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 3d stress vector transformation matrix from standard vector transformation matrix. More...
 
static void givePlaneStressVectorTranformationMtrx (FloatMatrix &answer, const FloatMatrix &base, bool transpose=false)
 Computes 2d stress vector transformation matrix from standard vector transformation matrix. More...
 
static void sortPrincDirAndValCloseTo (FloatArray *pVal, FloatMatrix *pDir, FloatMatrix *toPDir)
 Method for sorting newly computed principal values (pVal) and corresponding principal directions (pDir) to be closed to some (often previous) principal directions (toPDir). More...
 
static void computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode)
 Common functions for convenience. More...
 
static void computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
 Computes principal values and directions of stress or strain vector. More...
 
static double computeDeviatoricVolumetricSplit (FloatArray &dev, const FloatArray &s)
 Computes split of receiver into deviatoric and volumetric part. More...
 
static void computeDeviatoricVolumetricSum (FloatArray &s, const FloatArray &dev, double mean)
 
static void applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu)
 
static void applyDeviatoricElasticCompliance (FloatArray &strain, const FloatArray &stress, double GModulus)
 
static void applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu)
 
static void applyDeviatoricElasticStiffness (FloatArray &stress, const FloatArray &strain, double GModulus)
 
static void applyElasticStiffness (FloatArray &stress, const FloatArray &strain, double EModulus, double nu)
 
static void applyElasticCompliance (FloatArray &strain, const FloatArray &stress, double EModulus, double nu)
 
static double computeStressNorm (const FloatArray &stress)
 
static double computeFirstInvariant (const FloatArray &s)
 
static double computeSecondStressInvariant (const FloatArray &s)
 
static double computeThirdStressInvariant (const FloatArray &s)
 
static double computeFirstCoordinate (const FloatArray &s)
 
static double computeSecondCoordinate (const FloatArray &s)
 
static double computeThirdCoordinate (const FloatArray &s)
 

Protected Types

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
}
 Type characterizing the algorithm used to compute equivalent strain measure. More...
 
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
}
 Type characterizing the formula for the damage law. More...
 
- Protected Types inherited from oofem::IsotropicDamageMaterial
enum  loaUnloCriterium { idm_strainLevelCR, idm_damageLevelCR }
 Variable controlling type of loading/unloading law, default set to idm_strainLevel defines the two two possibilities: More...
 

Protected Member Functions

virtual void initDamaged (double kappa, FloatArray &totalStrainVector, GaussPoint *gp)
 Performs initialization, when damage first appear. More...
 
- Protected Member Functions inherited from oofem::IsotropicDamageMaterial
virtual void givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane stress stiffness matrix of receiver. More...
 
virtual void givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane strain stiffness matrix of receiver. More...
 
virtual void give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d stiffness matrix of receiver. More...
 
- Protected Member Functions inherited from oofem::RandomMaterialExtensionInterface
void _generateStatusVariables (GaussPoint *) const
 Sets up (generates) the variables identified in randVariables array using generators given in randomVariableGenerators and stores them in given status. More...
 

Protected Attributes

double e0
 Equivalent strain at stress peak (or a similar parameter). More...
 
double ef
 Determines ductility -> corresponds to fracturing strain. More...
 
double wf
 Determines ductility -> corresponds to crack opening in the cohesive crack model. More...
 
double gf
 Determines the softening -> corresponds to the initial fracture energy. More...
 
double gft
 Determines the softening for the bilinear law -> corresponds to the total fracture energy. More...
 
double ek
 Determines the softening for the bilinear law -> corresponds to the strain at the knee point. More...
 
double wk
 Determines the softening for the bilinear law -> corresponds to the crack opening at the knee point. More...
 
double sk
 Determines the softening for the bilinear law -> corresponds to the stress at the knee point. More...
 
double c1
 Parameters used in Hordijk's softening law. More...
 
double c2
 
EquivStrainType equivStrainType
 Parameter specifying the definition of equivalent strain. More...
 
double k
 Parameter used in Mises definition of equivalent strain. More...
 
double griff_n
 Parameter used in Griffith's criterion. More...
 
int damageLaw
 Temporary parameter reading type of softening law, used in other isotropic damage material models. More...
 
SofteningType softType
 Parameter specifying the type of softening (damage law). More...
 
double At
 Parameters used in Mazars damage law. More...
 
double Bt
 
double md
 Parameter used in "smooth damage law". More...
 
double e1
 Parameters used if softType = 7 (extended smooth damage law) More...
 
double e2
 
double s1
 
double nd
 
int checkSnapBack
 Check possible snap back flag. More...
 
double ep
 auxiliary input variablesfor softType == ST_SmoothExtended More...
 
double ft
 
double ps_alpha
 Parameters used by the model with permanent strain. More...
 
double ps_H
 
ElementCharSizeMethod ecsMethod
 Method used for evaluation of characteristic element size. More...
 
SetsourceElemSet
 Cached source element set used to map internal variables (adaptivity), created on demand. More...
 
- Protected Attributes inherited from oofem::IsotropicDamageMaterial
double tempDillatCoeff
 Coefficient of thermal dilatation. More...
 
double maxOmega
 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. More...
 
int permStrain
 Indicator of the type of permanent strain formulation (0 = standard damage with no permanent strain) More...
 
LinearElasticMateriallinearElasticMaterial
 Reference to bulk (undamaged) material. More...
 
enum oofem::IsotropicDamageMaterial::loaUnloCriterium llcriteria
 
- Protected Attributes inherited from oofem::StructuralMaterial
double referenceTemperature
 Reference temperature (temperature, when material has been built into structure). More...
 
- Protected Attributes inherited from oofem::Material
Dictionary propertyDictionary
 Property dictionary. More...
 
double castingTime
 Casting time. More...
 
- Protected Attributes inherited from oofem::FEMComponent
int number
 Component number. More...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 
- Protected Attributes inherited from oofem::RandomMaterialExtensionInterface
IntArray randVariables
 Array of randomized variables (identified by a key). More...
 
IntArray randomVariableGenerators
 Array of generators id's for corresponding randomized variables. More...
 

Static Protected Attributes

static MMAContainingElementProjection mapper
 Mapper used to map internal variables in adaptivity. More...
 

Additional Inherited Members

- Static Public Attributes inherited from oofem::StructuralMaterial
static std::vector< std::vector< int > > vIindex
 Voigt index map. More...
 
static std::vector< std::vector< int > > svIndex
 Symmetric Voigt index map. More...
 

Detailed Description

This class implements a simple local isotropic damage model for concrete in tension.

A model 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 137 of file idm1.h.

Member Enumeration Documentation

Type characterizing the algorithm used to compute equivalent strain measure.

Note that the assigned numbers to enum values have to correspond to values used in initializeFrom to resolve EquivStrainType. If not, the consistency between initializeFrom and giveInputRecord methods is lost.

Enumerator
EST_Mazars 
EST_Rankine_Smooth 
EST_ElasticEnergy 
EST_Mises 
EST_Rankine_Standard 
EST_ElasticEnergyPositiveStress 
EST_ElasticEnergyPositiveStrain 
EST_Griffith 
EST_Unknown 

Definition at line 177 of file idm1.h.

Type characterizing the formula for the damage law.

For example, linear softening can be specified with fracturing strain or crack opening.

Enumerator
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 

Definition at line 203 of file idm1.h.

Constructor & Destructor Documentation

oofem::IsotropicDamageMaterial1::IsotropicDamageMaterial1 ( int  n,
Domain d 
)
oofem::IsotropicDamageMaterial1::~IsotropicDamageMaterial1 ( )
virtual

Destructor.

Definition at line 101 of file idm1.C.

References sourceElemSet.

Member Function Documentation

double oofem::IsotropicDamageMaterial1::complianceFunction ( double  kappa,
GaussPoint gp 
)

Returns the value of compliance parameter corresponding to a given value of the damage-driving variable kappa, depending on the type of selected damage law, using a simple dependence (no adjustment for element size).

The compliance parameter gamma is defined as gamma = omega/(1-omega) where omega is the damage.

Parameters
kappaEquivalent strain measure.
gpIntegration point.

Definition at line 1125 of file idm1.C.

References damageFunction().

Referenced by oofem::IDNLMaterial::updateBeforeNonlocAverage().

void oofem::IsotropicDamageMaterial1::computeDamageParam ( double &  omega,
double  kappa,
const FloatArray strain,
GaussPoint gp 
)
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.

Implements oofem::IsotropicDamageMaterial.

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

Definition at line 831 of file idm1.C.

References computeDamageParamForCohesiveCrack(), damageFunction(), isCrackBandApproachUsed(), softType, and ST_Disable_Damage.

Referenced by oofem::IDGMaterial::giveRealStressVectorGrad().

void oofem::IsotropicDamageMaterial1::computeDamageParamForCohesiveCrack ( double &  omega,
double  kappa,
GaussPoint gp 
)

computes the value of damage parameter omega, based on a given value of equivalent strain, using iterations to achieve objectivity, based on the crack band concept (effective element size used)

Parameters
[out]omegaContains the resulting damage.
kappaEquivalent strain measure.
gpIntegration point.

Definition at line 843 of file idm1.C.

References checkSnapBack, E, e0, e0_ID, ef, ek, ek_ID, gf, gf_ID, gft, gft_ID, oofem::Material::give(), give(), oofem::GaussPoint::giveElement(), oofem::IsotropicDamageMaterialStatus::giveLe(), oofem::IsotropicDamageMaterial::giveLinearElasticMaterial(), oofem::FEMComponent::giveNumber(), giveStatus(), IDM1_ITERATION_LIMIT, oofem::IsotropicDamageMaterial::maxOmega, OOFEM_ERROR, OOFEM_WARNING, sk, softType, ST_BiLinear_Cohesive_Crack, ST_Exponential_Cohesive_Crack, ST_Linear_Cohesive_Crack, wf, wf_ID, and wk.

Referenced by computeDamageParam().

void oofem::IsotropicDamageMaterial1::computeEta ( FloatArray answer,
const FloatArray strain,
GaussPoint gp,
TimeStep tStep 
)
virtual
void oofem::IsotropicDamageMaterial1::computeStrainInvariants ( const FloatArray strainVector,
double &  I1e,
double &  J2e 
)
static

Computes invariants I1 and J2 of the strain tensor from the strain components stored in a vector.

Parameters
strainVectorInput strain components.
[out]I1eOutput value of strain invariant I1.
[out]J2eOutput value of strain invariant J2.

Definition at line 821 of file idm1.C.

References oofem::FloatArray::at(), and s1.

Referenced by computeEquivalentStrain().

MaterialStatus * oofem::IsotropicDamageMaterial1::CreateStatus ( GaussPoint gp) const
virtual

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::IDNLMaterial, oofem::MazarsNLMaterial, oofem::MazarsMaterial, and oofem::IDGMaterial.

Definition at line 1348 of file idm1.C.

References oofem::FEMComponent::domain.

Referenced by giveStatus().

double oofem::IsotropicDamageMaterial1::damageFunction ( double  kappa,
GaussPoint gp 
)

Returns the value of damage parameter corresponding to a given value of the damage-driving variable kappa, depending on the type of selected damage law, using a simple dependence (no adjustment for element size).

Parameters
kappaEquivalent strain measure.
gpIntegration point.

Definition at line 968 of file idm1.C.

References At, Bt, c2, e0, e0_ID, e1, e2, ef, ef_ID, give(), md, nd, OOFEM_WARNING, oofem::IsotropicDamageMaterial::permStrain, ps_alpha, s1, softType, ST_DoubleExponential, ST_Exponential, ST_Linear, ST_Mazars, ST_ModPowerExponential, ST_PowerExponential, ST_Smooth, and ST_SmoothExtended.

Referenced by complianceFunction(), oofem::IDNLMaterial::computeDamageParam(), computeDamageParam(), and oofem::IDNLMaterial::updateBeforeNonlocAverage().

double oofem::IsotropicDamageMaterial1::damageFunctionPrime ( double  kappa,
GaussPoint gp 
)
virtual

Returns the value of compliance parameter corresponding to a given value of the damage-driving variable kappa, depending on the type of selected damage law, using a simple dependence (no adjustment for element size).

The compliance parameter gamma is defined as gamma = omega/(1-omega) where omega is the damage.

Parameters
kappaEquivalent strain measure.
gpIntegration point. 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.
kappaEquivalent strain measure.
gpIntegration point.

Reimplemented from oofem::IsotropicDamageMaterial.

Definition at line 1039 of file idm1.C.

References At, Bt, E, e0, e0_ID, e1, e2, ef, ef_ID, gf, oofem::Material::give(), give(), oofem::IsotropicDamageMaterialStatus::giveLe(), oofem::IsotropicDamageMaterial::giveLinearElasticMaterial(), giveStatus(), oofem::IsotropicDamageMaterialStatus::giveTempDamage(), oofem::isnan(), md, nd, OOFEM_ERROR, s1, softType, ST_Exponential, ST_Exponential_Cohesive_Crack, ST_Linear, ST_Linear_Cohesive_Crack, ST_Mazars, ST_Smooth, ST_SmoothExtended, and wf.

Referenced by oofem::IDGMaterial::give1dGprime(), oofem::IDGMaterial::givePlaneStrainGprime(), and oofem::IDGMaterial::givePlaneStressGprime().

double oofem::IsotropicDamageMaterial1::evaluatePermanentStrain ( double  kappa,
double  omega 
)
virtual

Reimplemented from oofem::IsotropicDamageMaterial.

Definition at line 1132 of file idm1.C.

References oofem::IsotropicDamageMaterial::permStrain, and ps_alpha.

double oofem::IsotropicDamageMaterial1::give ( int  aProperty,
GaussPoint gp 
)
virtual

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::IsotropicDamageMaterial.

Definition at line 1313 of file idm1.C.

References e0, e0_ID, ef, ef_ID, ek, ek_ID, gf, gf_ID, gft, gft_ID, oofem::IsotropicDamageMaterial::give(), giveStatus(), wf, and wf_ID.

Referenced by computeDamageParamForCohesiveCrack(), damageFunction(), damageFunctionPrime(), oofem::IDNLMaterial::giveLocalNonlocalStiffnessContribution(), initDamaged(), and oofem::IDNLMaterial::NonlocalMaterialStiffnessInterface_showSparseMtrxStructure().

virtual const char* oofem::IsotropicDamageMaterial1::giveClassName ( ) const
inlinevirtual
Returns
Class name of the receiver.

Reimplemented from oofem::IsotropicDamageMaterial.

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

Definition at line 255 of file idm1.h.

void oofem::IsotropicDamageMaterial1::giveInputRecord ( DynamicInputRecord input)
virtual

Setups the input record string of receiver.

Parameters
inputDynamic input record to be filled by receiver.

Reimplemented from oofem::IsotropicDamageMaterial.

Reimplemented in oofem::IDNLMaterial.

Definition at line 395 of file idm1.C.

References _IFT_IsotropicDamageMaterial1_At, _IFT_IsotropicDamageMaterial1_Bt, _IFT_IsotropicDamageMaterial1_damageLaw, _IFT_IsotropicDamageMaterial1_e0, _IFT_IsotropicDamageMaterial1_e1, _IFT_IsotropicDamageMaterial1_e2, _IFT_IsotropicDamageMaterial1_ecsm, _IFT_IsotropicDamageMaterial1_ef, _IFT_IsotropicDamageMaterial1_ek, _IFT_IsotropicDamageMaterial1_ep, _IFT_IsotropicDamageMaterial1_equivstraintype, _IFT_IsotropicDamageMaterial1_ft, _IFT_IsotropicDamageMaterial1_gf, _IFT_IsotropicDamageMaterial1_gft, _IFT_IsotropicDamageMaterial1_md, _IFT_IsotropicDamageMaterial1_nd, _IFT_IsotropicDamageMaterial1_sk, _IFT_IsotropicDamageMaterial1_wf, _IFT_IsotropicDamageMaterial1_wk, At, Bt, damageLaw, e0, e1, e2, ecsMethod, ef, ek, ep, equivStrainType, ft, gf, gft, oofem::LinearElasticMaterial::giveInputRecord(), oofem::RandomMaterialExtensionInterface::giveInputRecord(), oofem::MaterialMappingAlgorithm::giveInputRecord(), oofem::IsotropicDamageMaterial::giveInputRecord(), oofem::IsotropicDamageMaterial::linearElasticMaterial, mapper, md, nd, oofem::DynamicInputRecord::setField(), sk, softType, ST_BiLinear_Cohesive_Crack, ST_Exponential, ST_Exponential_Cohesive_Crack, ST_Linear_Cohesive_Crack, wf, and wk.

Referenced by oofem::IDNLMaterial::giveInputRecord().

virtual const char* oofem::IsotropicDamageMaterial1::giveInputRecordName ( ) const
inlinevirtual
Interface * oofem::IsotropicDamageMaterial1::giveInterface ( InterfaceType  t)
virtual

Interface requesting service.

See also
InterfaceType
Returns
Requested interface if implemented, otherwise NULL.

Reimplemented from oofem::FEMComponent.

Reimplemented in oofem::IDNLMaterial, oofem::MazarsNLMaterial, and oofem::IDGMaterial.

Definition at line 1337 of file idm1.C.

References oofem::MaterialModelMapperInterfaceType.

MaterialStatus * oofem::IsotropicDamageMaterial1::giveStatus ( GaussPoint gp) const
virtual

Returns material status of receiver in given integration point.

If status does not exist yet, it is created using CreateStatus member function.

Parameters
gpReturns reference to material status belonging to integration point gp.
Returns
Material status associated with given integration point.

Reimplemented from oofem::Material.

Definition at line 1354 of file idm1.C.

References oofem::RandomMaterialExtensionInterface::_generateStatusVariables(), CreateStatus(), oofem::GaussPoint::giveMaterialStatus(), oofem::FEMComponent::giveNumber(), and oofem::GaussPoint::setMaterialStatus().

Referenced by oofem::IDNLMaterial::computeAngleAndSigmaRatio(), computeDamageParamForCohesiveCrack(), oofem::MazarsNLMaterial::computeEquivalentStrain(), oofem::IDNLMaterial::computeEquivalentStrain(), oofem::MazarsMaterial::computeGc(), oofem::MazarsMaterial::computeGt(), damageFunctionPrime(), give(), oofem::IDGMaterial::give1dGprime(), oofem::IDGMaterial::give1dKappaMatrix(), oofem::IDGMaterial::give1dStressStiffMtrx(), oofem::IDGMaterial::giveInternalLength(), oofem::IDGMaterial::giveInternalLengthDerivative(), oofem::IDNLMaterial::giveIPValue(), oofem::IDNLMaterial::giveLocalNonlocalStiffnessContribution(), oofem::IDNLMaterial::giveNonlocalMetricModifierAt(), oofem::IDGMaterial::givePlaneStrainGprime(), oofem::IDGMaterial::givePlaneStrainKappaMatrix(), oofem::IDGMaterial::givePlaneStrainStiffMtrx(), oofem::IDGMaterial::givePlaneStressGprime(), oofem::IDGMaterial::givePlaneStressKappaMatrix(), oofem::IDGMaterial::givePlaneStressStiffMtrx(), oofem::IDGMaterial::giveRealStressVectorGrad(), oofem::IDNLMaterial::giveRemoteNonlocalStiffnessContribution(), oofem::MazarsMaterial::initDamaged(), oofem::MazarsNLMaterial::initDamaged(), initDamaged(), MMI_map(), MMI_update(), oofem::IDNLMaterial::NonlocalMaterialStiffnessInterface_addIPContribution(), oofem::IDNLMaterial::NonlocalMaterialStiffnessInterface_giveIntegrationDomainList(), oofem::IDNLMaterial::NonlocalMaterialStiffnessInterface_showSparseMtrxStructure(), oofem::MazarsNLMaterial::updateBeforeNonlocAverage(), and oofem::IDNLMaterial::updateBeforeNonlocAverage().

void oofem::IsotropicDamageMaterial1::initDamaged ( double  kappa,
FloatArray totalStrainVector,
GaussPoint gp 
)
protectedvirtual

Performs initialization, when damage first appear.

The characteristic length is computed from the direction of largest positive principal strain and stored in corresponding status.

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

Reimplemented from oofem::IsotropicDamageMaterial.

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

Definition at line 1159 of file idm1.C.

References oofem::FloatMatrix::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beColumnOf(), oofem::FloatArray::beProductOf(), checkSnapBack, oofem::StructuralMaterial::computePrincipalValDir(), oofem::FEMComponent::domain, E, e0, e0_ID, ecsMethod, ef, ef_ID, equivStrainType, EST_Griffith, gf, gf_ID, oofem::Material::give(), give(), oofem::Element::giveCharacteristicSize(), oofem::EngngModel::giveCurrentStep(), oofem::IsotropicDamageMaterialStatus::giveDamage(), oofem::GaussPoint::giveElement(), oofem::Domain::giveEngngModel(), oofem::StructuralMaterial::giveFullSymVectorForm(), oofem::FloatArray::giveIndexMaxElem(), oofem::FloatArray::giveIndexMinElem(), oofem::Element::giveLabel(), oofem::IsotropicDamageMaterialStatus::giveLe(), oofem::IsotropicDamageMaterial::giveLinearElasticMaterial(), oofem::GaussPoint::giveMaterialMode(), giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), isCrackBandApproachUsed(), M_PI, OOFEM_ERROR, OOFEM_WARNING, oofem::principal_strain, oofem::principal_stress, oofem::FloatArray::rotatedWith(), oofem::IsotropicDamageMaterialStatus::setCrackAngle(), oofem::IsotropicDamageMaterialStatus::setCrackVector(), oofem::IsotropicDamageMaterialStatus::setLe(), sk, softType, ST_BiLinear_Cohesive_Crack, ST_Disable_Damage, ST_Exponential_Cohesive_Crack, ST_Linear_Cohesive_Crack, wf, wf_ID, wk, and oofem::FloatMatrix::zero().

Referenced by oofem::IDGMaterial::giveRealStressVectorGrad().

IRResultType oofem::IsotropicDamageMaterial1::initializeFrom ( InputRecord ir)
virtual

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.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
Returns
IRResultType

Reimplemented from oofem::IsotropicDamageMaterial.

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

Definition at line 112 of file idm1.C.

References _IFT_IsotropicDamageMaterial1_alphaps, _IFT_IsotropicDamageMaterial1_At, _IFT_IsotropicDamageMaterial1_Bt, _IFT_IsotropicDamageMaterial1_c1, _IFT_IsotropicDamageMaterial1_c2, _IFT_IsotropicDamageMaterial1_checkSnapBack, _IFT_IsotropicDamageMaterial1_damageLaw, _IFT_IsotropicDamageMaterial1_e0, _IFT_IsotropicDamageMaterial1_e1, _IFT_IsotropicDamageMaterial1_e2, _IFT_IsotropicDamageMaterial1_ecsm, _IFT_IsotropicDamageMaterial1_ef, _IFT_IsotropicDamageMaterial1_ek, _IFT_IsotropicDamageMaterial1_ep, _IFT_IsotropicDamageMaterial1_equivstraintype, _IFT_IsotropicDamageMaterial1_ft, _IFT_IsotropicDamageMaterial1_gf, _IFT_IsotropicDamageMaterial1_gft, _IFT_IsotropicDamageMaterial1_k, _IFT_IsotropicDamageMaterial1_md, _IFT_IsotropicDamageMaterial1_n, _IFT_IsotropicDamageMaterial1_nd, _IFT_IsotropicDamageMaterial1_sk, _IFT_IsotropicDamageMaterial1_skft, _IFT_IsotropicDamageMaterial1_wf, _IFT_IsotropicDamageMaterial1_wk, _IFT_IsotropicDamageMaterial1_wkwf, _IFT_IsotropicLinearElasticMaterial_e, At, Bt, c1, c2, checkSnapBack, damageLaw, E, e0, e1, e2, oofem::ECSM_Oliver1, oofem::ECSM_Oliver1modified, oofem::ECSM_Projection, oofem::ECSM_ProjectionCentered, oofem::ECSM_SquareRootOfArea, ecsMethod, ef, ek, ep, equivStrainType, EST_ElasticEnergy, EST_ElasticEnergyPositiveStrain, EST_ElasticEnergyPositiveStress, EST_Griffith, EST_Mazars, EST_Mises, EST_Rankine_Smooth, EST_Rankine_Standard, ft, gf, gft, griff_n, oofem::InputRecord::hasField(), oofem::LinearElasticMaterial::initializeFrom(), oofem::RandomMaterialExtensionInterface::initializeFrom(), oofem::MaterialMappingAlgorithm::initializeFrom(), oofem::IsotropicDamageMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_BAD_FORMAT, oofem::IRRT_OK, k, oofem::IsotropicDamageMaterial::linearElasticMaterial, mapper, md, nd, OOFEM_ERROR, OOFEM_WARNING, oofem::IsotropicDamageMaterial::permStrain, ps_alpha, s1, sk, softType, ST_BiLinear_Cohesive_Crack, ST_Disable_Damage, ST_DoubleExponential, ST_Exponential, ST_Exponential_Cohesive_Crack, ST_Hordijk_Cohesive_Crack, ST_Linear, ST_Linear_Cohesive_Crack, ST_Mazars, ST_ModPowerExponential, ST_PowerExponential, ST_Smooth, ST_SmoothExtended, wf, and wk.

Referenced by oofem::IDGMaterial::initializeFrom(), and oofem::IDNLMaterial::initializeFrom().

virtual bool oofem::IsotropicDamageMaterial1::isCharacteristicMtrxSymmetric ( MatResponseMode  rMode)
inlinevirtual

Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.

Reimplemented from oofem::Material.

Definition at line 341 of file idm1.h.

bool oofem::IsotropicDamageMaterial1::isCrackBandApproachUsed ( )
inline
int oofem::IsotropicDamageMaterial1::MMI_finish ( TimeStep tStep)
virtual

Finishes the mapping for given time step.

Used to perform cleanup. Typically some mappers require to compute some global mesh data related to current step, which are valid for example to all IPs - so they are computed only once for all IPs, stored and they need to be deallocated. These mappers are typically class variables, but their finish is invoked by all members.

Implements oofem::MaterialModelMapperInterface.

Definition at line 1454 of file idm1.C.

References oofem::MMAContainingElementProjection::finish(), and mapper.

int oofem::IsotropicDamageMaterial1::MMI_update ( GaussPoint gp,
TimeStep tStep,
FloatArray elemGPVec = NULL 
)
virtual

Updates the required internal state variables from previously mapped values.

The result is stored in gp status. This map and update splitting is necessary, for example for nonlocal models that local quantity to be averaged must be mapped in all integration points and then update can happen, because it may depend on nonlocal variable, which is computed from local values.

Parameters
gpIntegration point belonging to new domain which values will be mapped.
tStepTime step.
elemGPVecVector passed to MMI_update at material level, probably computed from primary unknowns (for structural elements this represent strain vector).
Returns
Nonzero if o.k.

Implements oofem::MaterialModelMapperInterface.

Definition at line 1435 of file idm1.C.

References oofem::IsotropicDamageMaterial::giveRealStressVector(), giveStatus(), oofem::StructuralMaterialStatus::giveStrainVector(), and oofem::GaussPoint::updateYourself().

Member Data Documentation

double oofem::IsotropicDamageMaterial1::At
protected

Parameters used in Mazars damage law.

Definition at line 209 of file idm1.h.

Referenced by damageFunction(), damageFunctionPrime(), giveInputRecord(), and initializeFrom().

double oofem::IsotropicDamageMaterial1::Bt
protected

Definition at line 209 of file idm1.h.

Referenced by damageFunction(), damageFunctionPrime(), giveInputRecord(), and initializeFrom().

double oofem::IsotropicDamageMaterial1::c1
protected

Parameters used in Hordijk's softening law.

Definition at line 170 of file idm1.h.

Referenced by initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::c2
protected

Definition at line 170 of file idm1.h.

Referenced by damageFunction(), initializeFrom(), and IsotropicDamageMaterial1().

int oofem::IsotropicDamageMaterial1::checkSnapBack
protected

Check possible snap back flag.

Definition at line 216 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), initDamaged(), and initializeFrom().

int oofem::IsotropicDamageMaterial1::damageLaw
protected

Temporary parameter reading type of softening law, used in other isotropic damage material models.

Definition at line 198 of file idm1.h.

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

double oofem::IsotropicDamageMaterial1::e1
protected

Parameters used if softType = 7 (extended smooth damage law)

Definition at line 214 of file idm1.h.

Referenced by oofem::IDNLMaterial::computeAngleAndSigmaRatio(), damageFunction(), damageFunctionPrime(), giveInputRecord(), and initializeFrom().

double oofem::IsotropicDamageMaterial1::e2
protected
ElementCharSizeMethod oofem::IsotropicDamageMaterial1::ecsMethod
protected

Method used for evaluation of characteristic element size.

Definition at line 225 of file idm1.h.

Referenced by giveInputRecord(), initDamaged(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::ek
protected

Determines the softening for the bilinear law -> corresponds to the strain at the knee point.

Definition at line 161 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), give(), giveInputRecord(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::ep
protected

auxiliary input variablesfor softType == ST_SmoothExtended

Definition at line 219 of file idm1.h.

Referenced by giveInputRecord(), and initializeFrom().

EquivStrainType oofem::IsotropicDamageMaterial1::equivStrainType
protected
double oofem::IsotropicDamageMaterial1::ft
protected

Definition at line 219 of file idm1.h.

Referenced by giveInputRecord(), and initializeFrom().

double oofem::IsotropicDamageMaterial1::gf
protected

Determines the softening -> corresponds to the initial fracture energy.

For a linear law, it is the area under the stress/strain curve. For an exponential law, it is the area bounded by the elastic range and a tangent to the softening part of the curve at the peak stress. For a bilinear law, gf corresponds to area bounded by elasticity and the first linear softening line projected to zero stress.

Definition at line 155 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), damageFunctionPrime(), give(), giveInputRecord(), initDamaged(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::gft
protected

Determines the softening for the bilinear law -> corresponds to the total fracture energy.

Definition at line 158 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), give(), giveInputRecord(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::griff_n
protected

Parameter used in Griffith's criterion.

Definition at line 195 of file idm1.h.

Referenced by computeEquivalentStrain(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::k
protected

Parameter used in Mises definition of equivalent strain.

Definition at line 192 of file idm1.h.

Referenced by computeEquivalentStrain(), initializeFrom(), and IsotropicDamageMaterial1().

MMAContainingElementProjection oofem::IsotropicDamageMaterial1::mapper
staticprotected

Mapper used to map internal variables in adaptivity.

Definition at line 236 of file idm1.h.

Referenced by giveInputRecord(), oofem::IDGMaterial::initializeFrom(), oofem::MazarsMaterial::initializeFrom(), initializeFrom(), MMI_finish(), and MMI_map().

double oofem::IsotropicDamageMaterial1::md
protected

Parameter used in "smooth damage law".

Definition at line 211 of file idm1.h.

Referenced by damageFunction(), damageFunctionPrime(), giveInputRecord(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::nd
protected

Definition at line 214 of file idm1.h.

Referenced by damageFunction(), damageFunctionPrime(), giveInputRecord(), and initializeFrom().

double oofem::IsotropicDamageMaterial1::ps_alpha
protected

Parameters used by the model with permanent strain.

Definition at line 222 of file idm1.h.

Referenced by damageFunction(), evaluatePermanentStrain(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::ps_H
protected

Definition at line 222 of file idm1.h.

Referenced by IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::s1
protected
double oofem::IsotropicDamageMaterial1::sk
protected

Determines the softening for the bilinear law -> corresponds to the stress at the knee point.

Definition at line 167 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), giveInputRecord(), initDamaged(), initializeFrom(), and IsotropicDamageMaterial1().

SofteningType oofem::IsotropicDamageMaterial1::softType
protected

Parameter specifying the type of softening (damage law).

Definition at line 206 of file idm1.h.

Referenced by computeDamageParam(), computeDamageParamForCohesiveCrack(), damageFunction(), damageFunctionPrime(), giveInputRecord(), initDamaged(), initializeFrom(), and IsotropicDamageMaterial1().

Set* oofem::IsotropicDamageMaterial1::sourceElemSet
protected

Cached source element set used to map internal variables (adaptivity), created on demand.

Definition at line 228 of file idm1.h.

Referenced by IsotropicDamageMaterial1(), MMI_map(), and ~IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::wf
protected

Determines ductility -> corresponds to crack opening in the cohesive crack model.

Definition at line 147 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), damageFunctionPrime(), give(), giveInputRecord(), initDamaged(), initializeFrom(), and IsotropicDamageMaterial1().

double oofem::IsotropicDamageMaterial1::wk
protected

Determines the softening for the bilinear law -> corresponds to the crack opening at the knee point.

Definition at line 164 of file idm1.h.

Referenced by computeDamageParamForCohesiveCrack(), giveInputRecord(), initDamaged(), initializeFrom(), and IsotropicDamageMaterial1().


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

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:37 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011