OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Class representing anisotropic damage model. More...
#include <anisodamagemodel.h>
Public Member Functions | |
AnisotropicDamageMaterial (int n, Domain *d) | |
Constructor. More... | |
virtual | ~AnisotropicDamageMaterial () |
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... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
IsotropicLinearElasticMaterial * | giveLinearElasticMaterial () |
Returns reference to undamaged (bulk) material. More... | |
virtual void | giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Plane-stress version of the stress evaluation algorithm. More... | |
void | computePrincValDir2D (double &D1, double &D2, double &c, double &s, double Dx, double Dy, double Dxy) |
bool | checkPrincVal2D (double Dx, double Dy, double Dxy) |
void | computeDamage (FloatMatrix &tempDamage, const FloatMatrix &damage, double kappa, double eps1, double eps2, double ceps, double seps, double epsZ) |
double | computeTraceD (double equivStrain) |
double | computeOutOfPlaneStrain (const FloatArray &inplaneStrain, const FloatMatrix &dam, bool tens_flag) |
double | computeDimensionlessOutOfPlaneStress (const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam) |
void | computeInplaneStress (FloatArray &inplaneStress, const FloatArray &inplaneStrain, double epsZ, const FloatMatrix &dam) |
double | obtainAlpha1 (FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, double damageThreshold) |
Obtains the proportion of the damage tensor that is needed to get the first eigenvalue equal to the damage threshold. More... | |
double | obtainAlpha2 (FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatMatrix ProjMatrix, double damageThreshold) |
Obtains the proportion of the damage tensor that is needed to get the second eigenvalue equal to the damage threshold. More... | |
double | obtainAlpha3 (FloatMatrix tempDamageTensor, double deltaLambda, FloatMatrix positiveStrainTensor, FloatArray vec3, double damageThreshold) |
Obtains the proportion of the damage tensor that is needed to get the third eigenvalue equal to the damage threshold. More... | |
double | checkSymmetry (FloatMatrix matrix) |
void | correctBigValues (FloatMatrix &matrix) |
double | computeTraceD (FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp) |
double | computeCorrectionFactor (FloatMatrix tempDamageTensor, FloatMatrix strainTensor, GaussPoint *gp) |
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_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 *atTime) |
Returns the integration point corresponding value in Reduced form. More... | |
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 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... | |
virtual void | computeDamageTensor (FloatMatrix &answer, GaussPoint *gp, const FloatArray &totalStrain, double equivStrain, TimeStep *atTime) |
MaterialStatus * | CreateStatus (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) |
virtual void | giveThermalDilatationVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep) |
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis. More... | |
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 | isCharacteristicMtrxSymmetric (MatResponseMode rMode) |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More... | |
virtual double | give (int aProperty, GaussPoint *gp) |
Returns the value of material property 'aProperty'. 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 MaterialStatus * | giveStatus (GaussPoint *gp) const |
Returns material status of receiver in given integration point. 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... | |
Domain * | giveDomain () 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... | |
virtual Interface * | giveInterface (InterfaceType t) |
Interface requesting service. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Types | |
enum | EquivStrainType { EST_Unknown, EST_Mazars, EST_Rankine_Smooth, EST_Rankine_Standard, EST_ElasticEnergy, EST_ElasticEnergyPositiveStress, EST_ElasticEnergyPositiveStrain, EST_Mises, EST_Griffith } |
Type characterizing the algorithm used to compute equivalent strain measure. More... | |
enum | DamLawType { DLT_Unknown, DLT_Desmorat1, DLT_Desmorat2, DLT_Linear, DLT_Exponential } |
Type characterizing the damage law. More... | |
Protected Member Functions | |
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... | |
virtual void | computePlaneStressStrain (FloatMatrix &answer, FloatMatrix damageTensor, FloatArray totalStrain, GaussPoint *gp, TimeStep *atTime) |
virtual void | computePlaneStressSigmaZ (double &answer, FloatMatrix damageTensor, FloatArray reducedTotalStrainVector, double epsZ, GaussPoint *gp, TimeStep *atTime) |
virtual void | computeSecantOperator (FloatMatrix &answer, FloatMatrix strainTensor, FloatMatrix damageTensor, GaussPoint *gp) |
double | computeK (GaussPoint *gp) |
double | computeKappa (FloatMatrix damageTensor) |
Protected Attributes | |
IsotropicLinearElasticMaterial * | linearElasticMaterial |
Reference to bulk (undamaged) material. More... | |
double | E |
Young's modulus. More... | |
double | nu |
Poisson's ratio. More... | |
double | kappa0 |
Damage threshold kappa0, as defined in the paper mentioned above. More... | |
double | kappaf |
Damage parameter kappa_f (in the paper denoted as "a") More... | |
double | aA |
Damage parameter a*A, needed to obtain Kappa(trD), according to eq. 33 in the paper mentioned above. More... | |
EquivStrainType | equivStrainType |
Parameter specifying the definition of equivalent strain. More... | |
DamLawType | damageLawType |
Parameter specifying the damage law. More... | |
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... | |
Domain * | domain |
Link to domain object, useful for communicating with other FEM components. More... | |
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 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) |
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... | |
Class representing anisotropic damage model.
Based on the paper : "Nonlocal anisotropic damage model and related computational aspects for quasi-brittle materials" by R. Desmorat, F. Gatuingt, and F. Ragueneau Plane stress case implemented by the algorithm described in a report by M. Jirasek & F. Suarez, 25 April 2014
Definition at line 179 of file anisodamagemodel.h.
|
protected |
Type characterizing the damage law.
Enumerator | |
---|---|
DLT_Unknown | |
DLT_Desmorat1 | |
DLT_Desmorat2 | |
DLT_Linear | |
DLT_Exponential |
Definition at line 200 of file anisodamagemodel.h.
|
protected |
Type characterizing the algorithm used to compute equivalent strain measure.
Enumerator | |
---|---|
EST_Unknown | |
EST_Mazars | |
EST_Rankine_Smooth | |
EST_Rankine_Standard | |
EST_ElasticEnergy | |
EST_ElasticEnergyPositiveStress | |
EST_ElasticEnergyPositiveStrain | |
EST_Mises | |
EST_Griffith |
Definition at line 196 of file anisodamagemodel.h.
oofem::AnisotropicDamageMaterial::AnisotropicDamageMaterial | ( | int | n, |
Domain * | d | ||
) |
Constructor.
Definition at line 53 of file anisodamagemodel.C.
References aA, damageLawType, DLT_Unknown, E, equivStrainType, EST_Unknown, kappa0, kappaf, linearElasticMaterial, and nu.
|
virtual |
bool oofem::AnisotropicDamageMaterial::checkPrincVal2D | ( | double | Dx, |
double | Dy, | ||
double | Dxy | ||
) |
Definition at line 234 of file anisodamagemodel.C.
Referenced by computeDamage().
double oofem::AnisotropicDamageMaterial::checkSymmetry | ( | FloatMatrix | matrix | ) |
Definition at line 1048 of file anisodamagemodel.C.
References oofem::FloatMatrix::at(), and oofem::FloatMatrix::giveNumberOfRows().
double oofem::AnisotropicDamageMaterial::computeCorrectionFactor | ( | FloatMatrix | tempDamageTensor, |
FloatMatrix | strainTensor, | ||
GaussPoint * | gp | ||
) |
Definition at line 1116 of file anisodamagemodel.C.
References oofem::FloatMatrix::at(), oofem::AnisotropicDamageMaterialStatus::giveFlag(), oofem::Material::giveStatus(), oofem::AnisotropicDamageMaterialStatus::giveStoredFactor(), oofem::AnisotropicDamageMaterialStatus::setTempFlag(), and oofem::AnisotropicDamageMaterialStatus::setTempStoredFactor().
Referenced by giveRealStressVector().
void oofem::AnisotropicDamageMaterial::computeDamage | ( | FloatMatrix & | tempDamage, |
const FloatMatrix & | damage, | ||
double | kappa, | ||
double | eps1, | ||
double | eps2, | ||
double | ceps, | ||
double | seps, | ||
double | epsZ | ||
) |
Definition at line 249 of file anisodamagemodel.C.
References oofem::FloatMatrix::at(), checkPrincVal2D(), computePrincValDir2D(), computeTraceD(), and oofem::macbra().
Referenced by giveRealStressVector_PlaneStress().
|
virtual |
Definition at line 1442 of file anisodamagemodel.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beProductOf(), computeKappa(), computePlaneStressStrain(), computeTraceD(), correctBigValues(), oofem::AnisotropicDamageMaterialStatus::giveDamage(), oofem::StructuralMaterial::giveFullSymVectorForm(), oofem::GaussPoint::giveMaterialMode(), oofem::Material::giveStatus(), oofem::AnisotropicDamageMaterialStatus::giveTempStrainZ(), oofem::FloatMatrix::jaco_(), obtainAlpha1(), obtainAlpha2(), obtainAlpha3(), oofem::FloatMatrix::resize(), oofem::FloatMatrix::subtract(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
Referenced by giveRealStressVector().
double oofem::AnisotropicDamageMaterial::computeDimensionlessOutOfPlaneStress | ( | const FloatArray & | inplaneStrain, |
double | epsZ, | ||
const FloatMatrix & | dam | ||
) |
Definition at line 393 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computePrincValDir2D(), oofem::macbra(), and nu.
Referenced by giveRealStressVector_PlaneStress().
|
virtual |
Computes the equivalent strain measure from given strain vector (full form).
[out] | kappa | Return parameter, containing the corresponding equivalent strain. |
strain | Total strain vector in full form. | |
gp | Integration point. | |
tStep | Time step. |
Definition at line 670 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::StructuralMaterial::computePrincipalValues(), equivStrainType, EST_Mazars, oofem::StructuralMaterial::giveFullSymVectorForm(), oofem::GaussPoint::giveMaterialMode(), oofem::Material::giveStatus(), oofem::AnisotropicDamageMaterialStatus::giveTempStrainZ(), oofem::FloatArray::isEmpty(), nu, OOFEM_ERROR, and oofem::principal_strain.
Referenced by giveRealStressVector().
void oofem::AnisotropicDamageMaterial::computeInplaneStress | ( | FloatArray & | inplaneStress, |
const FloatArray & | inplaneStrain, | ||
double | epsZ, | ||
const FloatMatrix & | dam | ||
) |
Definition at line 430 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computePrincValDir2D(), E, oofem::macbra(), nu, and oofem::FloatArray::resize().
Referenced by giveRealStressVector_PlaneStress().
|
protected |
|
protected |
Definition at line 814 of file anisodamagemodel.C.
References oofem::FloatMatrix::giveTrace(), kappa0, and kappaf.
Referenced by computeDamageTensor(), and giveRealStressVector().
double oofem::AnisotropicDamageMaterial::computeOutOfPlaneStrain | ( | const FloatArray & | inplaneStrain, |
const FloatMatrix & | dam, | ||
bool | tens_flag | ||
) |
Definition at line 357 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computePrincValDir2D(), oofem::macbra(), and nu.
Referenced by giveRealStressVector_PlaneStress().
|
protectedvirtual |
Definition at line 1364 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beProductOf(), E, giveLinearElasticMaterial(), oofem::FloatMatrix::giveTrace(), oofem::FloatMatrix::jaco_(), nu, oofem::FloatArray::resize(), and oofem::FloatMatrix::resize().
|
protectedvirtual |
Definition at line 1281 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatArray::beProductOf(), oofem::FloatMatrix::jaco_(), nu, oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by computeDamageTensor(), and givePlaneStressStiffMtrx().
void oofem::AnisotropicDamageMaterial::computePrincValDir2D | ( | double & | D1, |
double & | D2, | ||
double & | c, | ||
double & | s, | ||
double | Dx, | ||
double | Dy, | ||
double | Dxy | ||
) |
Definition at line 205 of file anisodamagemodel.C.
Referenced by computeDamage(), computeDimensionlessOutOfPlaneStress(), computeInplaneStress(), computeOutOfPlaneStrain(), and giveRealStressVector_PlaneStress().
|
protectedvirtual |
Definition at line 1753 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeTraceD(), E, oofem::FloatMatrix::jaco_(), nu, oofem::FloatMatrix::resize(), oofem::FloatMatrix::subtract(), and oofem::FloatMatrix::zero().
Referenced by give3dMaterialStiffnessMatrix(), and givePlaneStressStiffMtrx().
double oofem::AnisotropicDamageMaterial::computeTraceD | ( | double | equivStrain | ) |
Definition at line 328 of file anisodamagemodel.C.
References aA, damageLawType, DLT_Desmorat1, DLT_Desmorat2, DLT_Exponential, DLT_Linear, kappa0, kappaf, nu, and OOFEM_ERROR.
Referenced by computeDamage(), computeDamageTensor(), computeSecantOperator(), and giveRealStressVector().
double oofem::AnisotropicDamageMaterial::computeTraceD | ( | FloatMatrix | tempDamageTensor, |
FloatMatrix | strainTensor, | ||
GaussPoint * | gp | ||
) |
Definition at line 1083 of file anisodamagemodel.C.
References oofem::FloatMatrix::at(), oofem::AnisotropicDamageMaterialStatus::giveFlag(), oofem::Material::giveStatus(), and oofem::AnisotropicDamageMaterialStatus::setTempFlag().
void oofem::AnisotropicDamageMaterial::correctBigValues | ( | FloatMatrix & | matrix | ) |
Definition at line 1068 of file anisodamagemodel.C.
References oofem::FloatMatrix::at(), and oofem::FloatMatrix::giveNumberOfRows().
Referenced by computeDamageTensor(), giveRealStressVector(), obtainAlpha1(), obtainAlpha2(), and obtainAlpha3().
|
inlinevirtual |
Creates new copy of associated status and inserts it into given integration point.
gp | Integration point where newly created status will be stored. |
Reimplemented from oofem::Material.
Definition at line 283 of file anisodamagemodel.h.
References oofem::AnisotropicDamageMaterialStatus::AnisotropicDamageMaterialStatus(), and oofem::FEMComponent::domain.
|
protectedvirtual |
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.
answer | Stiffness matrix. |
mmode | Material response mode. |
gp | Integration point, which load history is used. |
tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 1263 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::IsotropicLinearElasticMaterial::give3dMaterialStiffnessMatrix(), giveLinearElasticMaterial(), oofem::Material::giveStatus(), and oofem::StructuralMaterialStatus::giveTempStrainVector().
|
virtual |
Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.
answer | Computed results. |
mode | Material response mode. |
gp | Integration point. |
tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 1156 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeSecantOperator(), giveLinearElasticMaterial(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::AnisotropicDamageMaterialStatus::giveTempDamage(), oofem::StructuralMaterialStatus::giveTempStrainVector(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
inlinevirtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 214 of file anisodamagemodel.h.
|
virtual |
Setups the input record string of receiver.
input | Dynamic input record to be filled by receiver. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 2082 of file anisodamagemodel.C.
References _IFT_AnisotropicDamageMaterial_kappa0, oofem::StructuralMaterial::giveInputRecord(), kappa0, and oofem::DynamicInputRecord::setField().
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 215 of file anisodamagemodel.h.
References _IFT_AnisotropicDamageMaterial_Name.
|
virtual |
Returns the integration point corresponding value in Reduced form.
answer | Contain corresponding ip value, zero sized if not available. |
gp | Integration point to which the value refers. |
type | Determines the type of internal variable. |
tStep | Determines the time step. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 1958 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::AnisotropicDamageMaterialStatus::giveDamage(), oofem::AnisotropicDamageMaterialStatus::giveDissWork(), oofem::StructuralMaterial::giveIPValue(), oofem::AnisotropicDamageMaterialStatus::giveKappa(), oofem::Material::giveStatus(), oofem::AnisotropicDamageMaterialStatus::giveStressWork(), oofem::AnisotropicDamageMaterialStatus::giveTempDamage(), oofem::FloatMatrix::jaco_(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
inline |
Returns reference to undamaged (bulk) material.
Definition at line 218 of file anisodamagemodel.h.
References oofem::AnisotropicDamageMaterialStatus::damage, oofem::IntegrationPointStatus::gp, oofem::AnisotropicDamageMaterialStatus::kappa, and oofem::AnisotropicDamageMaterialStatus::tempDamage.
Referenced by computePlaneStressSigmaZ(), give1dStressStiffMtrx(), give3dMaterialStiffnessMatrix(), givePlaneStressStiffMtrx(), and giveRealStressVector().
|
protectedvirtual |
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 , but relations for and are included).
answer | Stiffness matrix. |
mmode | Material response mode. |
gp | Integration point, which load history is used. |
tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 1259 of file anisodamagemodel.C.
|
protectedvirtual |
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.
answer | Stiffness matrix. |
mmode | Material response mode. |
gp | Integration point, which load history is used. |
tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 1197 of file anisodamagemodel.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computePlaneStressStrain(), computeSecantOperator(), E, giveLinearElasticMaterial(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::AnisotropicDamageMaterialStatus::giveTempDamage(), oofem::StructuralMaterialStatus::giveTempStrainVector(), oofem::StructuralMaterialStatus::giveTempStressVector(), nu, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
|
virtual |
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.
answer | Stress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress. |
gp | Integration point. |
reducedStrain | Strain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain. |
tStep | Current time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 508 of file anisodamagemodel.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beMatrixForm(), oofem::FloatArray::beProductOf(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beSymVectorForm(), computeCorrectionFactor(), computeDamageTensor(), computeEquivalentStrain(), computeKappa(), computeTraceD(), oofem::AnisotropicDamageMaterialStatus::computeWork(), correctBigValues(), oofem::AnisotropicDamageMaterialStatus::giveDamage(), oofem::StructuralMaterial::giveFullSymVectorForm(), giveLinearElasticMaterial(), oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveReducedSymVectorForm(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::Material::initTempStatus(), oofem::FloatMatrix::jaco_(), kappa0, oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::FloatMatrix::resize(), oofem::AnisotropicDamageMaterialStatus::setTempDamage(), oofem::AnisotropicDamageMaterialStatus::setTempKappa(), oofem::FloatMatrix::subtract(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
|
inlinevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 260 of file anisodamagemodel.h.
References oofem::FEMComponent::giveInputRecord(), and oofem::MaterialStatus::initializeFrom().
|
inlinevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 254 of file anisodamagemodel.h.
|
inlinevirtual |
Default implementation relies on giveRealStressVector_3d.
Reimplemented from oofem::StructuralMaterial.
Definition at line 256 of file anisodamagemodel.h.
|
virtual |
Plane-stress version of the stress evaluation algorithm.
Reimplemented from oofem::StructuralMaterial.
Definition at line 92 of file anisodamagemodel.C.
References AD_MAXITER, AD_TOLERANCE, oofem::FloatArray::at(), computeDamage(), computeDimensionlessOutOfPlaneStress(), computeInplaneStress(), computeOutOfPlaneStrain(), computePrincValDir2D(), oofem::AnisotropicDamageMaterialStatus::computeWork(), oofem::AnisotropicDamageMaterialStatus::giveDamage(), oofem::AnisotropicDamageMaterialStatus::giveKappa(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::Material::initTempStatus(), kappa0, oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::macbra(), OOFEM_ERROR, oofem::AnisotropicDamageMaterialStatus::setTempDamage(), oofem::AnisotropicDamageMaterialStatus::setTempKappa(), and oofem::AnisotropicDamageMaterialStatus::setTempStrainZ().
|
inlinevirtual |
Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero·
Reimplemented from oofem::StructuralMaterial.
Definition at line 258 of file anisodamagemodel.h.
|
virtual |
Tests if material supports material mode.
mode | Required material mode. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 77 of file anisodamagemodel.C.
|
inlinevirtual |
Returns nonzero if receiver is non linear.
Reimplemented from oofem::Material.
Definition at line 210 of file anisodamagemodel.h.
|
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.
ir | Input record to initialize from. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 2025 of file anisodamagemodel.C.
References _IFT_AnisotropicDamageMaterial_aA, _IFT_AnisotropicDamageMaterial_damageLawType, _IFT_AnisotropicDamageMaterial_kappa0, _IFT_AnisotropicDamageMaterial_kappaf, aA, damageLawType, DLT_Desmorat1, DLT_Desmorat2, DLT_Exponential, DLT_Linear, E, equivStrainType, EST_ElasticEnergy, EST_ElasticEnergyPositiveStrain, EST_ElasticEnergyPositiveStress, EST_Griffith, EST_Mazars, EST_Mises, EST_Rankine_Smooth, EST_Rankine_Standard, oofem::IsotropicLinearElasticMaterial::givePoissonsRatio(), oofem::IsotropicLinearElasticMaterial::giveYoungsModulus(), oofem::IsotropicLinearElasticMaterial::initializeFrom(), oofem::StructuralMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, kappa0, kappaf, linearElasticMaterial, and nu.
double oofem::AnisotropicDamageMaterial::obtainAlpha1 | ( | FloatMatrix | tempDamageTensor, |
double | deltaLambda, | ||
FloatMatrix | positiveStrainTensor, | ||
double | damageThreshold | ||
) |
Obtains the proportion of the damage tensor that is needed to get the first eigenvalue equal to the damage threshold.
Definition at line 826 of file anisodamagemodel.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::at(), oofem::FloatMatrix::beProductOf(), correctBigValues(), oofem::FloatArray::giveSize(), oofem::FloatMatrix::jaco_(), and oofem::FloatMatrix::times().
Referenced by computeDamageTensor().
double oofem::AnisotropicDamageMaterial::obtainAlpha2 | ( | FloatMatrix | tempDamageTensor, |
double | deltaLambda, | ||
FloatMatrix | positiveStrainTensor, | ||
FloatMatrix | ProjMatrix, | ||
double | damageThreshold | ||
) |
Obtains the proportion of the damage tensor that is needed to get the second eigenvalue equal to the damage threshold.
Definition at line 906 of file anisodamagemodel.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::at(), correctBigValues(), oofem::FloatArray::giveSize(), oofem::FloatMatrix::jaco_(), and oofem::FloatMatrix::times().
Referenced by computeDamageTensor().
double oofem::AnisotropicDamageMaterial::obtainAlpha3 | ( | FloatMatrix | tempDamageTensor, |
double | deltaLambda, | ||
FloatMatrix | positiveStrainTensor, | ||
FloatArray | vec3, | ||
double | damageThreshold | ||
) |
Obtains the proportion of the damage tensor that is needed to get the third eigenvalue equal to the damage threshold.
Definition at line 989 of file anisodamagemodel.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beDyadicProductOf(), oofem::FloatMatrix::beProductOf(), correctBigValues(), oofem::FloatMatrix::jaco_(), and oofem::FloatMatrix::times().
Referenced by computeDamageTensor().
|
protected |
Damage parameter a*A, needed to obtain Kappa(trD), according to eq. 33 in the paper mentioned above.
Definition at line 194 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeTraceD(), and initializeFrom().
|
protected |
Parameter specifying the damage law.
Definition at line 202 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeTraceD(), and initializeFrom().
|
protected |
Young's modulus.
Definition at line 186 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeInplaneStress(), computePlaneStressSigmaZ(), computeSecantOperator(), givePlaneStressStiffMtrx(), and initializeFrom().
|
protected |
Parameter specifying the definition of equivalent strain.
Definition at line 198 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeEquivalentStrain(), and initializeFrom().
|
protected |
Damage threshold kappa0, as defined in the paper mentioned above.
Definition at line 190 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeKappa(), computeTraceD(), giveInputRecord(), giveRealStressVector(), giveRealStressVector_PlaneStress(), and initializeFrom().
|
protected |
Damage parameter kappa_f (in the paper denoted as "a")
Definition at line 192 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeKappa(), computeTraceD(), and initializeFrom().
|
protected |
Reference to bulk (undamaged) material.
Definition at line 184 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), initializeFrom(), and ~AnisotropicDamageMaterial().
|
protected |
Poisson's ratio.
Definition at line 188 of file anisodamagemodel.h.
Referenced by AnisotropicDamageMaterial(), computeDimensionlessOutOfPlaneStress(), computeEquivalentStrain(), computeInplaneStress(), computeOutOfPlaneStrain(), computePlaneStressSigmaZ(), computePlaneStressStrain(), computeSecantOperator(), computeTraceD(), givePlaneStressStiffMtrx(), and initializeFrom().