OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::ConcreteDPM Class Reference

This class contains the combination of a local plasticity model for concrete with a local isotropic damage model. More...

#include <concretedpm.h>

+ Inheritance diagram for oofem::ConcreteDPM:
+ Collaboration diagram for oofem::ConcreteDPM:

Public Member Functions

 ConcreteDPM (int n, Domain *d)
 Constructor. More...
 
virtual ~ConcreteDPM ()
 Destructor. More...
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual const char * giveClassName () const
 
virtual const char * giveInputRecordName () const
 
virtual ConcreteDPMStatusgiveStatus (GaussPoint *gp) const
 Returns material status of receiver in given integration point. More...
 
LinearElasticMaterialgiveLinearElasticMaterial ()
 
virtual void giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More...
 
void performPlasticityReturn (GaussPoint *gp, FloatArray &strain)
 
bool checkForVertexCase (double &answer, double sig, double tempKappa)
 Check if the trial stress state falls within the vertex region of the plasticity model at the apex of triaxial extension or triaxial compression. More...
 
void performRegularReturn (FloatArray &stress, GaussPoint *gp)
 Perform regular stress return for the plasticity model, i.e. More...
 
void performVertexReturn (FloatArray &stress, double apexStress, GaussPoint *gp)
 Perform stress return for vertex case of the plasticity model, i.e. More...
 
double computeYieldValue (double sig, double rho, double theta, double tempKappa) const
 Compute the yield value based on stress and hardening variable. More...
 
double computeHardeningOne (double tempKappa) const
 Compute the value of the hardening function based on the hardening variable. More...
 
double computeHardeningOnePrime (double tempKappa) const
 Compute the derivative of the hardening function based on the hardening parameter. More...
 
double computeDFDKappa (double sig, double rho, double tempKappa)
 Compute the derivative of the yield surface with respect to the hardening variable based on the stress state and the hardening variable. More...
 
double computeDKappaDDeltaLambda (double sig, double rho, double tempKappa)
 Compute the derivative of kappa with respect of delta lambda based on the stress state and the hardening variable. More...
 
virtual double computeDuctilityMeasure (double sig, double rho, double theta)
 Compute the ductility measure based on the stress state. More...
 
void computeDDuctilityMeasureDInv (FloatArray &answer, double sig, double rho, double tempKappa)
 Computes the first derivative of the ductility measure with respect to the invariants sig and rho based on the stress state and the hardening parameter. More...
 
void computeAMatrix (FloatMatrix &answer, double sig, double rho, double tempKappa)
 This matrix is the core of the closest point projection and collects the derivative of the flow rule and the hardening parameters. More...
 
void computeDGDInv (FloatArray &answer, double sig, double rho, double tempKappa)
 Here, the first derivative of the plastic potential with respect to the invariants sig and rho are computed. More...
 
double computeRatioPotential (double sig, double tempKappa)
 This function computes the ratio of the volumetric and deviatoric component of the flow direction. More...
 
void computeDDGDDInv (FloatMatrix &answer, double sig, double rho, double tempKappa)
 Here, the second derivative of the plastic potential with respect to the invariants sig and rho are computed. More...
 
void computeDDGDInvDKappa (FloatArray &answer, double sig, double rho, double tempKappa)
 Here, the mixed derivative of the plastic potential with respect to the invariants and the hardening parameter are determined. More...
 
void computeDDKappaDDeltaLambdaDInv (FloatArray &answer, double sig, double rho, double tempKappa)
 Computes the mixed derivative of the hardening parameter kappa with respect to the plastic multiplier delta Lambda and the invariants sig and rho. More...
 
double computeDDKappaDDeltaLambdaDKappa (double sig, double rho, double tempKappa)
 Computes the derivative of the evolution law of the hardening parameter kappa with respect to the hardening variable kappa. More...
 
void computeDFDInv (FloatArray &answer, double sig, double rho, double tempKappa) const
 Computes the derivative of the yield surface with respect to the invariants sig and rho. More...
 
double computeTempKappa (double kappaInitial, double sigTrial, double rhoTrial, double sig)
 Compute temporary kappa. More...
 
double computeDamage (const FloatArray &strain, GaussPoint *gp, TimeStep *tStep)
 Perform stress return for the damage model, i.e. More...
 
virtual double computeDamageParam (double kappa, GaussPoint *gp)
 Compute damage parameter. More...
 
double computeInverseDamage (double dam, GaussPoint *gp)
 Compute the damage-driving variable from given damage. More...
 
