OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::MPlasticMaterial2 Class Referenceabstract

This class represents a base class for non-associated multisurface plasticity. More...

#include <mplasticmaterial2.h>

+ Inheritance diagram for oofem::MPlasticMaterial2:
+ Collaboration diagram for oofem::MPlasticMaterial2:

Public Member Functions

 MPlasticMaterial2 (int n, Domain *d)
 
virtual ~MPlasticMaterial2 ()
 
virtual int hasNonLinearBehaviour ()
 Returns nonzero if receiver is non linear. More...
 
virtual int hasMaterialModeCapability (MaterialMode mode)
 Tests if material supports material mode. More...
 
virtual const char * giveClassName () const
 
LinearElasticMaterialgiveLinearElasticMaterial ()
 Returns reference to undamaged (bulk) material. 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 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)
 
virtual int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. More...
 
virtual int giveSizeOfFullHardeningVarsVector ()
 
virtual int giveSizeOfReducedHardeningVarsVector (GaussPoint *) const
 
virtual MaterialStatusCreateStatus (GaussPoint *gp) const
 Creates new copy of associated status and inserts it into given integration point. More...
 
- Public Member Functions inherited from oofem::StructuralMaterial
 StructuralMaterial (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralMaterial ()
 Destructor. More...
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual void giveInputRecord (DynamicInputRecord &input)
 Setups the input record string of receiver. More...
 
virtual void 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 MaterialStatusgiveStatus (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...
 
virtual const char * giveInputRecordName () const =0
 
DomaingiveDomain () 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 InterfacegiveInterface (InterfaceType t)
 Interface requesting service. More...
 
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros). More...
 

Protected Types

enum  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 }
 

Protected Member Functions

virtual int giveMaxNumberOfActiveYieldConds (GaussPoint *gp)=0
 
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)
 
virtual double computeYieldValueAt (GaussPoint *gp, int isurf, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables)=0
 Computes the value of yield function. More...
 
virtual void computeStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables)=0
 Computes the stress gradient of yield/loading function (df/d_sigma). More...
 
void computeReducedStressGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, const FloatArray &stressVector, const FloatArray &strainSpaceHardeningVariables)
 
virtual void computeStrainHardeningVarsIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &stress, const FloatArray &dlambda, const FloatArray &dplasticStrain, const IntArray &activeConditionMap)=0
 Computes the increment of strain-space hardening variables. More...
 
virtual void computeKGradientVector (FloatArray &answer, functType ftype, int isurf, GaussPoint *gp, FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables)=0
 Computes the derivative of yield/loading function with respect to $ \kappa $ vector. More...
 
virtual void computeReducedHardeningVarsSigmaGradient (FloatMatrix &answer, GaussPoint *gp, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &gamma)=0
 Computes derivative of $ \kappa $ vector with respect to stress. More...
 
virtual void computeReducedHardeningVarsLamGradient (FloatMatrix &answer, GaussPoint *gp, int actSurf, const IntArray &activeConditionMap, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVars, const FloatArray &gamma)=0
 computes derivative of $ \kappa $ vector with respect to lambda vector More...
 
virtual int hasHardening ()=0
 Indicates, whether receiver model has hardening/softening behavior or behaves according to perfect plasticity theory. More...
 
virtual void computeReducedSSGradientMatrix (FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables)=0
 Computes second derivative of loading function with respect to stress. More...
 
virtual void computeReducedSKGradientMatrix (FloatMatrix &gradientMatrix, int i, GaussPoint *gp, const FloatArray &fullStressVector, const FloatArray &strainSpaceHardeningVariables)=0
 Computes second derivative of loading function with respect to stress and hardening vars. More...
 
virtual void computeTrialStressIncrement (FloatArray &answer, GaussPoint *gp, const FloatArray &strainIncrement, TimeStep *tStep)
 Computes full-space trial stress increment (elastic). More...
 
