OOFEM 3.0
Loading...
Searching...
No Matches
oofem::IsotropicLinearElasticMaterial Class Reference

#include <isolinearelasticmaterial.h>

Inheritance diagram for oofem::IsotropicLinearElasticMaterial:
Collaboration diagram for oofem::IsotropicLinearElasticMaterial:

Public Member Functions

 IsotropicLinearElasticMaterial (int n, Domain *d)
 IsotropicLinearElasticMaterial (int n, Domain *d, double E, double nu)
const char * giveClassName () const override
const char * giveInputRecordName () const override
void initializeFrom (InputRecord &ir) override
void giveInputRecord (DynamicInputRecord &input) override
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
void initTangents ()
 Initialized fixed size tangents. Called by ctor and initializeFrom.
double give (int aProperty, GaussPoint *gp) const override
double giveYoungsModulus () const
 Returns Young's modulus.
double givePoissonsRatio () const
 Returns Poisson's ratio.
double giveShearModulus () const override
 Returns the shear elastic modulus \( G = \frac{E}{2(1+\nu)} \).
double giveBulkModulus () const
 Returns the bulk elastic modulus \( K = \frac{E}{3(1-2\nu)} \).
FloatMatrixF< 3, 3 > givePlaneStressStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 4, 4 > givePlaneStrainStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
FloatMatrixF< 1, 1 > give1dStressStiffMtrx (MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
double giveQcElasticParamneter () override
double giveQcPlasticParamneter () override
InterfacegiveInterface (InterfaceType t) override
void giveCharacteristicMatrix (FloatMatrix &answer, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override
void giveCharacteristicVector (FloatArray &answer, FloatArray &flux, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override
 Returns characteristic vector of the receiver.
void giveDeviatoric3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
void giveDeviatoricPlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) const override
void giveInverseOfBulkModulus (double &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) override
void giveRealStressVectorUP_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, double pressure, TimeStep *tStep) const override
void giveRealStressVectorUP_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, double pressure, TimeStep *tStep) const override
double giveCharacteristicValue (MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override
 Returns characteristic value of the receiver.
Public Member Functions inherited from oofem::LinearElasticMaterial
 LinearElasticMaterial (int n, Domain *d)
 Constructor.
void initializeFrom (InputRecord &ir) override
void giveInputRecord (DynamicInputRecord &input) override
const FloatMatrixF< 6, 6 > & giveTangent () const
const FloatArrayF< 6 > & giveAlpha () const
void computesSubTangents ()
FloatMatrixF< 6, 6 > give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 6 > giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const override
FloatArrayF< 6 > giveRealStressVector_3d (const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
virtual FloatArrayF< 6 > giveRealStressVector_3dDegeneratedShell (const FloatArrayF< 6 > &reducedF, GaussPoint *gp, TimeStep *tStep) const
FloatArrayF< 3 > giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector_StressControl.
FloatArrayF< 1 > giveRealStressVector_1d (const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector_StressControl.
FloatArrayF< 2 > giveRealStressVector_Warping (const FloatArrayF< 2 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector_StressControl.
FloatArrayF< 2 > giveRealStressVector_2dBeamLayer (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector_StressControl.
FloatArrayF< 5 > giveRealStressVector_PlateLayer (const FloatArrayF< 5 > &reducedE, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector_StressControl.
FloatArrayF< 3 > giveRealStressVector_Fiber (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const override
 Default implementation relies on giveRealStressVector_StressControl.
void giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) override
double giveEnergyDensity (GaussPoint *gp, TimeStep *tStep)
int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
std::unique_ptr< MaterialStatusCreateStatus (GaussPoint *gp) const override
bool hasCastingTimeSupport () const override
const char * giveClassName () const override
Public Member Functions inherited from oofem::StructuralMaterial
 StructuralMaterial (int n, Domain *d)
bool hasMaterialModeCapability (MaterialMode mode) const override
const char * giveClassName () const override
void initializeFrom (InputRecord &ir) override
void giveInputRecord (DynamicInputRecord &input) override
void giveCharacteristicMatrix (FloatMatrix &answer, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override
 Returns characteristic matrix of the receiver.
void giveCharacteristicVector (FloatArray &answer, FloatArray &flux, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override
 Returns characteristic vector of the receiver.
virtual void giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
virtual void giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) const
virtual FloatArrayF< 4 > giveRealStressVector_PlaneStrain (const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const
 Default implementation relies on giveRealStressVector_3d.
virtual FloatArray giveRealStressVector_StressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const
 Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero·
virtual FloatArray giveRealStressVector_ShellStressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const
virtual FloatArrayF< 3 > giveRealStressVector_2dPlateSubSoil (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const
 Default implementation is not provided.
virtual FloatArrayF< 6 > giveRealStressVector_3dBeamSubSoil (const FloatArrayF< 6 > &reducedE, GaussPoint *gp, TimeStep *tStep) const
virtual FloatArrayF< 9 > giveFirstPKStressVector_3d (const FloatArrayF< 9 > &vF, GaussPoint *gp, TimeStep *tStep) const
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
virtual FloatArrayF< 5 > giveFirstPKStressVector_PlaneStrain (const FloatArrayF< 5 > &vF, GaussPoint *gp, TimeStep *tStep) const
 Default implementation relies on giveFirstPKStressVector_3d.
virtual FloatArray giveFirstPKStressVector_StressControl (const FloatArray &reducedvF, const IntArray &FControl, GaussPoint *gp, TimeStep *tStep) const
 Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero·
virtual FloatArrayF< 4 > giveFirstPKStressVector_PlaneStress (const FloatArrayF< 4 > &vF, GaussPoint *gp, TimeStep *tStep) const
 Default implementation relies on giveFirstPKStressVector_StressControl.
virtual FloatArrayF< 1 > giveFirstPKStressVector_1d (const FloatArrayF< 1 > &vF, GaussPoint *gp, TimeStep *tStep) const
 Default implementation relies on giveFirstPKStressVector_StressControl.
virtual void giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual void giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual void giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
virtual void giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
double giveReferenceTemperature ()
virtual FloatArray computeStressIndependentStrainVector (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const
FloatArrayF< 6 > computeStressIndependentStrainVector_3d (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const
virtual void giveStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
virtual FloatMatrixF< 9, 9 > give3dMaterialStiffnessMatrix_dPdF (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const
virtual void give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
void giveStressDependentPartOfStrainVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode) const
int setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) override
int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override
virtual FloatMatrixF< 4, 4 > givePlaneStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual void givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual FloatMatrixF< 5, 5 > givePlaneStrainStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual void givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual FloatMatrixF< 1, 1 > give1dStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual void give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
virtual FloatMatrixF< 2, 2 > give2dBeamLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 5, 5 > givePlateLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 3, 3 > giveFiberStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 3, 3 > give2dPlateSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
virtual FloatMatrixF< 6, 6 > give3dBeamSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const
Public Member Functions inherited from oofem::Material
 Material (int n, Domain *d)
virtual ~Material ()=default
 Destructor.
virtual bool isCharacteristicMtrxSymmetric (MatResponseMode rMode) const
virtual bool hasProperty (int aProperty, GaussPoint *gp) const
virtual void modifyProperty (int aProperty, double value, GaussPoint *gp)
double giveCastingTime () const
virtual bool isActivated (TimeStep *tStep) const
void printYourself () override
 Prints receiver state on stdout. Useful for debugging.
virtual void saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
virtual void restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
int checkConsistency () override
virtual void restoreConsistency (GaussPoint *gp)
virtual int initMaterial (Element *element)
virtual MaterialStatusgiveStatus (GaussPoint *gp) const
virtual int packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
virtual int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)
virtual int estimatePackSize (DataStream &buff, GaussPoint *ip)
virtual double predictRelativeComputationalCost (GaussPoint *gp)
virtual double predictRelativeRedistributionCost (GaussPoint *gp)
virtual void initTempStatus (GaussPoint *gp) const
void saveContext (DataStream &stream, ContextMode mode) override
void restoreContext (DataStream &stream, ContextMode mode) override
Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
virtual ~FEMComponent ()=default
 Virtual destructor.
DomaingiveDomain () const
virtual void setDomain (Domain *d)
int giveNumber () const
void setNumber (int num)
virtual void updateLocalNumbering (EntityRenumberingFunctor &f)
virtual void initializeFrom (InputRecord &ir, int priority)
virtual void initializeFinish ()
virtual void postInitialize ()
 Performs post initialization steps. Called after all components are created and initialized.
virtual void printOutputAt (FILE *file, TimeStep *tStep)
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros).
Public Member Functions inherited from oofem::QCMaterialExtensionInterface
 QCMaterialExtensionInterface ()
virtual ~QCMaterialExtensionInterface ()
 Destructor.
Public Member Functions inherited from oofem::Interface
 Interface ()
 Constructor.
virtual ~Interface ()
Public Member Functions inherited from oofem::MixedPressureMaterialExtensionInterface
 MixedPressureMaterialExtensionInterface (Domain *d)
virtual ~MixedPressureMaterialExtensionInterface ()
 Destructor.
virtual void giveDeviatoricConstitutiveMatrix (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) const
void giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, double pressure, TimeStep *tStep) const
virtual void giveFiniteStrainGeneralizedStressVectors (FloatArray &sigma, GaussPoint *gp, const FloatArray &devF, double pressure, TimeStep *tStep)

Static Public Member Functions

static double computeBulkModulusFromYoungAndPoisson (double young, double nu)
static double computeShearModulusFromYoungAndPoisson (double young, double nu)
Static Public Member Functions inherited from oofem::StructuralMaterial
static int giveSymVI (int ind1, int ind2)
static int giveVI (int ind1, int ind2)
static FloatMatrixF< 9, 9 > convert_dSdE_2_dPdF_3D (const FloatMatrixF< 6, 6 > &dSdE, const FloatArrayF< 6 > &S, const FloatArrayF< 9 > &F)
static FloatMatrixF< 5, 5 > convert_dSdE_2_dPdF_PlaneStrain (const FloatMatrixF< 4, 4 > &dSdE, const FloatArrayF< 4 > &S, const FloatArrayF< 5 > &F)
static FloatMatrixF< 4, 4 > convert_dSdE_2_dPdF_PlaneStress (const FloatMatrixF< 3, 3 > &dSdE, const FloatArrayF< 3 > &S, const FloatArrayF< 4 > &F)
static FloatMatrixF< 1, 1 > convert_dSdE_2_dPdF_1D (const FloatMatrixF< 1, 1 > &dSdE, const FloatArrayF< 1 > &S, const FloatArrayF< 1 > &F)
static void computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode)
 Common functions for convenience.
static FloatArrayF< 3 > computePrincipalValues (const FloatMatrixF< 3, 3 > &s)
static FloatArrayF< 3 > computePrincipalValues (double I1, double I2, double I3)
static void computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode)
static std::pair< FloatArrayF< 3 >, FloatMatrixF< 3, 3 > > computePrincipalValDir (const FloatMatrixF< 3, 3 > &s)
static FloatArrayF< 6 > computeDeviator (const FloatArrayF< 6 > &s)
static std::pair< FloatArrayF< 6 >, double > computeDeviatoricVolumetricSplit (const FloatArrayF< 6 > &s)
static FloatArrayF< 6 > computeDeviatoricVolumetricSum (const FloatArrayF< 6 > &dev, double mean)
static FloatArrayF< 6 > applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu)
static FloatArrayF< 6 > applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double GModulus)
static FloatArrayF< 6 > applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu)
static FloatArrayF< 6 > applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double GModulus)
static FloatArrayF< 6 > applyElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu)
static FloatArrayF< 6 > applyElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu)
static double computeStressNorm (const FloatArrayF< 6 > &stress)
static double computeFirstInvariant (const FloatArrayF< 6 > &s)
static double computeSecondStressInvariant (const FloatArrayF< 6 > &s)
static double computeThirdStressInvariant (const FloatArrayF< 6 > &s)
static double computeFirstCoordinate (const FloatArrayF< 6 > &s)
static double computeSecondCoordinate (const FloatArrayF< 6 > &s)
static double computeThirdCoordinate (const FloatArrayF< 6 > &s)
static int giveVoigtVectorMask (IntArray &answer, MaterialMode mmode)
static int giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode)
static void giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode)
static int giveSizeOfVoigtVector (MaterialMode mmode)
static int giveSizeOfVoigtSymVector (MaterialMode mmode)
static void giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
 Converts the reduced symmetric Voigt vector (2nd order tensor) to full form.