virtual void computeEquivalentStrain (double &kappaD, const FloatArray &elasticStrain, GaussPoint *gp, TimeStep *tStep)
 Compute equivalent strain value. More...
 
double computeDuctilityMeasureDamage (const FloatArray &strain, GaussPoint *gp)
 Compute the ductility measure for the damage model. More...
 
void initDamaged (double kappa, const FloatArray &elasticStrain, GaussPoint *gp)
 Initialize the characteristic length, if damage is not yet activated. More...
 
void computeTrialCoordinates (const FloatArray &stress, GaussPoint *gp)
 Compute the trial coordinates. More...
 
void assignStateFlag (GaussPoint *gp)
 Assign state flag. More...
 
void computeDRhoDStress (FloatArray &answer, const FloatArray &stress) const
 Computes the derivative of rho with respect to the stress. More...
 
void computeDSigDStress (FloatArray &answer) const
 Computes the derivative of sig with respect to the stress. More...
 
void computeDDRhoDDStress (FloatMatrix &answer, const FloatArray &stress) const
 Computes the second derivative of rho with the respect to the stress. More...
 
void computeDCosThetaDStress (FloatArray &answer, const FloatArray &stress) const
 Computes the derivative of costheta with respect to the stress. More...
 
double computeDRDCosTheta (double theta, double ecc) const
 Compute the derivative of R with respect to costheta. More...
 
virtual void give3dMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point. More...
 
virtual bool isCharacteristicMtrxSymmetric (MatResponseMode rMode)
 Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. More...
 
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 int giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. More...
 
