OOFEM
2.4
OOFEM.org - Object Oriented Finite Element Solver
|
This class implements an isotropic elasto-plasto-damage material with Drucker-Prager yield condition, tension cut-off, non-associated flow rule, linear isotropic hardening and isotropic damage. More...
#include <druckerpragercutmat.h>
Public Member Functions | |
DruckerPragerCutMat (int n, Domain *d) | |
virtual | ~DruckerPragerCutMat () |
virtual int | hasMaterialModeCapability (MaterialMode mode) |
Tests if material supports material mode. More... | |
virtual IRResultType | initializeFrom (InputRecord *ir) |
Initializes receiver according to object description stored in input record. More... | |
virtual MaterialStatus * | CreateStatus (GaussPoint *gp) const |
Creates new copy of associated status and inserts it into given integration point. More... | |
virtual int | hasNonLinearBehaviour () |
Returns nonzero if receiver is non linear. More... | |
virtual bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More... | |
virtual const char * | giveClassName () const |
virtual const char * | giveInputRecordName () const |
LinearElasticMaterial * | giveLinearElasticMaterial () |
Returns a reference to the basic elastic material. More... | |
virtual int | giveSizeOfFullHardeningVarsVector () |
virtual int | giveSizeOfReducedHardeningVarsVector (GaussPoint *) const |
Public Member Functions inherited from oofem::MPlasticMaterial2 | |
MPlasticMaterial2 (int n, Domain *d) | |
virtual | ~MPlasticMaterial2 () |
LinearElasticMaterial * | giveLinearElasticMaterial () |
Returns reference to undamaged (bulk) material. 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 void | give3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode, 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 *, const FloatArray &, TimeStep *) |
Computes the real stress vector for given total strain and integration point. More... | |
virtual void | giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
virtual void | giveRealStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep) |
Default implementation relies on giveRealStressVector_3d. More... | |
virtual void | giveRealStressVector_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 double | computeDamage (GaussPoint *gp, const FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) |
Public Member Functions inherited from oofem::StructuralMaterial | |
StructuralMaterial (int n, Domain *d) | |
Constructor. More... | |
virtual | ~StructuralMaterial () |
Destructor. 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_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_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 | give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d lattice stiffness matrix of receiver. More... | |
virtual void | give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 3d lattice stiffness matrix of receiver. More... | |
virtual void | give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing stiffness matrix of plate subsoil model. More... | |
virtual void | give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Method for computing stiffness matrix of beam3d subsoil model. More... | |
virtual void | giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More... | |
virtual void | giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
Default implementation relies on giveFirstPKStressVector_3d. More... | |
virtual void | giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
virtual void | givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give1dStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
Public Member Functions inherited from oofem::Material | |
Material (int n, Domain *d) | |
Constructor. More... | |
virtual | ~Material () |
Destructor. More... | |
virtual double | give (int aProperty, GaussPoint *gp) |
Returns the value of material property 'aProperty'. More... | |
virtual bool | hasProperty (int aProperty, GaussPoint *gp) |
Returns true if 'aProperty' exists on material. More... | |
virtual void | modifyProperty (int aProperty, double value, GaussPoint *gp) |
Modify 'aProperty', which already exists on material. More... | |
double | giveCastingTime () |
virtual bool | isActivated (TimeStep *tStep) |
virtual int | hasCastingTimeSupport () |
Tests if material supports casting time. More... | |
virtual void | printYourself () |
Prints receiver state on stdout. Useful for debugging. More... | |
virtual contextIOResultType | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
Stores integration point state to output stream. More... | |
virtual contextIOResultType | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
Reads integration point state to output stream. More... | |
virtual int | checkConsistency () |
Allows programmer to test some internal data, before computation begins. More... | |
virtual int | initMaterial (Element *element) |
Optional function to call specific procedures when initializing a material. More... | |
virtual MaterialStatus * | giveStatus (GaussPoint *gp) const |
Returns material status of receiver in given integration point. More... | |
virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More... | |
virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More... | |
virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip) |
Estimates the necessary pack size to hold all packed data of receiver. More... | |
virtual double | predictRelativeComputationalCost (GaussPoint *gp) |
Returns the weight representing relative computational cost of receiver The reference material model is linear isotropic material - its weight is set to 1.0 The other material models should compare to this reference model. More... | |
virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
Returns the relative redistribution cost of the receiver. More... | |
virtual void | initTempStatus (GaussPoint *gp) |
Initializes temporary variables stored in integration point status at the beginning of new time step. More... | |
Public Member Functions inherited from oofem::FEMComponent | |
FEMComponent (int n, Domain *d) | |
Regular constructor, creates component with given number and belonging to given domain. More... | |
virtual | ~FEMComponent () |
Virtual destructor. More... | |
Domain * | giveDomain () const |
virtual void | setDomain (Domain *d) |
Sets associated Domain. More... | |
int | giveNumber () const |
void | setNumber (int num) |
Sets number of receiver. More... | |
virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
Local renumbering support. More... | |
virtual contextIOResultType | saveContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Stores receiver state to output stream. More... | |
virtual contextIOResultType | restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL) |
Restores the receiver state previously written in stream. More... | |
virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
Prints output of receiver to stream, for given time step. More... | |
virtual Interface * | giveInterface (InterfaceType t) |
Interface requesting service. More... | |
std::string | errorInfo (const char *func) const |
Returns string for prepending output (used by error reporting macros). More... | |
Protected Member Functions | |
virtual int | giveMaxNumberOfActiveYieldConds (GaussPoint *gp) |
virtual double | computeYieldValueAt (GaussPoint *gp, int isurf, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) |
Computes the value of yield function. More... | |
virtual void | computeStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &stressSpaceHardeningVars) |
Computes the stress gradient of yield/loading function (df/d_sigma). More... | |
virtual void | computeReducedSSGradientMatrix (FloatMatrix &gradientMatrix, int isurf, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) |
Computes second derivative of yield/loading function with respect to stress. More... | |
virtual void | computeReducedElasticModuli (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep) |
virtual int | hasHardening () |
Functions related to hardening. More... | |
virtual void | computeStrainHardeningVarsIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &stress, const FloatArray &dlambda, const FloatArray &dplasticStrain, const IntArray &activeConditionMap) |
Compute dot(kappa_1), dot(kappa_2) etc. More... | |
virtual void | computeKGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) |
Computes the derivative of yield/loading function with respect to kappa_1, kappa_2 etc. More... | |
virtual void | computeReducedSKGradientMatrix (FloatMatrix &gradientMatrix, int isurf, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) |
computes mixed derivative of load function with respect to stress and hardening variables More... | |
virtual void | computeReducedHardeningVarsSigmaGradient (FloatMatrix &answer, GaussPoint *gp, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &dlambda) |
computes dk(i)/dsig(j) gradient matrix More... | |
virtual void | computeReducedHardeningVarsLamGradient (FloatMatrix &answer, GaussPoint *gp, int actSurf, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &dlambda) |
computes dKappa_i/dLambda_j More... | |
virtual int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) |
Returns the integration point corresponding value in Reduced form. More... | |
Protected Member Functions inherited from oofem::MPlasticMaterial2 | |
void | closestPointReturn (FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma, GaussPoint *gp, const FloatArray &totalStrain, FloatArray &plasticStrainR, FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) |
void | cuttingPlaneReturn (FloatArray &answer, IntArray &activeConditionMap, FloatArray &gamma, GaussPoint *gp, const FloatArray &totalStrain, FloatArray &plasticStrainR, FloatArray &strainSpaceHardeningVariables, TimeStep *tStep) |
void | computeResidualVector (FloatArray &answer, GaussPoint *gp, const FloatArray &gamma, const IntArray &activeConditionMap, const FloatArray &plasticStrainVectorR, const FloatArray &strainSpaceHardeningVariables, std::vector< FloatArray > &gradVec) |
virtual void | giveConsistentStiffnessMatrix (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
virtual void | giveElastoPlasticStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
void | computeAlgorithmicModuli (FloatMatrix &answer, GaussPoint *gp, const FloatMatrix &elasticModuliInverse, const FloatArray &gamma, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables) |
void | computeReducedStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables) |
virtual void | computeTrialStressIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &strainIncrement, TimeStep *tStep) |
Computes full-space trial stress increment (elastic). More... | |
virtual void | givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane stress stiffness matrix of receiver. More... | |
virtual void | givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
Method for computing plane strain stiffness matrix of receiver. More... | |
virtual void | give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d stiffness matrix of receiver. More... | |
virtual void | give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d beam layer stiffness matrix of receiver. More... | |
virtual void | givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 2d plate layer stiffness matrix of receiver. More... | |
virtual void | giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) |
Method for computing 1d fiber stiffness matrix of receiver. More... | |
long | getPopulationSignature (IntArray &mask) |
int | testPopulation (long pop) |
void | clearPopulationSet () |
void | addNewPopulation (IntArray &mask) |
int | getNewPopulation (IntArray &result, IntArray &candidateMask, int degree, int size) |
Protected Attributes | |
double | G |
Elastic shear modulus. More... | |
double | K |
Elastic bulk modulus. More... | |
double | H |
Hardening modulus. More... | |
double | sigT |
Uniaxial tensile strength for cut-off. More... | |
double | tau0 |
Initial yield stress under pure shear. More... | |
double | alpha |
Friction coefficient. More... | |
double | alphaPsi |
Dilatancy coefficient (allowing non-associated plasticity). More... | |
double | yieldTol |
Tolerance of the error in the yield criterion. More... | |
int | newtonIter |
Maximum number of iterations in lambda search. More... | |
double | omegaCrit |
Maximum damage value. More... | |
double | a |
Parameter for damage computation from cumulative plastic strain. More... | |
Protected Attributes inherited from oofem::MPlasticMaterial2 | |
LinearElasticMaterial * | linearElasticMaterial |
Reference to bulk (undamaged) material. More... | |
int | nsurf |
Number of yield surfaces. More... | |
enum oofem::MPlasticMaterial2::ReturnMappingAlgoType | rmType |
enum oofem::MPlasticMaterial2::plastType | plType |
bool | iterativeUpdateOfActiveConds |
Flag indicating whether iterative update of a set of active yield conditions takes place. More... | |
std::set< long > | populationSet |
Set for keeping record of generated populations of active yield conditions during 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... | |
Protected Types inherited from oofem::MPlasticMaterial2 | |
enum | ReturnMappingAlgoType { mpm_ClosestPoint, mpm_CuttingPlane } |
Protected type to determine the return mapping algorithm. More... | |
enum | functType { yieldFunction, loadFunction } |
Type that allows to distinguish between yield function and loading function. More... | |
enum | plastType { associatedPT, nonassociatedPT } |
This class implements an isotropic elasto-plasto-damage material with Drucker-Prager yield condition, tension cut-off, non-associated flow rule, linear isotropic hardening and isotropic damage.
Definition at line 63 of file druckerpragercutmat.h.
oofem::DruckerPragerCutMat::DruckerPragerCutMat | ( | int | n, |
Domain * | d | ||
) |
Definition at line 49 of file druckerpragercutmat.C.
References a, G, H, K, oofem::MPlasticMaterial2::linearElasticMaterial, oofem::MPlasticMaterial2::mpm_CuttingPlane, newtonIter, oofem::MPlasticMaterial2::nonassociatedPT, oofem::MPlasticMaterial2::nsurf, omegaCrit, oofem::MPlasticMaterial2::plType, oofem::MPlasticMaterial2::rmType, sigT, and yieldTol.
|
virtual |
Definition at line 69 of file druckerpragercutmat.C.
|
protectedvirtual |
Computes the derivative of yield/loading function with respect to kappa_1, kappa_2 etc.
Implements oofem::MPlasticMaterial2.
Definition at line 255 of file druckerpragercutmat.C.
References oofem::FloatArray::at(), H, oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by hasHardening().
|
protectedvirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 236 of file druckerpragercutmat.C.
References giveLinearElasticMaterial(), and oofem::StructuralMaterial::giveStiffnessMatrix().
Referenced by giveMaxNumberOfActiveYieldConds().
|
protectedvirtual |
computes dKappa_i/dLambda_j
Implements oofem::MPlasticMaterial2.
Definition at line 283 of file druckerpragercutmat.C.
References alphaPsi, oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by hasHardening().
|
protectedvirtual |
computes dk(i)/dsig(j) gradient matrix
Implements oofem::MPlasticMaterial2.
Definition at line 275 of file druckerpragercutmat.C.
References oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveSizeOfVoigtSymVector(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by hasHardening().
|
protectedvirtual |
computes mixed derivative of load function with respect to stress and hardening variables
Implements oofem::MPlasticMaterial2.
Definition at line 267 of file druckerpragercutmat.C.
References oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveSizeOfVoigtSymVector(), oofem::FloatMatrix::resize(), and oofem::FloatMatrix::zero().
Referenced by hasHardening().
|
protectedvirtual |
Computes second derivative of yield/loading function with respect to stress.
Implements oofem::MPlasticMaterial2.
Definition at line 172 of file druckerpragercutmat.C.
References oofem::__MaterialModeToString(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeSecondStressInvariant(), oofem::GaussPoint::giveMaterialMode(), OOFEM_ERROR, oofem::FloatMatrix::resize(), oofem::FloatMatrix::symmetrized(), and oofem::FloatMatrix::zero().
Referenced by giveMaxNumberOfActiveYieldConds().
|
protectedvirtual |
Compute dot(kappa_1), dot(kappa_2) etc.
Implements oofem::MPlasticMaterial2.
Definition at line 244 of file druckerpragercutmat.C.
References alphaPsi, oofem::FloatArray::at(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by hasHardening().
|
protectedvirtual |
Computes the stress gradient of yield/loading function (df/d_sigma).
Implements oofem::MPlasticMaterial2.
Definition at line 136 of file druckerpragercutmat.C.
References alphaPsi, oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondStressInvariant(), oofem::principal_stress, and oofem::FloatArray::resize().
Referenced by giveMaxNumberOfActiveYieldConds().
|
protectedvirtual |
Computes the value of yield function.
Implements oofem::MPlasticMaterial2.
Definition at line 115 of file druckerpragercutmat.C.
References alpha, oofem::FloatArray::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computePrincipalValues(), oofem::StructuralMaterial::computeSecondStressInvariant(), H, oofem::principal_stress, sigT, and tau0.
Referenced by giveMaxNumberOfActiveYieldConds().
|
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::MPlasticMaterial2.
Definition at line 109 of file druckerpragercutmat.C.
References oofem::FEMComponent::giveDomain(), and giveSizeOfReducedHardeningVarsVector().
|
inlinevirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 115 of file druckerpragercutmat.h.
|
inlinevirtual |
Implements oofem::FEMComponent.
Definition at line 116 of file druckerpragercutmat.h.
References _IFT_DruckerPragerCutMat_Name.
|
protectedvirtual |
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::MPlasticMaterial2.
Definition at line 298 of file druckerpragercutmat.C.
References oofem::MPlasticMaterial2::giveIPValue(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
Referenced by hasHardening().
|
inline |
Returns a reference to the basic elastic material.
Definition at line 119 of file druckerpragercutmat.h.
References oofem::MPlasticMaterial2::linearElasticMaterial.
Referenced by computeReducedElasticModuli().
|
inlineprotectedvirtual |
Implements oofem::MPlasticMaterial2.
Definition at line 125 of file druckerpragercutmat.h.
References computeReducedElasticModuli(), computeReducedSSGradientMatrix(), computeStressGradientVector(), and computeYieldValueAt().
|
inlinevirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 121 of file druckerpragercutmat.h.
|
inlinevirtual |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 122 of file druckerpragercutmat.h.
Referenced by CreateStatus().
|
inlineprotectedvirtual |
Functions related to hardening.
Implements oofem::MPlasticMaterial2.
Definition at line 137 of file druckerpragercutmat.h.
References computeKGradientVector(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedSKGradientMatrix(), computeStrainHardeningVarsIncrement(), and giveIPValue().
|
virtual |
Tests if material supports material mode.
mode | Required material mode. |
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 74 of file druckerpragercutmat.C.
|
inlinevirtual |
Returns nonzero if receiver is non linear.
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 112 of file druckerpragercutmat.h.
|
virtual |
Initializes receiver according to object description stored in input record.
This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.
ir | Input record to initialize from. |
Reimplemented from oofem::StructuralMaterial.
Definition at line 81 of file druckerpragercutmat.C.
References _IFT_DruckerPragerCutMat_a, _IFT_DruckerPragerCutMat_alpha, _IFT_DruckerPragerCutMat_alphapsi, _IFT_DruckerPragerCutMat_h, _IFT_DruckerPragerCutMat_newtonIter, _IFT_DruckerPragerCutMat_omegaCrit, _IFT_DruckerPragerCutMat_sigT, _IFT_DruckerPragerCutMat_tau0, _IFT_DruckerPragerCutMat_yieldTol, a, alpha, alphaPsi, G, H, oofem::LinearElasticMaterial::initializeFrom(), oofem::StructuralMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_OK, K, oofem::MPlasticMaterial2::linearElasticMaterial, newtonIter, omegaCrit, sigT, tau0, and yieldTol.
|
inlinevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::MPlasticMaterial2.
Definition at line 113 of file druckerpragercutmat.h.
|
protected |
Parameter for damage computation from cumulative plastic strain.
Definition at line 100 of file druckerpragercutmat.h.
Referenced by DruckerPragerCutMat(), and initializeFrom().
|
protected |
Friction coefficient.
Definition at line 85 of file druckerpragercutmat.h.
Referenced by computeYieldValueAt(), and initializeFrom().
|
protected |
Dilatancy coefficient (allowing non-associated plasticity).
Definition at line 88 of file druckerpragercutmat.h.
Referenced by computeReducedHardeningVarsLamGradient(), computeStrainHardeningVarsIncrement(), computeStressGradientVector(), and initializeFrom().
|
protected |
Elastic shear modulus.
Definition at line 70 of file druckerpragercutmat.h.
Referenced by DruckerPragerCutMat(), and initializeFrom().
|
protected |
Hardening modulus.
Definition at line 76 of file druckerpragercutmat.h.
Referenced by computeKGradientVector(), computeYieldValueAt(), DruckerPragerCutMat(), and initializeFrom().
|
protected |
Elastic bulk modulus.
Definition at line 73 of file druckerpragercutmat.h.
Referenced by DruckerPragerCutMat(), and initializeFrom().
|
protected |
Maximum number of iterations in lambda search.
Definition at line 94 of file druckerpragercutmat.h.
Referenced by DruckerPragerCutMat(), and initializeFrom().
|
protected |
Maximum damage value.
Definition at line 97 of file druckerpragercutmat.h.
Referenced by DruckerPragerCutMat(), and initializeFrom().
|
protected |
Uniaxial tensile strength for cut-off.
Definition at line 79 of file druckerpragercutmat.h.
Referenced by computeYieldValueAt(), DruckerPragerCutMat(), and initializeFrom().
|
protected |
Initial yield stress under pure shear.
Definition at line 82 of file druckerpragercutmat.h.
Referenced by computeYieldValueAt(), and initializeFrom().
|
protected |
Tolerance of the error in the yield criterion.
Definition at line 91 of file druckerpragercutmat.h.
Referenced by DruckerPragerCutMat(), and initializeFrom().