static void giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode)
 Converts the reduced deformation gradient Voigt vector (2nd order tensor).
static void giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form.
static void giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the full symmetric Voigt vector (2nd order tensor) to reduced form.
static void giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode)
 Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form.
static void giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode)
 Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form.
static void giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
 Converts the full symmetric Voigt matrix (4th order tensor) to reduced form.
static void giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode)
 Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form.
static FloatArrayF< 6 > transformStrainVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &strain, bool transpose=false)
static FloatArrayF< 6 > transformStressVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &stress, bool transpose=false)
static double computeVonMisesStress (const FloatArray &currentStress)
static double computeVonMisesStress_3D (const FloatArrayF< 6 > &stress)
static double computeVonMisesStress_PlaneStress (const FloatArrayF< 3 > &stress)
static FloatMatrixF< 6, 6 > giveStrainVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false)
static FloatMatrixF< 3, 3 > give2DStrainVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false)
static FloatMatrixF< 6, 6 > giveStressVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false)
static FloatMatrixF< 3, 3 > givePlaneStressVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false)
static void sortPrincDirAndValCloseTo (FloatArray &pVal, FloatMatrix &pDir, const FloatMatrix &toPDir)

Protected Attributes

double E = 0
 Young's modulus.
double nu = 0
 Poisson's ratio.