- Public Member Functions inherited from oofem::StructuralMaterial
 StructuralMaterial (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralMaterial ()
 Destructor. More...
 
virtual int hasMaterialModeCapability (MaterialMode mode)
 Tests if material supports material mode. More...
 
virtual void giveInputRecord (DynamicInputRecord &input)
 Setups the input record string of receiver. More...
 
virtual void giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 Computes the stiffness matrix for giveRealStressVector of receiver in given integration point, respecting its history. More...
 
virtual void giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
 Computes the real stress vector for given total strain and integration point. More...
 
virtual void giveRealStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_3d. More...
 
virtual void giveRealStressVector_StressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
 Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· More...
 
virtual void giveRealStressVector_ShellStressControl (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, const IntArray &strainControl, TimeStep *tStep)
 
virtual void giveRealStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_Warping (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_2dBeamLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_PlateLayer (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_Fiber (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_Lattice2d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 
virtual void giveRealStressVector_Lattice3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 
virtual void giveRealStressVector_2dPlateSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 Default implementation is not provided. More...
 
virtual void giveRealStressVector_3dBeamSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedE, TimeStep *tStep)
 
virtual void giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Prototype for computation of Eshelby stress. More...
 
void give_dPdF_from (const FloatMatrix &dSdE, FloatMatrix &answer, GaussPoint *gp)
 
void convert_dSdE_2_dPdF (FloatMatrix &answer, const FloatMatrix &dSdE, const FloatArray &S, const FloatArray &F, MaterialMode matMode)
 
virtual void giveThermalDilatationVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep)
 Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis. More...
 
double giveReferenceTemperature ()
 Returns the reference temperature of receiver. More...
 
virtual void computeStressIndependentStrainVector (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
 Computes reduced strain vector in given integration point, generated by internal processes in material, which are independent on loading in particular integration point. More...
 
virtual void computeStressIndependentStrainVector_3d (FloatArray &answer, GaussPoint *gp, TimeStep *tStep, ValueModeType mode)
 
virtual void give3dMaterialStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 
virtual void give3dMaterialStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 
void giveStressDependentPartOfStrainVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode)
 Method for subtracting from reduced space strain vector its stress-independent parts (caused by temperature, shrinkage, creep and possibly by other phenomena). More...
 
void giveStressDependentPartOfStrainVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrainVector, TimeStep *tStep, ValueModeType mode)
 
virtual void give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d beam layer stiffness matrix of receiver. More...
 
virtual void givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d plate layer stiffness matrix of receiver. More...
 
virtual void giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d fiber stiffness matrix of receiver. More...
 
virtual void give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d lattice stiffness matrix of receiver. More...
 
virtual void give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 3d lattice stiffness matrix of receiver. More...
 
virtual void give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of plate subsoil model. More...
 
virtual void give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of beam3d subsoil model. More...
 
virtual void giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More...
 
virtual void giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane stress stiffness matrix of receiver. More...
 
virtual void givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane strain stiffness matrix of receiver. More...
 
virtual void givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d stiffness matrix of receiver. More...
 
virtual void give1dStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
- Public Member Functions inherited from oofem::Material
 Material (int n, Domain *d)
 Constructor. More...
 
virtual ~Material ()
 Destructor. More...
 
virtual double give (int aProperty, GaussPoint *gp)
 Returns the value of material property 'aProperty'. More...
 
virtual bool hasProperty (int aProperty, GaussPoint *gp)
 Returns true if 'aProperty' exists on material. More...
 
virtual void modifyProperty (int aProperty, double value, GaussPoint *gp)
 Modify 'aProperty', which already exists on material. More...
 
double giveCastingTime ()
 
virtual bool isActivated (TimeStep *tStep)
 
virtual int hasNonLinearBehaviour ()
 Returns nonzero if receiver is non linear. More...
 
virtual int hasCastingTimeSupport ()
 Tests if material supports casting time. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual contextIOResultType saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Stores integration point state to output stream. More...
 
virtual contextIOResultType restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Reads integration point state to output stream. More...
 
virtual int checkConsistency ()
 Allows programmer to test some internal data, before computation begins. More...
 
virtual int initMaterial (Element *element)
 Optional function to call specific procedures when initializing a material. More...
 
virtual 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...
 
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  Concrete_VertexType { VT_Regular, VT_Tension, VT_Compression }
 

Protected Member Functions

virtual MaterialStatusCreateStatus (GaussPoint *gp) const
 Creates new copy of associated status and inserts it into given integration point. More...
 

Protected Attributes

Concrete_VertexType vertexType
 
double fc
 Parameters of the yield surface of the plasticity model. More...
 
double ft
 
double ecc
 
double AHard
 Parameter of the ductilityMeasure of the plasticity model. More...
 
double BHard
 Parameter of the ductilityMeasure of the plasticity model. More...
 
double CHard
 Parameter of the ductilityMeasure of the plasticity model. More...
 
double DHard
 Parameter of the ductilityMeasure of the plasticity model. More...
 
double ASoft
 Parameter of the ductilityMeasure of the damage model. More...
 
double yieldHardInitial
 Parameter of the hardening law of the plasticity model. More...
 
double dilationConst
 Control parameter for te volumetric plastic flow of the plastic potential. More...
 
double deltaLambda
 Plastic multiplier of the plasticity model. More...
 
double sig
 the volumetric stress. More...
 
double rho
 The length of the deviatoric stress. More...
 
double thetaTrial
 The lode angle of the trial stress. More...
 
double m
 The friction parameter of the yield surface. More...
 
double mQ
 The dilation parameter of the plastic potential. More...
 
double helem
 Element size (to be used in fracture energy approach (crack band). More...
 
LinearElasticMateriallinearElasticMaterial
 Pointer for linear elastic material. More...
 
double eM
 Elastic Young's modulus. More...
 
double gM
 Elastic shear modulus. More...
 
double kM
 Elastic bulk modulus. More...
 
double nu
 Elastic Poisson's ration. More...
 
double kappaP
 Hardening variable of plasticity model. More...
 
double tempKappaP
 
double kappaD
 Hardening variable of damage model. More...
 
double tempKappaD
 
double damage
 Damage variable of damage model. More...
 
double tempDamage
 
double ef
 Control parameter for the exponential softening law. More...
 
double yieldTol
 Yield tolerance for the plasticity model. More...
 
int newtonIter
 Maximum number of iterations for stress return. More...
 
FloatArray effectiveStress
 Stress and its deviatoric part. More...
 
double href
 Material parameter of the size-dependent adjustment (reference element size) 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 contains the combination of a local plasticity model for concrete with a local isotropic damage model.

The yield surface of the plasticity model is based on the extension of the Menetrey and Willam yield criterion. The flow rule is nonassociated. The evolution laws of the hardening variables depend on the stress state. The plasticity model describes only hardening and perfect plasticity. It is based on h effective stress. The damage parameter of the isotropic damage model is based on the total volumetric strain. An exponential softening law is implemented.

Author
Peter Grassl

Definition at line 365 of file concretedpm.h.

Member Enumeration Documentation

Enumerator
VT_Regular 
VT_Tension 
VT_Compression 

Definition at line 368 of file concretedpm.h.

Constructor & Destructor Documentation

oofem::ConcreteDPM::ConcreteDPM ( int  n,
Domain d 
)

Constructor.

Definition at line 316 of file concretedpm.C.

References damage, helem, href, kappaD, kappaP, linearElasticMaterial, newtonIter, tempDamage, tempKappaD, tempKappaP, and yieldTol.

oofem::ConcreteDPM::~ConcreteDPM ( )
virtual

Destructor.

Definition at line 332 of file concretedpm.C.

References linearElasticMaterial.

Member Function Documentation

bool oofem::ConcreteDPM::checkForVertexCase ( double &  answer,
double  sig,
double  tempKappa 
)

Check if the trial stress state falls within the vertex region of the plasticity model at the apex of triaxial extension or triaxial compression.

Returns
true for vertex case and false if regular stress return can be used.
Parameters
answerThe volumetric apex stress.
sigThe volumetric stress.
tempKappaThe hardening variable.

Definition at line 782 of file concretedpm.C.

References computeHardeningOne(), fc, m, newtonIter, vertexType, VT_Compression, VT_Regular, VT_Tension, and yieldTol.

Referenced by performPlasticityReturn().

void oofem::ConcreteDPM::computeAMatrix ( FloatMatrix answer,
double  sig,
double  rho,
double  tempKappa 
)

This matrix is the core of the closest point projection and collects the derivative of the flow rule and the hardening parameters.

Definition at line 1530 of file concretedpm.C.

References oofem::FloatMatrix::beInverseOf(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), deltaLambda, gM, kM, and oofem::FloatMatrix::zero().

Referenced by performRegularReturn().

double oofem::ConcreteDPM::computeDamage ( const FloatArray strain,
GaussPoint gp,
TimeStep tStep 
)

Perform stress return for the damage model, i.e.

if the trial stress state does not violate the plasticity surface.

Parameters
strainStrain.
gpGauss point.
tStepTime step.
Returns
Damage.

Definition at line 485 of file concretedpm.C.

References computeDamageParam(), computeEquivalentStrain(), oofem::ConcreteDPMStatus::giveDamage(), oofem::ConcreteDPMStatus::giveEpsLoc(), oofem::ConcreteDPMStatus::giveKappaD(), giveStatus(), href, initDamaged(), tempDamage, and tempKappaD.

Referenced by giveRealStressVector_3d().

double oofem::ConcreteDPM::computeDamageParam ( double  kappa,
GaussPoint gp 
)
virtual

Compute damage parameter.

Definition at line 584 of file concretedpm.C.

References DPM_DAMAGE_TOLERANCE, ef, eM, ft, oofem::ConcreteDPMStatus::giveEpsLoc(), oofem::ConcreteDPMStatus::giveLe(), giveStatus(), helem, href, and OOFEM_ERROR.

Referenced by computeDamage().

void oofem::ConcreteDPM::computeDCosThetaDStress ( FloatArray answer,
const FloatArray stress 
) const
void oofem::ConcreteDPM::computeDDGDDInv ( FloatMatrix answer,
double  sig,
double  rho,
double  tempKappa 
)

Here, the second derivative of the plastic potential with respect to the invariants sig and rho are computed.

Definition at line 1482 of file concretedpm.C.

References computeHardeningOne(), dilationConst, fc, ft, and m.

Referenced by computeAMatrix(), and computeDDKappaDDeltaLambdaDInv().

void oofem::ConcreteDPM::computeDDGDInvDKappa ( FloatArray answer,
double  sig,
double  rho,
double  tempKappa 
)

Here, the mixed derivative of the plastic potential with respect to the invariants and the hardening parameter are determined.

Definition at line 1446 of file concretedpm.C.

References computeHardeningOne(), computeHardeningOnePrime(), dilationConst, fc, ft, m, and mQ.

Referenced by computeAMatrix(), and computeDDKappaDDeltaLambdaDKappa().

void oofem::ConcreteDPM::computeDDKappaDDeltaLambdaDInv ( FloatArray answer,
double  sig,
double  rho,
double  tempKappa 
)

Computes the mixed derivative of the hardening parameter kappa with respect to the plastic multiplier delta Lambda and the invariants sig and rho.

Definition at line 1253 of file concretedpm.C.

References computeDDGDDInv(), computeDDuctilityMeasureDInv(), computeDGDInv(), computeDuctilityMeasure(), thetaTrial, and oofem::FloatArray::zero().

Referenced by computeAMatrix().

double oofem::ConcreteDPM::computeDDKappaDDeltaLambdaDKappa ( double  sig,
double  rho,
double  tempKappa 
)

Computes the derivative of the evolution law of the hardening parameter kappa with respect to the hardening variable kappa.

Definition at line 1293 of file concretedpm.C.

References computeDDGDInvDKappa(), computeDGDInv(), computeDuctilityMeasure(), and thetaTrial.

Referenced by computeAMatrix().

void oofem::ConcreteDPM::computeDDRhoDDStress ( FloatMatrix answer,
const FloatArray stress 
) const
void oofem::ConcreteDPM::computeDDuctilityMeasureDInv ( FloatArray answer,
double  sig,
double  rho,
double  tempKappa 
)

Computes the first derivative of the ductility measure with respect to the invariants sig and rho based on the stress state and the hardening parameter.

Definition at line 1356 of file concretedpm.C.

References AHard, BHard, CHard, DHard, fc, and thetaTrial.

Referenced by computeDDKappaDDeltaLambdaDInv().

void oofem::ConcreteDPM::computeDFDInv ( FloatArray answer,
double  sig,
double  rho,
double  tempKappa 
) const

Computes the derivative of the yield surface with respect to the invariants sig and rho.

Definition at line 1205 of file concretedpm.C.

References oofem::AL, computeHardeningOne(), ecc, fc, m, and thetaTrial.

Referenced by performRegularReturn().

double oofem::ConcreteDPM::computeDFDKappa ( double  sig,
double  rho,
double  tempKappa 
)

Compute the derivative of the yield surface with respect to the hardening variable based on the stress state and the hardening variable.

Parameters
sigVolumetric stress.
rhoDeviatoric length.
tempKappaHardening variable.
Returns
The derivative of the yield surface.

Definition at line 1159 of file concretedpm.C.

References computeHardeningOne(), computeHardeningOnePrime(), ecc, fc, m, and thetaTrial.

Referenced by performRegularReturn().

void oofem::ConcreteDPM::computeDGDInv ( FloatArray answer,
double  sig,
double  rho,
double  tempKappa 
)

Here, the first derivative of the plastic potential with respect to the invariants sig and rho are computed.

Definition at line 1385 of file concretedpm.C.

References computeHardeningOne(), dilationConst, fc, ft, m, and mQ.

Referenced by computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDKappaDDeltaLambda(), and performRegularReturn().

double oofem::ConcreteDPM::computeDKappaDDeltaLambda ( double  sig,
double  rho,
double  tempKappa 
)

Compute the derivative of kappa with respect of delta lambda based on the stress state and the hardening variable.

Parameters
sigVolumetric stress.
rhoLength of the deviatoric stress.
tempKappaHardening variable.
Returns
Derivative of kappa with respect to delta lambda.

Definition at line 1234 of file concretedpm.C.

References computeDGDInv(), computeDuctilityMeasure(), and thetaTrial.

Referenced by performRegularReturn().

double oofem::ConcreteDPM::computeDRDCosTheta ( double  theta,
double  ecc 
) const

Compute the derivative of R with respect to costheta.

Definition at line 1797 of file concretedpm.C.

References ecc.

void oofem::ConcreteDPM::computeDRhoDStress ( FloatArray answer,
const FloatArray stress 
) const

Computes the derivative of rho with respect to the stress.

Definition at line 1661 of file concretedpm.C.

References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeSecondCoordinate(), rho, and oofem::FloatArray::times().

Referenced by computeDCosThetaDStress().

void oofem::ConcreteDPM::computeDSigDStress ( FloatArray answer) const

Computes the derivative of sig with respect to the stress.

Definition at line 1686 of file concretedpm.C.

double oofem::ConcreteDPM::computeDuctilityMeasure ( double  sig,
double  rho,
double  theta 
)
virtual

Compute the ductility measure based on the stress state.

Parameters
sigVolumetric stress.
rhoLength of the deviatoric strength.
thetaLode angle of stress state.
Returns
Ductility measure.

Definition at line 1328 of file concretedpm.C.

References AHard, BHard, CHard, DHard, and fc.

Referenced by computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDKappaDDeltaLambda(), and computeTempKappa().

double oofem::ConcreteDPM::computeDuctilityMeasureDamage ( const FloatArray strain,
GaussPoint gp 
)
double oofem::ConcreteDPM::computeHardeningOne ( double  tempKappa) const

Compute the value of the hardening function based on the hardening variable.

Parameters
tempKappaHardening variable.
Returns
Value of the hardening function.

Definition at line 1567 of file concretedpm.C.

References yieldHardInitial.

Referenced by checkForVertexCase(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeRatioPotential(), and computeYieldValue().

double oofem::ConcreteDPM::computeHardeningOnePrime ( double  tempKappa) const

Compute the derivative of the hardening function based on the hardening parameter.

Parameters
tempKappaHardening variable.
Returns
The derivative of the hardening function.

Definition at line 1581 of file concretedpm.C.

References yieldHardInitial.

Referenced by computeDDGDInvDKappa(), and computeDFDKappa().

double oofem::ConcreteDPM::computeInverseDamage ( double  dam,
GaussPoint gp 
)
double oofem::ConcreteDPM::computeRatioPotential ( double  sig,
double  tempKappa 
)

This function computes the ratio of the volumetric and deviatoric component of the flow direction.

It is used within the vertex return to check, if the vertex return is admissible.

Definition at line 1416 of file concretedpm.C.

References computeHardeningOne(), dilationConst, fc, ft, m, mQ, and nu.

Referenced by performVertexReturn().

double oofem::ConcreteDPM::computeTempKappa ( double  kappaInitial,
double  sigTrial,
double  rhoTrial,
double  sig 
)

Compute temporary kappa.

Definition at line 1112 of file concretedpm.C.

References computeDuctilityMeasure(), gM, kM, and rho.

Referenced by performVertexReturn().

void oofem::ConcreteDPM::computeTrialCoordinates ( const FloatArray stress,
GaussPoint gp 
)
double oofem::ConcreteDPM::computeYieldValue ( double  sig,
double  rho,
double  theta,
double  tempKappa 
) const

Compute the yield value based on stress and hardening variable.

Parameters
sigVolumetric stress.
rhoLength of the deviatoric stress.
thetaLode angle of the stress state.
tempKappaHardening variable.
Returns
Yield value.

Definition at line 1131 of file concretedpm.C.

References computeHardeningOne(), ecc, fc, and m.

Referenced by performPlasticityReturn(), performRegularReturn(), and performVertexReturn().

MaterialStatus * oofem::ConcreteDPM::CreateStatus ( GaussPoint gp) const
protectedvirtual

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.

Definition at line 1870 of file concretedpm.C.

References oofem::FEMComponent::giveDomain().

void oofem::ConcreteDPM::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 1593 of file concretedpm.C.

References oofem::StructuralMaterial::give3dMaterialStiffnessMatrix(), giveLinearElasticMaterial(), giveStatus(), oofem::ConcreteDPMStatus::giveTempDamage(), and oofem::FloatMatrix::times().

Referenced by oofem::ConcreteDPMStatus::restoreConsistency().

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

Reimplemented from oofem::StructuralMaterial.

Definition at line 464 of file concretedpm.h.

virtual const char* oofem::ConcreteDPM::giveInputRecordName ( ) const
inlinevirtual
Returns
Input record name of the receiver.

Implements oofem::FEMComponent.

Definition at line 465 of file concretedpm.h.

References _IFT_ConcreteDPM_Name.

int oofem::ConcreteDPM::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.

Reimplemented from oofem::StructuralMaterial.

Definition at line 1825 of file concretedpm.C.

References oofem::FloatArray::at(), oofem::ConcreteDPMStatus::giveDamage(), oofem::StructuralMaterial::giveIPValue(), oofem::ConcreteDPMStatus::giveKappaD(), oofem::ConcreteDPMStatus::giveKappaP(), oofem::ConcreteDPMStatus::givePlasticStrain(), giveStatus(), oofem::ConcreteDPMStatus::giveTempDamage(), oofem::ConcreteDPMStatus::giveVolumetricPlasticStrain(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().

LinearElasticMaterial* oofem::ConcreteDPM::giveLinearElasticMaterial ( )
inline

Definition at line 470 of file concretedpm.h.

References oofem::IntegrationPointStatus::gp.

Referenced by give3dMaterialStiffnessMatrix().

virtual ConcreteDPMStatus* oofem::ConcreteDPM::giveStatus ( GaussPoint gp) const
inlinevirtual

Returns material status of receiver in given integration point.

If status does not exist yet, it is created using CreateStatus member function.

Parameters
gpReturns reference to material status belonging to integration point gp.
Returns
Material status associated with given integration point.

Reimplemented from oofem::Material.

Definition at line 467 of file concretedpm.h.

References oofem::Material::giveStatus().

Referenced by assignStateFlag(), computeDamage(), computeDamageParam(), computeDuctilityMeasureDamage(), computeEquivalentStrain(), computeInverseDamage(), give3dMaterialStiffnessMatrix(), giveIPValue(), giveRealStressVector_3d(), initDamaged(), performPlasticityReturn(), performRegularReturn(), performVertexReturn(), and setIPValue().

IRResultType oofem::ConcreteDPM::initializeFrom ( InputRecord ir)
virtual
virtual bool oofem::ConcreteDPM::isCharacteristicMtrxSymmetric ( MatResponseMode  rMode)
inlinevirtual

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

Reimplemented from oofem::Material.

Definition at line 727 of file concretedpm.h.

References oofem::ConcreteDPMStatus::setIPValue().

void oofem::ConcreteDPM::performVertexReturn ( FloatArray stress,
double  apexStress,
GaussPoint gp 
)

Perform stress return for vertex case of the plasticity model, i.e.

if the trial stress state lies within the vertex region.

Parameters
stressStress vector of this Gauss point.
apexStressVolumetric stress at the apex of the yield surface.
gpGauss point.

Definition at line 1016 of file concretedpm.C.

References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeRatioPotential(), oofem::StructuralMaterial::computeSecondCoordinate(), computeTempKappa(), computeYieldValue(), effectiveStress, ft, oofem::FloatArray::giveSize(), giveStatus(), oofem::ConcreteDPMStatus::giveTempKappaP(), tempKappaP, vertexType, VT_Compression, VT_Regular, VT_Tension, and yieldTol.

Referenced by performPlasticityReturn().

int oofem::ConcreteDPM::setIPValue ( const FloatArray value,
GaussPoint gp,
InternalStateType  type 
)
virtual

Sets the value of a certain variable at a given integration point to the given value.

Parameters
valueContains the value(s) to be set (in reduced form).
gpIntegration point.
typeDetermines the type of internal variable.
typeDetermines the type of internal variable.
Returns
Nonzero if ok, zero if var not supported.

Reimplemented from oofem::StructuralMaterial.

Definition at line 1814 of file concretedpm.C.

References giveStatus(), oofem::ConcreteDPMStatus::setIPValue(), and oofem::StructuralMaterial::setIPValue().

Member Data Documentation

double oofem::ConcreteDPM::AHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 377 of file concretedpm.h.

Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), and initializeFrom().

double oofem::ConcreteDPM::ASoft
protected

Parameter of the ductilityMeasure of the damage model.

Definition at line 386 of file concretedpm.h.

Referenced by computeDuctilityMeasureDamage(), and initializeFrom().

double oofem::ConcreteDPM::BHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 379 of file concretedpm.h.

Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), and initializeFrom().

double oofem::ConcreteDPM::CHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 381 of file concretedpm.h.

Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), and initializeFrom().

double oofem::ConcreteDPM::damage
protected

Damage variable of damage model.

Definition at line 435 of file concretedpm.h.

Referenced by assignStateFlag(), ConcreteDPM(), and giveRealStressVector_3d().

double oofem::ConcreteDPM::deltaLambda
protected

Plastic multiplier of the plasticity model.

Definition at line 395 of file concretedpm.h.

Referenced by computeAMatrix(), and performRegularReturn().

double oofem::ConcreteDPM::DHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 383 of file concretedpm.h.

Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), and initializeFrom().

