OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
Material damage model for transversely orthotropic material. More...
#include <compodamagemat.h>
Public Member Functions | |
CompoDamageMat (int n, Domain *d) | |
Constructor. More... | |
virtual | ~CompoDamageMat () |
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... | |
virtual MaterialStatus * | CreateStatus (GaussPoint *gp) const |
Creates new copy of associated status and inserts it into given integration point. More... | |
virtual void | give3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mmode, 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 &, 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_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... | |
Public Member Functions inherited from oofem::StructuralMaterial | |
StructuralMaterial (int n, Domain *d) | |
Constructor. More... | |
virtual | ~StructuralMaterial () |
Destructor. More... | |
virtual int | hasMaterialModeCapability (MaterialMode mode) |
Tests if material supports material mode. 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_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_ShellStressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep) |
virtual void | giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_StressControl. More... | |
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 (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane stress stiffness matrix of receiver. More... | |
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 (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane strain stiffness matrix of receiver. More... | |
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 (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d stiffness matrix of receiver. More... | |
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 | hasNonLinearBehaviour () |
Returns nonzero if receiver is non linear. More... | |
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... | |
Public Attributes | |
int | afterIter |
Optional parameter determining after how many iterations within the time step the damage is calculated. More... | |
Protected Member Functions | |
void | giveUnrotated3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp) |
Returns 3D material stiffness matrix [6x6] in unrotated form. More... | |
int | giveMatStiffRotationMatrix (FloatMatrix &answer, GaussPoint *gp) |
Returns [6x6] rotation matrix in the global coordinate system. More... | |
void | giveCharLength (CompoDamageMatStatus *status, GaussPoint *gp, FloatMatrix &elementCs) |
Fills array elemCharLength with characteristic length related to three perpendicular planes. More... | |
void | giveCharLengthForModes (FloatArray &charLenModes, GaussPoint *gp) |
Computes characteristic length for fixed planes of material orientation. More... | |
void | checkSnapBack (GaussPoint *gp, MaterialMode mMode) |
Check that element is small enough or Gf is large enough to prevent the snap-back. More... | |
Protected Attributes | |
FloatArray | inputTension |
Six stress components of tension components read from the input file. More... | |
FloatArray | inputCompression |
Six stress components of compression components read from the input file. More... | |
IntArray | allowSnapBack |
Stress components which are allowed for snap back [6 tension, 6 compression]. 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... | |
Material damage model for transversely orthotropic material.
Fixed cracks are induced in principal material coordinates, always perpendicular to material axis. Six cracking modes are implemented in three orthogonal directions - three for tension/compression failure and three for shear failure. Material orientation can be specified on each element with "lmcs" keyword, otherwise global material orientation is assumed.
Evolution of cracks is determined separately in compression and tension. Linear softening is assumed (can be easily changed to exponential but requires then inner iterations)
For derivation of the model see the book Bazant and Planas, Fracture and Size Effect in Concrete and Other Quasibrittle Materials, pp.236 or article Bazant and Oh: Crack band theory for fracture of concrete, Materials and Structures, 1983.
The model is aimed for 3D problems but extension for 1D truss works. In this particular case, only the first array component is used in all involved variables.
Definition at line 139 of file compodamagemat.h.
oofem::CompoDamageMat::CompoDamageMat | ( | int | n, |
Domain * | d | ||
) |
Constructor.
Definition at line 48 of file compodamagemat.C.
|
virtual |
Destructor.
Definition at line 53 of file compodamagemat.C.
|
protected |
Check that element is small enough or Gf is large enough to prevent the snap-back.
gp | Integration point. |
mMode | Type of material (_1dMat, _3dMat supported). |
Definition at line 478 of file compodamagemat.C.
References oofem::__MaterialModeToString(), allowSnapBack, oofem::FloatArray::at(), oofem::IntArray::contains(), oofem::CompoDamageMatStatus::elemCharLength, Ex, Ey, Ez, oofem::Material::give(), giveCharLengthForModes(), oofem::GaussPoint::giveElement(), oofem::FEMComponent::giveNumber(), oofem::GaussPoint::giveNumber(), oofem::Material::giveStatus(), Gxy, Gxz, Gyz, inputCompression, inputTension, OOFEM_ERROR, and OOFEM_LOG_INFO.
Referenced by giveRealStressVector().
|
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 153 of file compodamagemat.h.
References oofem::CompoDamageMatStatus::CompoDamageMatStatus(), oofem::FEMComponent::domain, and oofem::IntegrationPointStatus::gp.
|
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 125 of file compodamagemat.C.
References giveMatStiffRotationMatrix(), giveUnrotated3dMaterialStiffnessMatrix(), and oofem::FloatMatrix::rotatedWith().
|
protected |
Fills array elemCharLength with characteristic length related to three perpendicular planes.
The planes are of the same orientation as material.
status | Pointer to integration point's status. |
gp | Integration point. |
elementCs | Material orientation matrix. |
Definition at line 445 of file compodamagemat.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::CompoDamageMatStatus::elemCharLength, oofem::Element::giveCharacteristicLength(), and oofem::GaussPoint::giveElement().
Referenced by giveRealStressVector().
|
protected |
Computes characteristic length for fixed planes of material orientation.
charLenModes | Returns six lengths. |
gp | Integration point. |
Definition at line 464 of file compodamagemat.C.
References oofem::FloatArray::at(), oofem::CompoDamageMatStatus::elemCharLength, oofem::Material::giveStatus(), and oofem::FloatArray::resize().
Referenced by checkSnapBack(), and giveRealStressVector().
|
inlinevirtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 147 of file compodamagemat.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 117 of file compodamagemat.C.
References oofem::StructuralMaterial::giveInputRecord(), and OOFEM_ERROR.
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 148 of file compodamagemat.h.
References _IFT_CompoDamageMat_Name, oofem::FEMComponent::giveInputRecord(), and oofem::MaterialStatus::initializeFrom().
|
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 352 of file compodamagemat.C.
References oofem::StructuralMaterial::giveIPValue(), oofem::Material::giveStatus(), and oofem::CompoDamageMatStatus::omega.
|
protected |
Returns [6x6] rotation matrix in the global coordinate system.
The matrix relates local c.s. to global c.s. Local c.s. can be specified with 'mcs' flag defined on element.
answer | Rotation matrix [3x3]. |
gp | Integration point. |
Definition at line 415 of file compodamagemat.C.
References oofem::__MaterialModeToString(), oofem::GaussPoint::giveElement(), oofem::Element::giveLocalCoordinateSystem(), oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveStrainVectorTranformationMtrx(), and OOFEM_ERROR.
Referenced by give3dMaterialStiffnessMatrix().
|
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 137 of file compodamagemat.C.
References oofem::__MaterialModeToString(), oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::beUnitMatrix(), checkSnapBack(), oofem::CompoDamageMatStatus::elemCharLength, Ex, Ey, Ez, oofem::Material::give(), oofem::Element::giveCharacteristicLength(), giveCharLength(), giveCharLengthForModes(), oofem::GaussPoint::giveElement(), oofem::Element::giveLocalCoordinateSystem(), oofem::GaussPoint::giveMaterialMode(), oofem::FEMComponent::giveNumber(), oofem::Material::giveStatus(), oofem::StructuralMaterialStatus::giveStrainVector(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::StructuralMaterialStatus::giveStressVector(), giveUnrotated3dMaterialStiffnessMatrix(), Gxy, Gxz, Gyz, oofem::CompoDamageMatStatus::hasSnapBack, oofem::CompoDamageMatStatus::initDamageStress, inputCompression, inputTension, oofem::CompoDamageMatStatus::Iteration, oofem::CompoDamageMatStatus::kappa, oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::max(), oofem::CompoDamageMatStatus::maxStrainAtZeroStress, oofem::min(), NYxy, NYxz, NYyx, NYyz, NYzx, NYzy, oofem::CompoDamageMatStatus::omega, OOFEM_ERROR, OOFEM_WARNING, oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), oofem::CompoDamageMatStatus::strainAtMaxStress, oofem::CompoDamageMatStatus::tempKappa, oofem::CompoDamageMatStatus::tempOmega, oofem::CompoDamageMatStatus::tempStressMLCS, oofem::StructuralMaterial::transformStrainVectorTo(), and oofem::StructuralMaterial::transformStressVectorTo().
|
inlinevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 165 of file compodamagemat.h.
|
inlinevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 163 of file compodamagemat.h.
|
protected |
Returns 3D material stiffness matrix [6x6] in unrotated form.
The matrix is reduced by by omega variables.
answer | Full symmetric matrix. |
mode | Material mode of stiffness matrix (elastic, secant). |
gp | Integration point. |
Definition at line 364 of file compodamagemat.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), Ex, Ey, Ez, oofem::Material::give(), oofem::Material::giveStatus(), Gxy, Gxz, Gyz, NYxy, NYyz, oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), oofem::CompoDamageMatStatus::tempOmega, and oofem::FloatMatrix::zero().
Referenced by give3dMaterialStiffnessMatrix(), and giveRealStressVector().
|
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 58 of file compodamagemat.C.
References _IFT_CompoDamageMat_afteriter, _IFT_CompoDamageMat_allowSnapBack, _IFT_CompoDamageMat_compres_f0_gf, _IFT_CompoDamageMat_exx, _IFT_CompoDamageMat_eyyezz, _IFT_CompoDamageMat_Gxy, _IFT_CompoDamageMat_nuxynuxz, _IFT_CompoDamageMat_nuyz, _IFT_CompoDamageMat_tension_f0_gf, oofem::Dictionary::add(), afterIter, allowSnapBack, oofem::FloatArray::at(), Ex, Ey, Ez, oofem::Material::give(), oofem::FloatArray::giveSize(), Gxy, Gxz, Gyz, oofem::Material::initializeFrom(), inputCompression, inputTension, IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, NYxy, NYxz, NYyx, NYyz, NYzx, NYzy, OOFEM_ERROR, and oofem::Material::propertyDictionary.
int oofem::CompoDamageMat::afterIter |
Optional parameter determining after how many iterations within the time step the damage is calculated.
This is important for stress evaluation which is unequilibrated in the beginning. Variables strainAtMaxStress, initDamageStress, maxStrainAtZeroStress are evaluated afterIter.
Definition at line 175 of file compodamagemat.h.
Referenced by initializeFrom().
|
protected |
Stress components which are allowed for snap back [6 tension, 6 compression].
Definition at line 201 of file compodamagemat.h.
Referenced by checkSnapBack(), and initializeFrom().
|
protected |
Six stress components of compression components read from the input file.
Definition at line 198 of file compodamagemat.h.
Referenced by checkSnapBack(), giveRealStressVector(), and initializeFrom().
|
protected |
Six stress components of tension components read from the input file.
Definition at line 196 of file compodamagemat.h.
Referenced by checkSnapBack(), giveRealStressVector(), and initializeFrom().