double G = 0
 Shear modulus.
double a = 0
 Alpha.
Protected Attributes inherited from oofem::LinearElasticMaterial
double preCastStiffnessReduction = 0.
 artificial isotropic damage to reflect reduction in stiffness for time < castingTime.
FloatMatrixF< 6, 6 > tangent
 Preconstructed 3d tangent.
FloatMatrixF< 4, 4 > tangentPlaneStrain
FloatMatrixF< 3, 3 > tangentPlaneStress
FloatArrayF< 6 > alpha
 Thermal expansion.
Protected Attributes inherited from oofem::StructuralMaterial
double referenceTemperature = 0.
 Reference temperature (temperature, when material has been built into structure).
MatResponseMode SCStiffMode = TangentStiffness
 stifness mode used in stress control
double SCRelTol = 1.e-3
 relative tolerance for stress control
double SCAbsTol = 1.e-12
 absolute stress tolerance for stress control
int SCMaxiter = 100000
 maximum iterations for stress-control
Protected Attributes inherited from oofem::Material
Dictionary propertyDictionary
double castingTime
int preCastingTimeMat
 Material existing before casting time - optional parameter, zero by default.
Protected Attributes inherited from oofem::FEMComponent
int number
 Component number.
Domaindomain
 Link to domain object, useful for communicating with other FEM components.