double oofem::ConcreteDPM::dilationConst
protected

Control parameter for te volumetric plastic flow of the plastic potential.

Definition at line 392 of file concretedpm.h.

Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDGDInv(), computeRatioPotential(), and initializeFrom().

double oofem::ConcreteDPM::ecc
protected
double oofem::ConcreteDPM::ef
protected

Control parameter for the exponential softening law.

Definition at line 439 of file concretedpm.h.

Referenced by computeDamageParam(), computeInverseDamage(), and initializeFrom().

FloatArray oofem::ConcreteDPM::effectiveStress
protected

Stress and its deviatoric part.

Definition at line 448 of file concretedpm.h.

Referenced by computeTrialCoordinates(), giveRealStressVector_3d(), performPlasticityReturn(), and performVertexReturn().

double oofem::ConcreteDPM::eM
protected

Elastic Young's modulus.

Definition at line 418 of file concretedpm.h.

Referenced by computeDamageParam(), computeInverseDamage(), giveRealStressVector_3d(), initializeFrom(), and performPlasticityReturn().

double oofem::ConcreteDPM::fc
protected

Parameters of the yield surface of the plasticity model.

fc is the uniaxial compressive strength, ft the uniaxial tensile strength and ecc controls the out of roundness of the deviatoric section.

Definition at line 374 of file concretedpm.h.