virtual void computeReducedElasticModuli (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
 
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

LinearElasticMateriallinearElasticMaterial
 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...
 
Domaindomain
 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...
 

Detailed Description

This class represents a base class for non-associated multisurface plasticity.

The Multisurface plasticity is characterized by the following: Let $\sigma, \varepsilon $, and $\varepsilon^p $ be the stress, total strain, and plastic strain vectors, respectively. It is assumed that the total strain is decomposed into reversible elastic and irreversible plastic parts $\varepsilon=\varepsilon^e+\varepsilon^p $. The elastic response is characterized in terms of elastic constitutive matrix $ D^e $ as $ \sigma=D^e(\varepsilon-\varepsilon^e) $ As long as the stress remains inside the elastic domain, the deformation process is purely elastic and the plastic strain does not change.

It is assumed that the elastic domain, denoted as $ IE $ is bounded by a composite yield surface. It is defined as

\[ IE=\{(\sigma,\kappa)|f_i(\sigma,\kappa)<0,\; i\in\{1,\cdots,m\}\} \]

where $ f_i(\sigma,\kappa) $ are $ m\ge1 $ yield functions intersecting in a possibly non-smooth fashion. The vector $ \kappa $ contains internal variables controlling the evolution of yield surfaces (amount of hardening or softening). The evolution of plastic strain $ \varepsilon^p $ is expressed in Koiter's form. Assuming the non-associated plasticity, this reads

\[ \label{epe} \varepsilon^p=\sum^{m}_{i=1} \lambda^i \partial_{\sigma}g_i(\sigma,\kappa) \]

where $ g_i $ are plastic potential functions. The $ \lambda^i $ are referred as plastic consistency parameters, which satisfy the following Kuhn-Tucker conditions

\[ \label{ktc} \lambda^i\ge0,\;f_i\le0,\;{\rm and}\ \lambda^i f_i=0 \]

These conditions imply that in the elastic regime the yield function must remain negative and the rate of the plastic multiplier is zero (plastic strain remains constant) while in the plastic regime the yield function must be equal to zero (stress remains on the surface) and the rate of the plastic multiplier is positive. The evolution of vector of internal hardening/softening variables $ \kappa $ is expressed in terms of a general hardening/softening law of the form

\[ \dot{\kappa} = \dot{\kappa}(\sigma, \lambda) \]

where $ \lambda $ is the vector of plastic consistency parameters $ \lambda_i $.

Definition at line 175 of file mplasticmaterial2.h.

Member Enumeration Documentation

Type that allows to distinguish between yield function and loading function.

Enumerator
yieldFunction 
loadFunction 

Definition at line 185 of file mplasticmaterial2.h.

Enumerator
associatedPT 
nonassociatedPT 

Definition at line 186 of file mplasticmaterial2.h.

Protected type to determine the return mapping algorithm.

Enumerator
mpm_ClosestPoint 
mpm_CuttingPlane 

Definition at line 183 of file mplasticmaterial2.h.

Constructor & Destructor Documentation

oofem::MPlasticMaterial2::MPlasticMaterial2 ( int  n,
Domain d 
)
oofem::MPlasticMaterial2::~MPlasticMaterial2 ( )
virtual

Definition at line 63 of file mplasticmaterial2.C.

References linearElasticMaterial.

Member Function Documentation

void oofem::MPlasticMaterial2::addNewPopulation ( IntArray mask)
protected

Definition at line 1920 of file mplasticmaterial2.C.

References getPopulationSignature(), and populationSet.

Referenced by closestPointReturn().

void oofem::MPlasticMaterial2::clearPopulationSet ( )
protected

Definition at line 1914 of file mplasticmaterial2.C.

References populationSet.

Referenced by closestPointReturn().

void oofem::MPlasticMaterial2::closestPointReturn ( FloatArray answer,
IntArray activeConditionMap,
FloatArray gamma,
GaussPoint gp,
const FloatArray totalStrain,
FloatArray plasticStrainR,
FloatArray strainSpaceHardeningVariables,
TimeStep tStep 
)
protected

Definition at line 184 of file mplasticmaterial2.C.

References oofem::FloatArray::add(), oofem::FloatMatrix::add(), addNewPopulation(), associatedPT, oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), oofem::FloatArray::beProductOf(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), clearPopulationSet(), computeAlgorithmicModuli(), computeKGradientVector(), oofem::FloatArray::computeNorm(), computeReducedElasticModuli(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedSKGradientMatrix(), computeReducedStressGradientVector(), computeResidualVector(), computeStrainHardeningVarsIncrement(), computeTrialStressIncrement(), computeYieldValueAt(), getNewPopulation(), oofem::GaussPoint::giveElement(), giveMaxNumberOfActiveYieldConds(), oofem::FEMComponent::giveNumber(), oofem::GaussPoint::giveNumber(), oofem::MPlasticMaterial2Status::givePlasticStrainVector(), oofem::FloatArray::giveSize(), oofem::Material::giveStatus(), oofem::MPlasticMaterial2Status::giveStrainSpaceHardeningVars(), hasHardening(), iterativeUpdateOfActiveConds, loadFunction, nonassociatedPT, nsurf, OOFEM_ERROR, OOFEM_WARNING, PLASTIC_MATERIAL_MAX_ITERATIONS, plType, oofem::FloatArray::printYourself(), oofem::IntArray::printYourself(), RES_TOL, oofem::IntArray::resize(), oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), oofem::FloatMatrix::solveForRhs(), oofem::FloatArray::subtract(), oofem::FloatMatrix::times(), YIELD_TOL, YIELD_TOL_SECONDARY, yieldFunction, oofem::FloatArray::zero(), oofem::FloatMatrix::zero(), and oofem::IntArray::zero().