Protected Attributes inherited from oofem::MixedPressureMaterialExtensionInterface
Domaindom = nullptr

Additional Inherited Members

Static Public Attributes inherited from oofem::StructuralMaterial
static std::array< std::array< int, 3 >, 3 > vIindex
 Voigt index map.
static std::array< std::array< int, 3 >, 3 > svIndex
 Symmetric Voigt index map.

Detailed Description

This class implements an isotropic linear elastic material in a finite element problem. For large deformation analysis it becomes the St. Venant-Kirchoff hyperelasticity model.

Tasks:

  • Returning standard material stiffness matrix for 3d-case. according to current state determined by using data stored in Gausspoint.
  • methods give2dPlaneStressMtrx, givePlaneStrainMtrx, give1dStressMtrx are introduced since form of this matrices is well known, and for faster response mainly in linear elastic problems.
  • Returning a material property (method 'give'). Only for non-standard elements.
  • Returning real stress state vector(tensor) at gauss point for 3d - case.

Definition at line 74 of file isolinearelasticmaterial.h.

Constructor & Destructor Documentation

◆ IsotropicLinearElasticMaterial() [1/2]

oofem::IsotropicLinearElasticMaterial::IsotropicLinearElasticMaterial ( int n,
Domain * d )

Creates a new IsotropicLinearElasticMaterial class instance with given number belonging to domain d.