Referenced by checkForVertexCase(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDDuctilityMeasureDInv(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeDuctilityMeasure(), computeRatioPotential(), computeYieldValue(), and initializeFrom().

double oofem::ConcreteDPM::gM
protected

Elastic shear modulus.

Definition at line 420 of file concretedpm.h.

Referenced by computeAMatrix(), computeTempKappa(), initializeFrom(), and performRegularReturn().

double oofem::ConcreteDPM::helem
protected

Element size (to be used in fracture energy approach (crack band).

Definition at line 412 of file concretedpm.h.

Referenced by computeDamageParam(), computeInverseDamage(), ConcreteDPM(), initDamaged(), and initializeFrom().

double oofem::ConcreteDPM::href
protected

Material parameter of the size-dependent adjustment (reference element size)

Definition at line 453 of file concretedpm.h.

Referenced by computeDamage(), computeDamageParam(), ConcreteDPM(), and initializeFrom().

double oofem::ConcreteDPM::kappaD
protected

Hardening variable of damage model.

Definition at line 431 of file concretedpm.h.

Referenced by ConcreteDPM(), and giveRealStressVector_3d().

double oofem::ConcreteDPM::kappaP
protected

Hardening variable of plasticity model.

Definition at line 427 of file concretedpm.h.

Referenced by assignStateFlag(), computeEquivalentStrain(), ConcreteDPM(), and performRegularReturn().

double oofem::ConcreteDPM::kM
protected

Elastic bulk modulus.

Definition at line 422 of file concretedpm.h.

Referenced by computeAMatrix(), computeTempKappa(), initializeFrom(), and performRegularReturn().

LinearElasticMaterial* oofem::ConcreteDPM::linearElasticMaterial
protected

Pointer for linear elastic material.

Definition at line 415 of file concretedpm.h.

Referenced by ConcreteDPM(), initializeFrom(), and ~ConcreteDPM().

double oofem::ConcreteDPM::m
protected
double oofem::ConcreteDPM::mQ
protected

The dilation parameter of the plastic potential.

Definition at line 409 of file concretedpm.h.

Referenced by computeDDGDInvDKappa(), computeDGDInv(), and computeRatioPotential().

int oofem::ConcreteDPM::newtonIter
protected

Maximum number of iterations for stress return.

Definition at line 445 of file concretedpm.h.

Referenced by checkForVertexCase(), ConcreteDPM(), initializeFrom(), and performRegularReturn().

double oofem::ConcreteDPM::nu
protected

Elastic Poisson's ration.

Definition at line 424 of file concretedpm.h.

Referenced by computeRatioPotential(), giveRealStressVector_3d(), initializeFrom(), and performPlasticityReturn().

double oofem::ConcreteDPM::rho
protected
double oofem::ConcreteDPM::sig
protected

the volumetric stress.

Definition at line 398 of file concretedpm.h.

Referenced by performPlasticityReturn(), and performRegularReturn().

double oofem::ConcreteDPM::tempDamage
protected
double oofem::ConcreteDPM::tempKappaD
protected

Definition at line 432 of file concretedpm.h.

Referenced by computeDamage(), ConcreteDPM(), and giveRealStressVector_3d().

double oofem::ConcreteDPM::tempKappaP
protected
Concrete_VertexType oofem::ConcreteDPM::vertexType
protected

Definition at line 369 of file concretedpm.h.

Referenced by checkForVertexCase(), performPlasticityReturn(), and performVertexReturn().

double oofem::ConcreteDPM::yieldHardInitial
protected

Parameter of the hardening law of the plasticity model.

Definition at line 389 of file concretedpm.h.

Referenced by computeHardeningOne(), computeHardeningOnePrime(), and initializeFrom().

double oofem::ConcreteDPM::yieldTol
protected

Yield tolerance for the plasticity model.

Definition at line 442 of file concretedpm.h.

Referenced by checkForVertexCase(), ConcreteDPM(), initializeFrom(), performPlasticityReturn(), performRegularReturn(), and performVertexReturn().


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:34 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011