OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class implements a (local) nonassociated plasticity model based on the Drucker-Prager yield criterion with hardening and softening. More...
#include <druckerPragerPlasticitySM.h>
Public Types | |
enum | state_flag_values { DP_Elastic, DP_Unloading, DP_Yielding, DP_Vertex } |
Values of history variable state_flag. More... | |
Public Member Functions | |
DruckerPragerPlasticitySM (int n, Domain *d) | |
Constructor. More... | |
virtual | ~DruckerPragerPlasticitySM () |
Destructor. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
virtual void | giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &strainVector, TimeStep *tStep) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. 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... | |
void | performLocalStressReturn (GaussPoint *gp, const FloatArray &strain) |
Perform a standard local stress return using the function computeYieldValue at the specified Gauss point. More... | |
bool | checkForVertexCase (double eM, double gM, double kM, double trialStressJTwo, double volumetricStress, double tempKappa) |
Check if the trial stress state falls within the vertex region. More... | |
void | performRegularReturn (double eM, double gM, double kM, double trialStressJTwo, FloatArray &stressDeviator, double &volumetricStress, double &tempKappa) |
Perform stress return for regular case, i.e. More... | |
void | performVertexReturn (double eM, double gM, double kM, double trialStressJTwo, FloatArray &stressDeviator, double &volumetricStress, double &tempKappa, double volumetricElasticTrialStrain, double kappa) |
Perform stress return for vertex case, i.e. More... | |
double | computeYieldValue (double meanStress, double JTwo, double kappa, double eM) const |
Compute the yield value based on stress and hardening variable. More... | |
virtual double | computeYieldStressInShear (double kappa, double eM) const |
Compute the current yield stress in pure shear of the Drucker-Prager model according to the used hardening law. More... | |
virtual double | computeYieldStressPrime (double kappa, double eM) const |
Compute derivative of yield stress with respect to the hardening variable kappa. More... | |
void | giveRegAlgorithmicStiffMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
Compute and give back algorithmic stiffness matrix for the regular case (no vertex). More... | |
void | giveVertexAlgorithmicStiffMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
Compute consistent stiffness matrix for the vertex case. More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in Reduced form. More... | |
virtual bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More... | |
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... | |
virtual MaterialStatus * | CreateStatus (GaussPoint *gp) const |
Creates new copy of associated status and inserts it into given integration point. 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... | |
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 | giveInputRecord (DynamicInputRecord &input) |
Setups the input record string of receiver. 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 (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_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_1d (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) |
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 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 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 Attributes | |
int | hardeningType |
Controls the hardening function in the yield stress: 1: linear hardening/softening with cutoff at zero stress. More... | |
double | kappaC |
Parameter of the exponential laws. More... | |
double | hardeningModulus |
Hardening modulus normalized with the elastic modulus, parameter of the linear hardening/softening law. More... | |
double | limitYieldStress |
Parameter of the exponential hardening law. More... | |
double | initialYieldStress |
Parameter of all three laws, this is the initial value of the yield stress in pure shear. More... | |
double | alpha |
Friction coefficient, parameter of the yield criterion. More... | |
double | alphaPsi |
Dilatancy coefficient, parameter of the flow rule. More... | |
double | kFactor |
Scalar factor between rate of plastic multiplier and rate of hardening variable. More... | |
IsotropicLinearElasticMaterial * | LEMaterial |
Associated linear elastic material. More... | |
double | yieldTol |
Yield tolerance. More... | |
int | newtonIter |
Maximum number of iterations for stress return. 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... | |
This class implements a (local) nonassociated plasticity model based on the Drucker-Prager yield criterion with hardening and softening.
Definition at line 194 of file druckerPragerPlasticitySM.h.
Values of history variable state_flag.
Enumerator | |
---|---|
DP_Elastic | |
DP_Unloading | |
DP_Yielding | |
DP_Vertex |
Definition at line 198 of file druckerPragerPlasticitySM.h.
oofem::DruckerPragerPlasticitySM::DruckerPragerPlasticitySM | ( | int | n, |
Domain * | d | ||
) |
Constructor.
Definition at line 201 of file druckerPragerPlasticitySM.C.
References kFactor, LEMaterial, newtonIter, and yieldTol.
|
virtual |
bool oofem::DruckerPragerPlasticitySM::checkForVertexCase | ( | double | eM, |
double | gM, | ||
double | kM, | ||
double | trialStressJTwo, | ||
double | volumetricStress, | ||
double | tempKappa | ||
) |
Check if the trial stress state falls within the vertex region.
eM | Elasticity modulus. |
gM | Shear modulus. |
kM | Bulk modulus. |
Definition at line 368 of file druckerPragerPlasticitySM.C.
References alphaPsi, computeYieldValue(), kFactor, and yieldTol.
Referenced by performLocalStressReturn().
|
virtual |
Compute the current yield stress in pure shear of the Drucker-Prager model according to the used hardening law.
The yield stress is tauY in f(sigma, kappa) = F(sigma) - tauY(kappa).
kappa | Hardening variable. |
eM | Elasticity modulus. |
Definition at line 525 of file druckerPragerPlasticitySM.C.
References hardeningModulus, hardeningType, initialYieldStress, kappaC, limitYieldStress, and OOFEM_ERROR.
Referenced by computeYieldValue().
|
virtual |
Compute derivative of yield stress with respect to the hardening variable kappa.
kappa | Hardening variable. |
eM | Elasticity modulus. |
Definition at line 553 of file druckerPragerPlasticitySM.C.
References hardeningModulus, hardeningType, initialYieldStress, kappaC, limitYieldStress, and OOFEM_ERROR.
Referenced by giveRegAlgorithmicStiffMatrix(), giveVertexAlgorithmicStiffMatrix(), performRegularReturn(), and performVertexReturn().
double oofem::DruckerPragerPlasticitySM::computeYieldValue | ( | double | meanStress, |
double | JTwo, | ||
double | kappa, | ||
double | eM | ||
) | const |
Compute the yield value based on stress and hardening variable.
meanStress | 1/3 of trace of sigma. |
JTwo | Second deviatoric invariant. |
kappa | Hardening variable. |
eM | Elasticity modulus. |
Definition at line 516 of file druckerPragerPlasticitySM.C.
References alpha, and computeYieldStressInShear().
Referenced by checkForVertexCase(), performLocalStressReturn(), performRegularReturn(), and performVertexReturn().
|
virtual |
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 801 of file druckerPragerPlasticitySM.C.
References oofem::FEMComponent::giveDomain().
|
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 579 of file druckerPragerPlasticitySM.C.
References oofem::DruckerPragerPlasticitySMStatus::DP_Elastic, oofem::DruckerPragerPlasticitySMStatus::DP_Unloading, oofem::DruckerPragerPlasticitySMStatus::DP_Vertex, oofem::DruckerPragerPlasticitySMStatus::DP_Yielding, oofem::IsotropicLinearElasticMaterial::give3dMaterialStiffnessMatrix(), giveRegAlgorithmicStiffMatrix(), oofem::Material::giveStatus(), giveVertexAlgorithmicStiffMatrix(), LEMaterial, and OOFEM_ERROR.
|
inlinevirtual |
Reimplemented from oofem::StructuralMaterial.
Definition at line 238 of file druckerPragerPlasticitySM.h.
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 239 of file druckerPragerPlasticitySM.h.
References _IFT_DruckerPragerPlasticitySM_Name, oofem::IntegrationPointStatus::gp, oofem::DruckerPragerPlasticitySMStatus::kappa, oofem::StructuralMaterialStatus::strainVector, and oofem::DruckerPragerPlasticitySMStatus::tempKappa.
|
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 776 of file druckerPragerPlasticitySM.C.
References oofem::FloatArray::at(), oofem::StructuralMaterial::giveIPValue(), oofem::DruckerPragerPlasticitySMStatus::giveKappa(), oofem::DruckerPragerPlasticitySMStatus::givePlasticStrainVector(), oofem::Material::giveStatus(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
virtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 265 of file druckerPragerPlasticitySM.C.
References oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector_3d(), oofem::StructuralMaterialStatus::giveTempStressVector(), oofem::Material::initTempStatus(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), and performLocalStressReturn().
void oofem::DruckerPragerPlasticitySM::giveRegAlgorithmicStiffMatrix | ( | FloatMatrix & | answer, |
MatResponseMode | mode, | ||
GaussPoint * | gp, | ||
TimeStep * | tStep | ||
) |
Compute and give back algorithmic stiffness matrix for the regular case (no vertex).
answer | Consistent stiffness matrix. |
mode | Material reponse mode. |
gp | Gauss point. |
tStep | Time step. |
Definition at line 625 of file druckerPragerPlasticitySM.C.
References alpha, alphaPsi, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeSecondStressInvariant(), computeYieldStressPrime(), Ex, oofem::IsotropicLinearElasticMaterial::give(), oofem::IsotropicLinearElasticMaterial::give3dMaterialStiffnessMatrix(), oofem::DruckerPragerPlasticitySMStatus::giveKappa(), oofem::Material::giveStatus(), oofem::DruckerPragerPlasticitySMStatus::giveTempKappa(), oofem::StructuralMaterialStatus::giveTempStressVector(), kFactor, LEMaterial, NYxz, OOFEM_ERROR, oofem::FloatMatrix::solveForRhs(), oofem::FloatArray::times(), and oofem::FloatMatrix::zero().
Referenced by give3dMaterialStiffnessMatrix().
|
inlinevirtual |
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis.
answer | Vector of thermal dilatation coefficients. |
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 336 of file druckerPragerPlasticitySM.h.
References oofem::IsotropicLinearElasticMaterial::giveThermalDilatationVector().
void oofem::DruckerPragerPlasticitySM::giveVertexAlgorithmicStiffMatrix | ( | FloatMatrix & | answer, |
MatResponseMode | mode, | ||
GaussPoint * | gp, | ||
TimeStep * | tStep | ||
) |
Compute consistent stiffness matrix for the vertex case.
answer | Consistent stiffness matrix. |
mode | Material reponse mode. |
gp | Gauss point. |
tStep | Time step. |
Definition at line 713 of file druckerPragerPlasticitySM.C.
References alpha, oofem::FloatMatrix::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeYieldStressPrime(), Ex, oofem::IsotropicLinearElasticMaterial::give(), oofem::IsotropicLinearElasticMaterial::give3dMaterialStiffnessMatrix(), oofem::DruckerPragerPlasticitySMStatus::giveKappa(), oofem::FEMComponent::giveNumber(), oofem::DruckerPragerPlasticitySMStatus::givePlasticStrainDeviator(), oofem::Material::giveStatus(), oofem::DruckerPragerPlasticitySMStatus::giveTempKappa(), oofem::StructuralMaterialStatus::giveTempStrainVector(), oofem::DruckerPragerPlasticitySMStatus::giveTempVolumetricPlasticStrain(), oofem::DruckerPragerPlasticitySMStatus::giveVolumetricPlasticStrain(), LEMaterial, NYxz, OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatArray::subtract(), oofem::FloatMatrix::times(), and oofem::FloatMatrix::zero().
Referenced by give3dMaterialStiffnessMatrix().
|
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 216 of file druckerPragerPlasticitySM.C.
References _IFT_DruckerPragerPlasticitySM_alpha, _IFT_DruckerPragerPlasticitySM_alphapsi, _IFT_DruckerPragerPlasticitySM_hm, _IFT_DruckerPragerPlasticitySM_ht, _IFT_DruckerPragerPlasticitySM_iys, _IFT_DruckerPragerPlasticitySM_kc, _IFT_DruckerPragerPlasticitySM_lys, _IFT_DruckerPragerPlasticitySM_newtoniter, _IFT_DruckerPragerPlasticitySM_yieldtol, alpha, alphaPsi, hardeningModulus, hardeningType, oofem::IsotropicLinearElasticMaterial::initializeFrom(), oofem::StructuralMaterial::initializeFrom(), initialYieldStress, IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_BAD_FORMAT, oofem::IRRT_OK, kappaC, kFactor, LEMaterial, limitYieldStress, newtonIter, OOFEM_WARNING, and yieldTol.
|
inlinevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 334 of file druckerPragerPlasticitySM.h.
void oofem::DruckerPragerPlasticitySM::performLocalStressReturn | ( | GaussPoint * | gp, |
const FloatArray & | strain | ||
) |
Perform a standard local stress return using the function computeYieldValue at the specified Gauss point.
This function computes the history variables, i.e. the temp variable of plastic strain, hardening variable, state flag, and the temp stress, and stores them in the temp-status.
gp | Gauss point. |
strain | Strain vector of this Gauss point. |
Definition at line 295 of file druckerPragerPlasticitySM.C.
References oofem::StructuralMaterial::applyDeviatoricElasticCompliance(), oofem::StructuralMaterial::applyDeviatoricElasticStiffness(), checkForVertexCase(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeDeviatoricVolumetricSum(), oofem::StructuralMaterial::computeSecondStressInvariant(), computeYieldValue(), oofem::DruckerPragerPlasticitySMStatus::DP_Elastic, oofem::DruckerPragerPlasticitySMStatus::DP_Unloading, oofem::DruckerPragerPlasticitySMStatus::DP_Vertex, oofem::DruckerPragerPlasticitySMStatus::DP_Yielding, Ex, oofem::IsotropicLinearElasticMaterial::give(), oofem::DruckerPragerPlasticitySMStatus::giveKappa(), oofem::DruckerPragerPlasticitySMStatus::givePlasticStrainDeviator(), oofem::DruckerPragerPlasticitySMStatus::giveStateFlag(), oofem::Material::giveStatus(), oofem::DruckerPragerPlasticitySMStatus::giveVolumetricPlasticStrain(), LEMaterial, oofem::DruckerPragerPlasticitySMStatus::letTempKappaBe(), oofem::DruckerPragerPlasticitySMStatus::letTempPlasticStrainDeviatorBe(), oofem::DruckerPragerPlasticitySMStatus::letTempStateFlagBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::DruckerPragerPlasticitySMStatus::letTempVolumetricPlasticStrainBe(), NYxz, performRegularReturn(), performVertexReturn(), oofem::FloatArray::subtract(), and yieldTol.
Referenced by giveRealStressVector_3d().
void oofem::DruckerPragerPlasticitySM::performRegularReturn | ( | double | eM, |
double | gM, | ||
double | kM, | ||
double | trialStressJTwo, | ||
FloatArray & | stressDeviator, | ||
double & | volumetricStress, | ||
double & | tempKappa | ||
) |
Perform stress return for regular case, i.e.
if the trial stress state does not lie within the vertex region.
eM | Elasticity modulus. |
gM | Shear modulus. |
kM | Bulk modulus. |
Definition at line 390 of file druckerPragerPlasticitySM.C.
References oofem::FloatArray::add(), alpha, alphaPsi, oofem::StructuralMaterial::computeSecondStressInvariant(), computeYieldStressPrime(), computeYieldValue(), kFactor, newtonIter, OOFEM_ERROR, OOFEM_LOG_DEBUG, oofem::FloatArray::times(), and yieldTol.
Referenced by performLocalStressReturn().
void oofem::DruckerPragerPlasticitySM::performVertexReturn | ( | double | eM, |
double | gM, | ||
double | kM, | ||
double | trialStressJTwo, | ||
FloatArray & | stressDeviator, | ||
double & | volumetricStress, | ||
double & | tempKappa, | ||
double | volumetricElasticTrialStrain, | ||
double | kappa | ||
) |
Perform stress return for vertex case, i.e.
if the trial stress state lies within the vertex region.
eM | Elasticity modulus. |
gM | Shear modulus. |
kM | Bulk modulus. |
Definition at line 458 of file druckerPragerPlasticitySM.C.
References alpha, computeYieldStressPrime(), computeYieldValue(), newtonIter, OOFEM_ERROR, OOFEM_LOG_DEBUG, yieldTol, and oofem::FloatArray::zero().
Referenced by performLocalStressReturn().
|
virtual |
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.
Reimplemented from oofem::Material.
Definition at line 808 of file druckerPragerPlasticitySM.C.
References oofem::DruckerPragerPlasticitySMStatus::DP_Vertex, oofem::DruckerPragerPlasticitySMStatus::DP_Yielding, oofem::DruckerPragerPlasticitySMStatus::giveStateFlag(), and oofem::Material::giveStatus().
|
inlinevirtual |
Returns the relative redistribution cost of the receiver.
Reimplemented from oofem::Material.
Definition at line 344 of file druckerPragerPlasticitySM.h.
|
protected |
Friction coefficient, parameter of the yield criterion.
Definition at line 216 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldValue(), giveRegAlgorithmicStiffMatrix(), giveVertexAlgorithmicStiffMatrix(), initializeFrom(), performRegularReturn(), and performVertexReturn().
|
protected |
Dilatancy coefficient, parameter of the flow rule.
Definition at line 218 of file druckerPragerPlasticitySM.h.
Referenced by checkForVertexCase(), giveRegAlgorithmicStiffMatrix(), initializeFrom(), and performRegularReturn().
|
protected |
Hardening modulus normalized with the elastic modulus, parameter of the linear hardening/softening law.
Definition at line 210 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Controls the hardening function in the yield stress: 1: linear hardening/softening with cutoff at zero stress.
2: exponential hardening/softening to limitYieldStress.
Definition at line 206 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Parameter of all three laws, this is the initial value of the yield stress in pure shear.
Definition at line 214 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Parameter of the exponential laws.
Definition at line 208 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Scalar factor between rate of plastic multiplier and rate of hardening variable.
Definition at line 220 of file druckerPragerPlasticitySM.h.
Referenced by checkForVertexCase(), DruckerPragerPlasticitySM(), giveRegAlgorithmicStiffMatrix(), initializeFrom(), and performRegularReturn().
|
protected |
Associated linear elastic material.
Definition at line 223 of file druckerPragerPlasticitySM.h.
Referenced by DruckerPragerPlasticitySM(), give3dMaterialStiffnessMatrix(), giveRegAlgorithmicStiffMatrix(), giveVertexAlgorithmicStiffMatrix(), initializeFrom(), performLocalStressReturn(), and ~DruckerPragerPlasticitySM().
|
protected |
Parameter of the exponential hardening law.
Definition at line 212 of file druckerPragerPlasticitySM.h.
Referenced by computeYieldStressInShear(), computeYieldStressPrime(), and initializeFrom().
|
protected |
Maximum number of iterations for stress return.
Definition at line 228 of file druckerPragerPlasticitySM.h.
Referenced by DruckerPragerPlasticitySM(), initializeFrom(), performRegularReturn(), and performVertexReturn().
|
protected |
Yield tolerance.
Definition at line 226 of file druckerPragerPlasticitySM.h.
Referenced by checkForVertexCase(), DruckerPragerPlasticitySM(), initializeFrom(), performLocalStressReturn(), performRegularReturn(), and performVertexReturn().