Parameters
nmaterial model number in domain
ddomain which receiver belongs to

Definition at line 48 of file isolinearelasticmaterial.C.

References oofem::LinearElasticMaterial::LinearElasticMaterial(), and oofem::MixedPressureMaterialExtensionInterface::MixedPressureMaterialExtensionInterface().

◆ IsotropicLinearElasticMaterial() [2/2]

oofem::IsotropicLinearElasticMaterial::IsotropicLinearElasticMaterial ( int n,
Domain * d,
double E,
double nu )

Creates a new IsotropicLinearElasticMaterial class instance with given number belonging to domain d.

Parameters
nMaterial model number in domain.
dDomain which receiver belongs to.
EYoung modulus.
nuPoisson ratio.

Definition at line 52 of file isolinearelasticmaterial.C.

References a, E, G, initTangents(), oofem::LinearElasticMaterial::LinearElasticMaterial(), oofem::MixedPressureMaterialExtensionInterface::MixedPressureMaterialExtensionInterface(), and nu.

Member Function Documentation

◆ computeBulkModulusFromYoungAndPoisson()

double oofem::IsotropicLinearElasticMaterial::computeBulkModulusFromYoungAndPoisson ( double young,
double nu )
inlinestatic

Computes bulk modulus from given Young's modulus and Poisson's ratio.