Referenced by giveRealStressVector().

double oofem::MPlasticMaterial2::computeDamage ( GaussPoint gp,
const FloatArray strainSpaceHardeningVariables,
TimeStep tStep 
)
virtual

Definition at line 178 of file mplasticmaterial2.C.

Referenced by giveRealStressVector().

virtual void oofem::MPlasticMaterial2::computeKGradientVector ( FloatArray answer,
functType  ftype,
int  isurf,
GaussPoint gp,
FloatArray fullStressVector,
const FloatArray strainSpaceHardeningVariables 
)
protectedpure virtual

Computes the derivative of yield/loading function with respect to $ \kappa $ vector.

Implemented in oofem::DruckerPragerCutMat, oofem::Masonry02, and oofem::J2Mat.

Referenced by closestPointReturn(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().

void oofem::MPlasticMaterial2::computeReducedElasticModuli ( FloatMatrix answer,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual
virtual void oofem::MPlasticMaterial2::computeReducedHardeningVarsLamGradient ( FloatMatrix answer,
GaussPoint gp,
int  actSurf,
const IntArray activeConditionMap,
const FloatArray fullStressVector,
const FloatArray strainSpaceHardeningVars,
const FloatArray gamma 
)
protectedpure virtual

computes derivative of $ \kappa $ vector with respect to lambda vector

Implemented in oofem::DruckerPragerCutMat, oofem::Masonry02, and oofem::J2Mat.

Referenced by closestPointReturn(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().

virtual void oofem::MPlasticMaterial2::computeReducedHardeningVarsSigmaGradient ( FloatMatrix answer,
GaussPoint gp,
const IntArray activeConditionMap,
const FloatArray fullStressVector,
const FloatArray strainSpaceHardeningVars,
const FloatArray gamma 
)
protectedpure virtual
virtual void oofem::MPlasticMaterial2::computeReducedSKGradientMatrix ( FloatMatrix gradientMatrix,
int  i,
GaussPoint gp,
const FloatArray fullStressVector,
const FloatArray strainSpaceHardeningVariables 
)
protectedpure virtual

Computes second derivative of loading function with respect to stress and hardening vars.

Implemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Referenced by closestPointReturn(), computeAlgorithmicModuli(), and giveConsistentStiffnessMatrix().

virtual void oofem::MPlasticMaterial2::computeReducedSSGradientMatrix ( FloatMatrix gradientMatrix,
int  i,
GaussPoint gp,
const FloatArray fullStressVector,
const FloatArray strainSpaceHardeningVariables 
)
protectedpure virtual

Computes second derivative of loading function with respect to stress.

Implemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Referenced by computeAlgorithmicModuli().

void oofem::MPlasticMaterial2::computeReducedStressGradientVector ( FloatArray answer,
functType  ftype,
int  isurf,
GaussPoint gp,
const FloatArray stressVector,
const FloatArray strainSpaceHardeningVariables 
)
protected
void oofem::MPlasticMaterial2::computeResidualVector ( FloatArray answer,
GaussPoint gp,
const FloatArray gamma,
const IntArray activeConditionMap,
const FloatArray plasticStrainVectorR,
const FloatArray strainSpaceHardeningVariables,
std::vector< FloatArray > &  gradVec 
)
protected
virtual void oofem::MPlasticMaterial2::computeStrainHardeningVarsIncrement ( FloatArray answer,
GaussPoint gp,
const FloatArray stress,
const FloatArray dlambda,
const FloatArray dplasticStrain,
const IntArray activeConditionMap 
)
protectedpure virtual

Computes the increment of strain-space hardening variables.

Parameters
answerResult.
gpGauss point to compute at.
stressUpdated stress (corresponds to newly reached state).
dlambdaIncrement of consistency parameters.
dplasticStrainActual plastic strain increment.
activeConditionMapArray of active yield conditions.

Implemented in oofem::DruckerPragerCutMat, oofem::Masonry02, and oofem::J2Mat.

Referenced by closestPointReturn(), and cuttingPlaneReturn().

virtual void oofem::MPlasticMaterial2::computeStressGradientVector ( FloatArray answer,
functType  ftype,
int  isurf,
GaussPoint gp,
const FloatArray stressVector,
const FloatArray strainSpaceHardeningVariables 
)
protectedpure virtual

Computes the stress gradient of yield/loading function (df/d_sigma).

Implemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Referenced by computeReducedStressGradientVector().

void oofem::MPlasticMaterial2::computeTrialStressIncrement ( FloatArray answer,
GaussPoint gp,
const FloatArray strainIncrement,
TimeStep tStep 
)
protectedvirtual

Computes full-space trial stress increment (elastic).

Parameters
answerContains result in full space.
gpIntegration point.
strainIncrementStrain increment vector.
tStepSolution step.

Definition at line 1260 of file mplasticmaterial2.C.

References oofem::FloatArray::beProductOf(), computeReducedElasticModuli(), oofem::StructuralMaterial::giveFullSymVectorForm(), and oofem::GaussPoint::giveMaterialMode().

Referenced by closestPointReturn(), and cuttingPlaneReturn().

virtual double oofem::MPlasticMaterial2::computeYieldValueAt ( GaussPoint gp,
int  isurf,
const FloatArray stressVector,
const FloatArray strainSpaceHardeningVariables 
)
protectedpure virtual

Computes the value of yield function.

Implemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Referenced by closestPointReturn(), and cuttingPlaneReturn().

MaterialStatus * oofem::MPlasticMaterial2::CreateStatus ( GaussPoint gp) const
virtual

Creates new copy of associated status and inserts it into given integration point.

Parameters
gpIntegration point where newly created status will be stored.
Returns
Reference to new status.

Reimplemented from oofem::Material.

Reimplemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Definition at line 88 of file mplasticmaterial2.C.

References oofem::FEMComponent::giveDomain(), and giveSizeOfReducedHardeningVarsVector().

void oofem::MPlasticMaterial2::cuttingPlaneReturn ( FloatArray answer,
IntArray activeConditionMap,
FloatArray gamma,
GaussPoint gp,
const FloatArray totalStrain,
FloatArray plasticStrainR,
FloatArray strainSpaceHardeningVariables,
TimeStep tStep 
)
protected
int oofem::MPlasticMaterial2::getNewPopulation ( IntArray result,
IntArray candidateMask,
int  degree,
int  size 
)
protected
long oofem::MPlasticMaterial2::getPopulationSignature ( IntArray mask)
protected

Definition at line 1886 of file mplasticmaterial2.C.

References oofem::IntArray::at(), and oofem::IntArray::giveSize().

Referenced by addNewPopulation(), and getNewPopulation().

void oofem::MPlasticMaterial2::give1dStressStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Method for computing 1d stiffness matrix of receiver.

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

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 1764 of file mplasticmaterial2.C.

References oofem::StructuralMaterial::give1dStressStiffMtrx(), giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveLinearElasticMaterial(), mpm_ClosestPoint, and rmType.

void oofem::MPlasticMaterial2::give2dBeamLayerStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Method for computing 2d beam layer stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d beam layer 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 1787 of file mplasticmaterial2.C.

References oofem::StructuralMaterial::give2dBeamLayerStiffMtrx(), giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveLinearElasticMaterial(), mpm_ClosestPoint, and rmType.

void oofem::MPlasticMaterial2::give3dMaterialStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.

Parameters
answerComputed results.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented from oofem::StructuralMaterial.

Definition at line 1668 of file mplasticmaterial2.C.

References oofem::StructuralMaterial::give3dMaterialStiffnessMatrix(), giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveLinearElasticMaterial(), oofem::GaussPoint::giveMaterialMode(), mpm_ClosestPoint, OOFEM_ERROR, and rmType.

virtual const char* oofem::MPlasticMaterial2::giveClassName ( ) const
inlinevirtual
Returns
Class name of the receiver.

Reimplemented from oofem::StructuralMaterial.

Reimplemented in oofem::DruckerPragerCutMat, oofem::Masonry02, and oofem::J2Mat.

Definition at line 199 of file mplasticmaterial2.h.

void oofem::MPlasticMaterial2::giveConsistentStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Definition at line 1339 of file mplasticmaterial2.C.

References oofem::FloatArray::add(), oofem::FloatMatrix::add(), associatedPT, oofem::IntArray::at(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), oofem::FloatMatrix::beProductOf(), oofem::FloatArray::beTProductOf(), computeAlgorithmicModuli(), computeKGradientVector(), computeReducedElasticModuli(), computeReducedHardeningVarsLamGradient(), computeReducedHardeningVarsSigmaGradient(), computeReducedSKGradientMatrix(), computeReducedStressGradientVector(), oofem::StructuralMaterial::giveFullSymVectorForm(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatMatrix::giveNumberOfRows(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::MPlasticMaterial2Status::giveTempActiveConditionMap(), oofem::MPlasticMaterial2Status::giveTempGamma(), oofem::MPlasticMaterial2Status::giveTempStateFlag(), oofem::MPlasticMaterial2Status::giveTempStrainSpaceHardeningVarsVector(), oofem::StructuralMaterialStatus::giveTempStressVector(), hasHardening(), loadFunction, oofem::FloatMatrix::negated(), nonassociatedPT, nsurf, plType, oofem::MPlasticMaterial2Status::PM_Elastic, oofem::MPlasticMaterial2Status::PM_Unloading, oofem::FloatMatrix::resize(), oofem::FloatMatrix::times(), yieldFunction, and oofem::FloatMatrix::zero().

Referenced by give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), oofem::Masonry02::give2dInterfaceMaterialStiffnessMatrix(), give3dMaterialStiffnessMatrix(), giveFiberStiffMtrx(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), and givePlateLayerStiffMtrx().

void oofem::MPlasticMaterial2::giveElastoPlasticStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual
void oofem::MPlasticMaterial2::giveFiberStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Method for computing 1d fiber stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d fiber 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 1837 of file mplasticmaterial2.C.

References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), oofem::StructuralMaterial::giveFiberStiffMtrx(), giveLinearElasticMaterial(), mpm_ClosestPoint, and rmType.

int oofem::MPlasticMaterial2::giveIPValue ( FloatArray answer,
GaussPoint gp,
InternalStateType  type,
TimeStep tStep 
)
virtual

Returns the integration point corresponding value in Reduced form.

Parameters
answerContain corresponding ip value, zero sized if not available.
gpIntegration point to which the value refers.
typeDetermines the type of internal variable.
tStepDetermines the time step.
Returns
Nonzero if the assignment can be done, zero if this type of variable is not supported.
Todo:
Fill in correct full form values here! This just adds zeros!
Todo:
Fill in correct full form values here! This just adds zeros!

Reimplemented from oofem::StructuralMaterial.

Reimplemented in oofem::DruckerPragerCutMat.

Definition at line 1862 of file mplasticmaterial2.C.

References oofem::StructuralMaterial::computePrincipalValues(), oofem::StructuralMaterial::giveFullSymVectorForm(), oofem::StructuralMaterial::giveIPValue(), oofem::GaussPoint::giveMaterialMode(), oofem::MPlasticMaterial2Status::givePlasticStrainVector(), oofem::Material::giveStatus(), and oofem::principal_strain.

Referenced by oofem::DruckerPragerCutMat::giveIPValue().

virtual int oofem::MPlasticMaterial2::giveMaxNumberOfActiveYieldConds ( GaussPoint gp)
protectedpure virtual
void oofem::MPlasticMaterial2::givePlaneStrainStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Method for computing plane strain stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane strain stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix. Note: as already described, if zero strain component is imposed (Plane strain, ..) this condition must be taken into account in geometrical relations, and corresponding component has to be included in reduced vector. (So plane strain conditions are $ \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 1740 of file mplasticmaterial2.C.

References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveLinearElasticMaterial(), oofem::StructuralMaterial::givePlaneStrainStiffMtrx(), mpm_ClosestPoint, and rmType.

void oofem::MPlasticMaterial2::givePlaneStressStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Method for computing plane stress stiffness matrix of receiver.

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

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 1714 of file mplasticmaterial2.C.

References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveLinearElasticMaterial(), oofem::StructuralMaterial::givePlaneStressStiffMtrx(), mpm_ClosestPoint, and rmType.

void oofem::MPlasticMaterial2::givePlateLayerStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
protectedvirtual

Method for computing 2d plate layer stiffness matrix of receiver.

Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d plate layer 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 1813 of file mplasticmaterial2.C.

References giveConsistentStiffnessMatrix(), giveElastoPlasticStiffnessMatrix(), giveLinearElasticMaterial(), oofem::StructuralMaterial::givePlateLayerStiffMtrx(), mpm_ClosestPoint, and rmType.

void oofem::MPlasticMaterial2::giveRealStressVector ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
virtual

Computes the real stress vector for given total strain and integration point.

The total strain is defined as strain computed directly from displacement field at given time. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process.

Parameters
answerStress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress.
gpIntegration point.
reducedStrainStrain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain.
tStepCurrent time step (most models are able to respond only when tStep is current time step).
Todo:
Move this to StructuralCrossSection ?

Reimplemented from oofem::StructuralMaterial.

Definition at line 95 of file mplasticmaterial2.C.

References oofem::FloatArray::at(), closestPointReturn(), computeDamage(), cuttingPlaneReturn(), oofem::GaussPoint::giveMaterialMode(), oofem::StructuralMaterial::giveReducedSymVectorForm(), oofem::MPlasticMaterial2Status::giveStateFlag(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::Material::initTempStatus(), oofem::MPlasticMaterial2Status::letTempDamageBe(), oofem::MPlasticMaterial2Status::letTempPlasticStrainVectorBe(), oofem::MPlasticMaterial2Status::letTempStateFlagBe(), oofem::MPlasticMaterial2Status::letTempStrainSpaceHardeningVarsVectorBe(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), mpm_ClosestPoint, nsurf, oofem::MPlasticMaterial2Status::PM_Elastic, oofem::MPlasticMaterial2Status::PM_Unloading, oofem::MPlasticMaterial2Status::PM_Yielding, rmType, oofem::MPlasticMaterial2Status::setTempActiveConditionMap(), oofem::MPlasticMaterial2Status::setTempGamma(), and oofem::FloatArray::times().

virtual void oofem::MPlasticMaterial2::giveRealStressVector_1d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
inlinevirtual

Default implementation relies on giveRealStressVector_StressControl.

Reimplemented from oofem::StructuralMaterial.

Definition at line 227 of file mplasticmaterial2.h.

virtual void oofem::MPlasticMaterial2::giveRealStressVector_3d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
inlinevirtual

Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.

Reimplemented from oofem::StructuralMaterial.

Definition at line 221 of file mplasticmaterial2.h.

virtual void oofem::MPlasticMaterial2::giveRealStressVector_PlaneStrain ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
inlinevirtual

Default implementation relies on giveRealStressVector_3d.

Reimplemented from oofem::StructuralMaterial.

Definition at line 223 of file mplasticmaterial2.h.

virtual void oofem::MPlasticMaterial2::giveRealStressVector_PlaneStress ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedE,
TimeStep tStep 
)
inlinevirtual

Default implementation relies on giveRealStressVector_StressControl.

Reimplemented from oofem::StructuralMaterial.

Definition at line 225 of file mplasticmaterial2.h.

virtual int oofem::MPlasticMaterial2::giveSizeOfFullHardeningVarsVector ( )
inlinevirtual

Reimplemented in oofem::DruckerPragerCutMat, oofem::Masonry02, and oofem::J2Mat.

Definition at line 235 of file mplasticmaterial2.h.

virtual int oofem::MPlasticMaterial2::giveSizeOfReducedHardeningVarsVector ( GaussPoint ) const
inlinevirtual
virtual void oofem::MPlasticMaterial2::giveThermalDilatationVector ( FloatArray answer,
GaussPoint gp,
TimeStep tStep 
)
inlinevirtual

Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis.

Parameters
answerVector of thermal dilatation coefficients.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented from oofem::StructuralMaterial.

Definition at line 206 of file mplasticmaterial2.h.

References oofem::Material::give(), oofem::IntegrationPointStatus::gp, and tAlpha.

virtual int oofem::MPlasticMaterial2::hasHardening ( )
protectedpure virtual

Indicates, whether receiver model has hardening/softening behavior or behaves according to perfect plasticity theory.

Implemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Referenced by closestPointReturn(), computeAlgorithmicModuli(), cuttingPlaneReturn(), giveConsistentStiffnessMatrix(), and giveElastoPlasticStiffnessMatrix().

int oofem::MPlasticMaterial2::hasMaterialModeCapability ( MaterialMode  mode)
virtual

Tests if material supports material mode.

Parameters
modeRequired material mode.
Returns
Nonzero if supported, zero otherwise.

Reimplemented from oofem::StructuralMaterial.

Reimplemented in oofem::Masonry02, and oofem::DruckerPragerCutMat.

Definition at line 73 of file mplasticmaterial2.C.

virtual int oofem::MPlasticMaterial2::hasNonLinearBehaviour ( )
inlinevirtual

Returns nonzero if receiver is non linear.

Reimplemented from oofem::Material.

Reimplemented in oofem::DruckerPragerCutMat.

Definition at line 197 of file mplasticmaterial2.h.

virtual bool oofem::MPlasticMaterial2::isCharacteristicMtrxSymmetric ( MatResponseMode  rMode)
inlinevirtual

Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.

Reimplemented from oofem::Material.

Reimplemented in oofem::Masonry02, oofem::DruckerPragerCutMat, and oofem::J2Mat.

Definition at line 204 of file mplasticmaterial2.h.

int oofem::MPlasticMaterial2::testPopulation ( long  pop)
protected

Definition at line 1904 of file mplasticmaterial2.C.

References populationSet.

Referenced by getNewPopulation().

Member Data Documentation

bool oofem::MPlasticMaterial2::iterativeUpdateOfActiveConds
protected

Flag indicating whether iterative update of a set of active yield conditions takes place.

Definition at line 188 of file mplasticmaterial2.h.

Referenced by closestPointReturn(), cuttingPlaneReturn(), and MPlasticMaterial2().

std :: set< long > oofem::MPlasticMaterial2::populationSet
protected

Set for keeping record of generated populations of active yield conditions during return.

Definition at line 190 of file mplasticmaterial2.h.

Referenced by addNewPopulation(), clearPopulationSet(), getNewPopulation(), and testPopulation().


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

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:38 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011