Parameters
youngYoung's modulus ( \( E \)).
nuPoisson's ratio ( \( \nu \)).
Returns
Bulk modulus ( \( K = E/(3*(1-2*nu) \)).

Definition at line 150 of file isolinearelasticmaterial.h.

References nu.

◆ computeShearModulusFromYoungAndPoisson()

double oofem::IsotropicLinearElasticMaterial::computeShearModulusFromYoungAndPoisson ( double young,
double nu )
inlinestatic

Computes shear modulus from given Young's modulus and Poisson's ratio.

Parameters
youngYoung's modulus ( \( E \)).
nuPoisson's ratio ( \( \nu \)).
Returns
Shear modulus ( \( G = \frac{E}{2 (1+\nu)} \)).

Definition at line 161 of file isolinearelasticmaterial.h.

References nu.

◆ give()

double oofem::IsotropicLinearElasticMaterial::give ( int aProperty,
GaussPoint * gp ) const
overridevirtual

Returns the value of material property 'aProperty'. Property must be identified by unique int id. Integration point also passed to allow for materials with spatially varying properties

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

Reimplemented from oofem::Material.

Definition at line 142 of file isolinearelasticmaterial.C.

References E, Ex, Ey, Ez, G, Gxy, Gxz, Gyz, nu, NYxy, NYxz, NYyx, NYyz, NYzx, and NYzy.

◆ give1dStressStiffMtrx()

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

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

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

Reimplemented from oofem::StructuralMaterial.

Definition at line 194 of file isolinearelasticmaterial.C.

References E, oofem::TimeStep::giveIntrinsicTime(), and oofem::LinearElasticMaterial::preCastStiffnessReduction.

◆ giveBulkModulus()

double oofem::IsotropicLinearElasticMaterial::giveBulkModulus ( ) const
inline

Returns the bulk elastic modulus \( K = \frac{E}{3(1-2\nu)} \).

Definition at line 133 of file isolinearelasticmaterial.h.

References E, and nu.

Referenced by giveCharacteristicValue().

◆ giveCharacteristicMatrix()

void oofem::IsotropicLinearElasticMaterial::giveCharacteristicMatrix ( FloatMatrix & answer,
MatResponseMode type,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

◆ giveCharacteristicValue()

double oofem::IsotropicLinearElasticMaterial::giveCharacteristicValue ( MatResponseMode type,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

Returns characteristic value of the receiver.

Reimplemented from oofem::Material.

Definition at line 305 of file isolinearelasticmaterial.C.

References giveBulkModulus(), and oofem::Material::giveCharacteristicValue().

◆ giveCharacteristicVector()

void oofem::IsotropicLinearElasticMaterial::giveCharacteristicVector ( FloatArray & answer,
FloatArray & flux,
MatResponseMode type,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

◆ giveClassName()

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

Implements oofem::FEMComponent.

Definition at line 104 of file isolinearelasticmaterial.h.

◆ giveDeviatoric3dMaterialStiffnessMatrix()

void oofem::IsotropicLinearElasticMaterial::giveDeviatoric3dMaterialStiffnessMatrix ( FloatMatrix & answer,
MatResponseMode mode,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

Reimplemented from oofem::MixedPressureMaterialExtensionInterface.

Definition at line 213 of file isolinearelasticmaterial.C.

References E, and nu.

Referenced by giveRealStressVectorUP_3d().

◆ giveDeviatoricPlaneStrainStiffMtrx()

void oofem::IsotropicLinearElasticMaterial::giveDeviatoricPlaneStrainStiffMtrx ( FloatMatrix & answer,
MatResponseMode mode,
GaussPoint * gp,
TimeStep * tStep ) const
overridevirtual

◆ giveInputRecord()

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

Setups the input record string of receiver.

Parameters
inputDynamic input record to be filled by receiver.

Reimplemented from oofem::FEMComponent.

Definition at line 91 of file isolinearelasticmaterial.C.

References _IFT_IsotropicLinearElasticMaterial_e, _IFT_IsotropicLinearElasticMaterial_n, _IFT_IsotropicLinearElasticMaterial_talpha, a, E, nu, and oofem::DynamicInputRecord::setField().

◆ giveInputRecordName()

const char * oofem::IsotropicLinearElasticMaterial::giveInputRecordName ( ) const
inlineoverridevirtual
Returns
Input record name of the receiver.

Implements oofem::FEMComponent.

Definition at line 105 of file isolinearelasticmaterial.h.

References _IFT_IsotropicLinearElasticMaterial_Name.

◆ giveInterface()

Interface * oofem::IsotropicLinearElasticMaterial::giveInterface ( InterfaceType t)
inlineoverridevirtual

◆ giveInverseOfBulkModulus()

void oofem::IsotropicLinearElasticMaterial::giveInverseOfBulkModulus ( double & answer,
MatResponseMode mode,
GaussPoint * gp,
TimeStep * tStep )
inlineoverridevirtual

Implements oofem::MixedPressureMaterialExtensionInterface.

Definition at line 200 of file isolinearelasticmaterial.h.

References E, and nu.

◆ givePlaneStrainStiffMtrx()

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

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

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

Reimplemented from oofem::StructuralMaterial.

Definition at line 181 of file isolinearelasticmaterial.C.

References oofem::TimeStep::giveIntrinsicTime(), oofem::LinearElasticMaterial::preCastStiffnessReduction, and oofem::LinearElasticMaterial::tangentPlaneStrain.

◆ givePlaneStressStiffMtrx()

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

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

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

Reimplemented from oofem::StructuralMaterial.

Definition at line 168 of file isolinearelasticmaterial.C.

References oofem::TimeStep::giveIntrinsicTime(), oofem::LinearElasticMaterial::preCastStiffnessReduction, and oofem::LinearElasticMaterial::tangentPlaneStress.

◆ givePoissonsRatio()

double oofem::IsotropicLinearElasticMaterial::givePoissonsRatio ( ) const
inline

Returns Poisson's ratio.

Definition at line 127 of file isolinearelasticmaterial.h.

References nu.

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

◆ giveQcElasticParamneter()

double oofem::IsotropicLinearElasticMaterial::giveQcElasticParamneter ( )
inlineoverridevirtual

Implements oofem::QCMaterialExtensionInterface.

Definition at line 166 of file isolinearelasticmaterial.h.

References E.

◆ giveQcPlasticParamneter()

double oofem::IsotropicLinearElasticMaterial::giveQcPlasticParamneter ( )
inlineoverridevirtual

Implements oofem::QCMaterialExtensionInterface.

Definition at line 167 of file isolinearelasticmaterial.h.

◆ giveRealStressVectorUP_3d()

◆ giveRealStressVectorUP_PlaneStrain()

◆ giveShearModulus()

double oofem::IsotropicLinearElasticMaterial::giveShearModulus ( ) const
inlineoverridevirtual

Returns the shear elastic modulus \( G = \frac{E}{2(1+\nu)} \).

Reimplemented from oofem::LinearElasticMaterial.

Definition at line 130 of file isolinearelasticmaterial.h.

References G.

◆ giveYoungsModulus()

double oofem::IsotropicLinearElasticMaterial::giveYoungsModulus ( ) const
inline

Returns Young's modulus.

Definition at line 124 of file isolinearelasticmaterial.h.

References E.

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

◆ initializeFrom()

void oofem::IsotropicLinearElasticMaterial::initializeFrom ( InputRecord & ir)
overridevirtual

Initializes receiver according to object description stored in input record. The E modulus (keyword "E"), Poisson ratio ("nu") and coefficient of thermal dilatation alpha ("talpha") are read. The parent class instanciateFrom method is called.

Reimplemented from oofem::FEMComponent.

Definition at line 65 of file isolinearelasticmaterial.C.

References _IFT_IsotropicLinearElasticMaterial_e, _IFT_IsotropicLinearElasticMaterial_n, _IFT_IsotropicLinearElasticMaterial_talpha, a, E, G, initTangents(), IR_GIVE_FIELD, and nu.

◆ initTangents()

void oofem::IsotropicLinearElasticMaterial::initTangents ( )

◆ restoreContext()

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

Restores the receiver state previously written in stream.

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

Reimplemented from oofem::FEMComponent.

Definition at line 120 of file isolinearelasticmaterial.C.

References a, oofem::CIO_IOERR, CM_Definition, E, G, initTangents(), nu, oofem::DataStream::read(), and THROW_CIOERR.

◆ saveContext()

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

Stores receiver state to output stream.

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

Reimplemented from oofem::FEMComponent.

Definition at line 102 of file isolinearelasticmaterial.C.

References a, oofem::CIO_IOERR, CM_Definition, E, nu, THROW_CIOERR, and oofem::DataStream::write().

Member Data Documentation

◆ a

double oofem::IsotropicLinearElasticMaterial::a = 0
protected

◆ E

◆ G

double oofem::IsotropicLinearElasticMaterial::G = 0
protected

◆ nu


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

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