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

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

#include <concretedpm2.h>

+ Inheritance diagram for oofem::ConcreteDPM2:
+ Collaboration diagram for oofem::ConcreteDPM2:

Public Member Functions

 ConcreteDPM2 (int n, Domain *d)
 Constructor. More...
 
virtual ~ConcreteDPM2 ()
 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
 
LinearElasticMaterialgiveLinearElasticMaterial ()
 
virtual void giveRealStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &totalStrain, TimeStep *tStep)
 Default implementation relies on giveRealStressVector_StressControl. More...
 
virtual void giveRealStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &strainVector, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More...
 
void performPlasticityReturn (GaussPoint *gp, const FloatMatrix &D, const FloatArray &strain)
 Perform stress return of the plasticity model and compute history variables. More...
 
bool checkForVertexCase (double &answer, double sig, double tempKappa, bool mode1d)
 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...
 
double performRegularReturn (FloatArray &stress, double kappaP, GaussPoint *gp)
 Perform regular stress return for the plasticity model, i.e. More...
 
void compute1dJacobian (FloatMatrix &answer, double totalsigma, double tempKappa, double deltaLambda, GaussPoint *gp)
 Compute jacobian for 1D case. More...
 
void computeJacobian (FloatMatrix &answer, double sig, double rho, double tempKappa, double deltaLambda, GaussPoint *gp)
 Compute jacobian for 2D(plane strain) and 3d cases. More...
 
double performVertexReturn (FloatArray &stress, double apexStress, double tempKappaP, 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 computeHardeningTwo (double tempKappa) const
 Compute the value of the hardening function based on the hardening variable. More...
 
double computeHardeningTwoPrime (double tempKappa) const
 Compute the derivative of the hardening function based on the hardening parameter. More...
 
double computeDFDKappa (double sig, double rho, double tempKappa, bool mode1d)
 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...
 
double computeDKappaDDeltaLambda1d (double sig, double tempKappa)
 
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)
 Compute derivative the ductility measure with respect to the stress state. More...
 
double computeDDuctilityMeasureDInv1d (double sigma, double tempKappa)
 Compute derivative the ductility measure with respect to the stress state. More...
 
void computeDGDInv (FloatArray &answer, double sig, double rho, double tempKappa)
 Compute derivative the palstic potential function with respect to the stress state. More...
 
double computeDGDInv1d (double sig, double tempKappa)
 Compute derivative the palstic potential function with respect to the stress state. More...
 
double computeRatioPotential (double sig, double tempKappa)
 This function computes the ratio of the volumetric and deviatoric component of the flow direction. More...
 
double computeRateFactor (double alpha, double timeFactor, GaussPoint *gp, TimeStep *deltaTime)
 This function computes the rate factor which is used to take into account the strain rate dependence of the material. 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...
 
double computeDDGDDInv1d (double sigma, 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...
 
double computeDDGDInvDKappa1d (double sigma, double tempKappa)
 
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 computeDDKappaDDeltaLambdaDInv1d (double sigma, double tempKappa)
 
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...
 
double computeDDKappaDDeltaLambdaDKappa1d (double sig, double tempKappa)
 
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 computeDFDInv1d (double sigma, double tempKappa) const
 
double computeTempKappa (double kappaInitial, double sigTrial, double rhoTrial, double sig)
 Compute tempKappa. More...
 
void computeDamage (FloatArray &answer, const FloatArray &strain, const FloatMatrix &D, double timeFactor, GaussPoint *gp, TimeStep *tStep, double alpha)
 Compute damage parameters. More...
 
int checkForUnAndReloading (double &tempEquivStrain, double &minEquivStrain, const FloatMatrix &D, GaussPoint *gp)
 Check for un- and reloading in the damage part. More...
 
double computeAlpha (FloatArray &effectiveStressTension, FloatArray &effectiveStressCompression, FloatArray &effectiveStress)
 
virtual double computeDamageParamTension (double equivStrain, double kappaOne, double kappaTwo, double le, double omegaOld)
 Compute damage parameter. More...
 
virtual double computeDamageParamCompression (double equivStrain, double kappaOne, double kappaTwo, double omegaOld)
 
double computeDeltaPlasticStrainNormTension (double tempKappaD, double kappaD, GaussPoint *gp)
 Compute equivalent strain value. More...
 
double computeDeltaPlasticStrainNormCompression (double tempAlpha, double tempKappaD, double kappaD, GaussPoint *gp)
 
virtual double computeEquivalentStrain (double sig, double rho, double theta)
 
double computeDuctilityMeasureDamage (const FloatArray &strain, GaussPoint *gp)
 Compute the ductility measure for the damage model. More...
 
void initDamaged (double kappa, const FloatArray &strain, GaussPoint *gp)
 Initialize the characteristic length, if damage is not yet activated Set the increase factor for the strain rate dependence. More...
 
void computeTrialCoordinates (const FloatArray &stress, double &sig, double &rho, double &theta)
 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 seconfd 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 give1dStressStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d stiffness matrix of receiver. 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...
 
void compute3dSecantStiffness (FloatMatrix &answer, GaussPoint *gp, TimeStep *tStep)
 
virtual bool isCharacteristicMtrxSymmetric (MatResponseMode rMode)
 Returns true if stiffness matrix of receiver is symmetric Default implementation returns true. 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_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 int setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type)
 Sets the value of a certain variable at a given integration point to the given value. More...
 
virtual void give2dBeamLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d beam layer stiffness matrix of receiver. More...
 
virtual void givePlateLayerStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d plate layer stiffness matrix of receiver. More...
 
virtual void giveFiberStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 1d fiber stiffness matrix of receiver. More...
 
virtual void give2dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 2d lattice stiffness matrix of receiver. More...
 
virtual void give3dLatticeStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 3d lattice stiffness matrix of receiver. More...
 
virtual void give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of plate subsoil model. More...
 
virtual void give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing stiffness matrix of beam3d subsoil model. More...
 
virtual void giveFirstPKStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. More...
 
virtual void giveFirstPKStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveFirstPKStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveFirstPKStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 Default implementation relies on giveFirstPKStressVector_3d. More...
 
virtual void giveCauchyStressVector_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void giveCauchyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void giveCauchyStressVector_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void giveCauchyStressVector_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep)
 
virtual void givePlaneStressStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane stress stiffness matrix of receiver. More...
 
virtual void givePlaneStressStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStrainStiffMtrx (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 Method for computing plane strain stiffness matrix of receiver. More...
 
virtual void givePlaneStrainStiffMtrx_dPdF (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep)
 
virtual void give1dStressStiffMtrx_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 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...
 
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  ConcreteDPM2_ReturnType { RT_Regular, RT_Tension, RT_Compression, RT_Auxiliary }
 
enum  ConcreteDPM2_ReturnResult { RR_NotConverged, RR_Converged }
 

Protected Member Functions

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

Protected Attributes

ConcreteDPM2_ReturnType returnType
 
ConcreteDPM2_ReturnResult returnResult
 
double fc
 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. More...
 
double ft
 
double ecc
 
int isotropicFlag
 
double e0
 
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 hardeningModulus
 Hardening modulus. More...
 
double ASoft
 Parameter of the ductilityMeasure of the damage model. More...
 
double yieldHardPrimePeak
 Parameter of the hardening law of the plasticity 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 sig
 Volumetric stress. More...
 
double rho
 Length of the deviatoric stress. More...
 
double thetaTrial
 Lode angle of the trial stress.. More...
 
double m
 Friction parameter of the yield surface. More...
 
double mQ
 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 efCompression
 Control parameter for the exponential softening law. More...
 
double wf
 Control parameter for the linear/bilinear softening law in tension. More...
 
double wfOne
 Control parameter for the bilinear softening law in tension. More...
 
double ftOne
 Control parameter for the bilinear softening law. More...
 
double yieldTol
 yield tolerance for the plasticity model. More...
 
double yieldTolDamage
 yield tolerance for the damage model. More...
 
int newtonIter
 Maximum number of iterations for stress return. More...
 
int softeningType
 Type of softening function used. More...
 
double timeFactor
 Input parameter which simulates a loading rate. Only for debugging purposes. More...
 
double fcZero
 This parameter is needed for the rate dependence. It should be read in if rate dependence is considered. More...
 
int strainRateFlag
 Flag which signals if strainRate effects should be considered. 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.

This is an extension of concretedpm2. 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 the effective stress. The damage parameter of the isotropic damage model is based on the total volumetric strain. An exponential softening law is implemented.

Definition at line 588 of file concretedpm2.h.

Member Enumeration Documentation

Enumerator
RR_NotConverged 
RR_Converged 

Definition at line 596 of file concretedpm2.h.

Enumerator
RT_Regular 
RT_Tension 
RT_Compression 
RT_Auxiliary 

Definition at line 593 of file concretedpm2.h.

Constructor & Destructor Documentation

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

Constructor.

Definition at line 503 of file concretedpm2.C.

References linearElasticMaterial, newtonIter, yieldTol, and yieldTolDamage.

oofem::ConcreteDPM2::~ConcreteDPM2 ( )
virtual

Destructor.

Definition at line 513 of file concretedpm2.C.

References linearElasticMaterial.

Member Function Documentation

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

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
answerVolumetric apex stress.
sigVolumetric stress.
tempKappaHardening variable.

Definition at line 1520 of file concretedpm2.C.

References returnType, RT_Compression, RT_Regular, and RT_Tension.

Referenced by performPlasticityReturn().

void oofem::ConcreteDPM2::compute1dJacobian ( FloatMatrix answer,
double  totalsigma,
double  tempKappa,
double  deltaLambda,
GaussPoint gp 
)

Compute jacobian for 1D case.

Parameters
totalsigmastress value
tempKappaplastic strain
deltaLambdaplastic multiplier
gpGauss point

Definition at line 1907 of file concretedpm2.C.

References oofem::FloatMatrix::at(), computeDDGDDInv1d(), computeDDGDInvDKappa1d(), computeDDKappaDDeltaLambdaDInv1d(), computeDDKappaDDeltaLambdaDKappa1d(), computeDFDInv1d(), computeDFDKappa(), computeDGDInv1d(), computeDKappaDDeltaLambda(), eM, and oofem::FloatMatrix::resize().

Referenced by performRegularReturn().

double oofem::ConcreteDPM2::computeAlpha ( FloatArray effectiveStressTension,
FloatArray effectiveStressCompression,
FloatArray effectiveStress 
)
void oofem::ConcreteDPM2::computeDamage ( FloatArray answer,
const FloatArray strain,
const FloatMatrix D,
double  timeFactor,
GaussPoint gp,
TimeStep tStep,
double  alpha 
)

Compute damage parameters.

Definition at line 794 of file concretedpm2.C.

References oofem::FloatArray::at(), checkForUnAndReloading(), computeDamageParamCompression(), computeDamageParamTension(), computeDeltaPlasticStrainNormCompression(), computeDeltaPlasticStrainNormTension(), computeDuctilityMeasureDamage(), computeRateFactor(), e0, oofem::ConcreteDPM2Status::giveAlpha(), oofem::ConcreteDPM2Status::giveDamageCompression(), oofem::ConcreteDPM2Status::giveDamageTension(), oofem::ConcreteDPM2Status::giveEquivStrain(), oofem::ConcreteDPM2Status::giveEquivStrainCompression(), oofem::ConcreteDPM2Status::giveEquivStrainTension(), oofem::ConcreteDPM2Status::giveKappaDCompression(), oofem::ConcreteDPM2Status::giveKappaDCompressionOne(), oofem::ConcreteDPM2Status::giveKappaDCompressionTwo(), oofem::ConcreteDPM2Status::giveKappaDTension(), oofem::ConcreteDPM2Status::giveKappaDTensionOne(), oofem::ConcreteDPM2Status::giveKappaDTensionTwo(), oofem::ConcreteDPM2Status::giveLe(), oofem::ConcreteDPM2Status::giveRateFactor(), oofem::Material::giveStatus(), initDamaged(), oofem::TimeStep::isTheFirstStep(), oofem::ConcreteDPM2Status::letTempDamageCompressionBe(), oofem::ConcreteDPM2Status::letTempDamageTensionBe(), oofem::ConcreteDPM2Status::letTempEquivStrainBe(), oofem::ConcreteDPM2Status::letTempEquivStrainCompressionBe(), oofem::ConcreteDPM2Status::letTempEquivStrainTensionBe(), oofem::ConcreteDPM2Status::letTempKappaDCompressionBe(), oofem::ConcreteDPM2Status::letTempKappaDCompressionOneBe(), oofem::ConcreteDPM2Status::letTempKappaDCompressionTwoBe(), oofem::ConcreteDPM2Status::letTempKappaDTensionBe(), oofem::ConcreteDPM2Status::letTempKappaDTensionOneBe(), oofem::ConcreteDPM2Status::letTempKappaDTensionTwoBe(), oofem::ConcreteDPM2Status::letTempRateFactorBe(), oofem::FloatArray::resize(), and yieldTolDamage.

Referenced by giveRealStressVector_1d(), and giveRealStressVector_3d().

double oofem::ConcreteDPM2::computeDamageParamCompression ( double  equivStrain,
double  kappaOne,
double  kappaTwo,
double  omegaOld 
)
virtual

Definition at line 1278 of file concretedpm2.C.

References e0, efCompression, eM, ft, isotropicFlag, newtonIter, OOFEM_ERROR, and yieldTolDamage.

Referenced by computeDamage().

double oofem::ConcreteDPM2::computeDamageParamTension ( double  equivStrain,
double  kappaOne,
double  kappaTwo,
double  le,
double  omegaOld 
)
virtual

Compute damage parameter.

Definition at line 1215 of file concretedpm2.C.

References e0, eM, ft, ftOne, newtonIter, OOFEM_ERROR, softeningType, wf, wfOne, and yieldTolDamage.

Referenced by computeDamage().

void oofem::ConcreteDPM2::computeDCosThetaDStress ( FloatArray answer,
const FloatArray stress 
) const

Computes the derivative of costheta with respect to the stress.

void oofem::ConcreteDPM2::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 2483 of file concretedpm2.C.

References computeHardeningOne(), computeHardeningTwo(), dilationConst, fc, ft, m, and oofem::FloatMatrix::resize().

Referenced by computeDDKappaDDeltaLambdaDInv(), and computeJacobian().

double oofem::ConcreteDPM2::computeDDGDDInv1d ( double  sigma,
double  tempKappa 
)

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

Definition at line 2467 of file concretedpm2.C.

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

Referenced by compute1dJacobian(), and computeDDKappaDDeltaLambdaDInv1d().

void oofem::ConcreteDPM2::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 2404 of file concretedpm2.C.

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

Referenced by computeDDKappaDDeltaLambdaDKappa(), and computeJacobian().

double oofem::ConcreteDPM2::computeDDGDInvDKappa1d ( double  sigma,
double  tempKappa 
)
void oofem::ConcreteDPM2::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 2183 of file concretedpm2.C.

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

Referenced by computeJacobian().

double oofem::ConcreteDPM2::computeDDKappaDDeltaLambdaDInv1d ( double  sigma,
double  tempKappa 
)
double oofem::ConcreteDPM2::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 2222 of file concretedpm2.C.

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

Referenced by computeJacobian().

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

Compute derivative the ductility measure with respect to the stress state.

Parameters
answerarray of the derivative of the ductility measure with respect to volumetric and deviatoric stress
sigVolumetric stress.
rhoLength of the deviatoric strength.
tempKappaplastic strain

Definition at line 2293 of file concretedpm2.C.

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

Referenced by computeDDKappaDDeltaLambdaDInv().

double oofem::ConcreteDPM2::computeDDuctilityMeasureDInv1d ( double  sigma,
double  tempKappa 
)

Compute derivative the ductility measure with respect to the stress state.

Returns
The derivative of the ductility measure with respect to stress
Parameters
sigstress.
tempKappaplastic strain

Definition at line 2272 of file concretedpm2.C.

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

Referenced by computeDDKappaDDeltaLambdaDInv1d().

double oofem::ConcreteDPM2::computeDeltaPlasticStrainNormTension ( double  tempKappaD,
double  kappaD,
GaussPoint gp 
)
void oofem::ConcreteDPM2::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 2097 of file concretedpm2.C.

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

Referenced by computeJacobian().

double oofem::ConcreteDPM2::computeDFDInv1d ( double  sigma,
double  tempKappa 
) const

Definition at line 2082 of file concretedpm2.C.

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

Referenced by compute1dJacobian().

double oofem::ConcreteDPM2::computeDFDKappa ( double  sig,
double  rho,
double  tempKappa,
bool  mode1d 
)

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
Derivative of the yield surface.

Definition at line 2028 of file concretedpm2.C.

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

Referenced by compute1dJacobian(), and computeJacobian().

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

Compute derivative the palstic potential function with respect to the stress state.

Parameters
answerarray of the derivative of the plastic potential with respect to volumetric and deviatoric stress
sigvolumetric stress.
rhodeviatoric stress.
tempKappaplastic strain

Definition at line 2334 of file concretedpm2.C.

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

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

double oofem::ConcreteDPM2::computeDGDInv1d ( double  sig,
double  tempKappa 
)

Compute derivative the palstic potential function with respect to the stress state.

Parameters
returnThe derivative of the plastic potential with respect to stress
sigstress.
tempKappaplastic strain

Definition at line 2319 of file concretedpm2.C.

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

Referenced by compute1dJacobian(), computeDDKappaDDeltaLambdaDInv1d(), computeDDKappaDDeltaLambdaDKappa1d(), computeDKappaDDeltaLambda1d(), and performRegularReturn().

double oofem::ConcreteDPM2::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 2139 of file concretedpm2.C.

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

Referenced by compute1dJacobian(), computeJacobian(), and performRegularReturn().

double oofem::ConcreteDPM2::computeDKappaDDeltaLambda1d ( double  sig,
double  tempKappa 
)

Definition at line 2130 of file concretedpm2.C.

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

Referenced by performRegularReturn().

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

Compute the derivative of R with respect to costheta.

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

Computes the derivative of rho with respect to the stress.

Definition at line 2789 of file concretedpm2.C.

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

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

Computes the derivative of sig with respect to the stress.

Definition at line 2811 of file concretedpm2.C.

double oofem::ConcreteDPM2::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 1668 of file concretedpm2.C.

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

Referenced by computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDInv1d(), computeDDKappaDDeltaLambdaDKappa(), computeDDKappaDDeltaLambdaDKappa1d(), computeDKappaDDeltaLambda(), computeDKappaDDeltaLambda1d(), and computeTempKappa().

double oofem::ConcreteDPM2::computeDuctilityMeasureDamage ( const FloatArray strain,
GaussPoint gp 
)

Compute the ductility measure for the damage model.

Definition at line 1369 of file concretedpm2.C.

References ASoft, rho, and sig.

Referenced by computeDamage().

double oofem::ConcreteDPM2::computeEquivalentStrain ( double  sig,
double  rho,
double  theta 
)
virtual

Definition at line 1194 of file concretedpm2.C.

References e0, ecc, fc, and m.

Referenced by checkForUnAndReloading().

double oofem::ConcreteDPM2::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 2587 of file concretedpm2.C.

References yieldHardInitial, and yieldHardPrimePeak.

Referenced by computeDDGDDInv(), computeDDGDDInv1d(), computeDDGDInvDKappa(), computeDDGDInvDKappa1d(), computeDFDInv(), computeDFDInv1d(), computeDFDKappa(), computeDGDInv(), computeDGDInv1d(), computeRatioPotential(), and computeYieldValue().

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

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

Parameters
tempKappaHardening variable.
Returns
Derivative of the hardening function.

Definition at line 2604 of file concretedpm2.C.

References yieldHardInitial, and yieldHardPrimePeak.

Referenced by computeDDGDInvDKappa(), computeDDGDInvDKappa1d(), and computeDFDKappa().

double oofem::ConcreteDPM2::computeHardeningTwo ( 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 2620 of file concretedpm2.C.

References yieldHardPrimePeak.

Referenced by computeDDGDDInv(), computeDDGDDInv1d(), computeDDGDInvDKappa(), computeDDGDInvDKappa1d(), computeDeltaPlasticStrainNormCompression(), computeDFDInv(), computeDFDInv1d(), computeDFDKappa(), computeDGDInv(), computeDGDInv1d(), computeRatioPotential(), and computeYieldValue().

double oofem::ConcreteDPM2::computeHardeningTwoPrime ( double  tempKappa) const

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

Parameters
tempKappaHardening variable.
Returns
Derivative of the hardening function.

Definition at line 2633 of file concretedpm2.C.

References yieldHardPrimePeak.

Referenced by computeDDGDInvDKappa(), computeDDGDInvDKappa1d(), and computeDFDKappa().

void oofem::ConcreteDPM2::computeJacobian ( FloatMatrix answer,
double  sig,
double  rho,
double  tempKappa,
double  deltaLambda,
GaussPoint gp 
)

Compute jacobian for 2D(plane strain) and 3d cases.

Parameters
sigvolumetric strain
rhodeviatoric
tempKappaplastic strain
deltaLambdaplastic multiplier
gpGauss point

Definition at line 1940 of file concretedpm2.C.

References oofem::FloatArray::at(), oofem::FloatMatrix::at(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeDKappaDDeltaLambda(), gM, kM, and oofem::FloatMatrix::resize().

Referenced by performRegularReturn().

double oofem::ConcreteDPM2::computeRateFactor ( double  alpha,
double  timeFactor,
GaussPoint gp,
TimeStep deltaTime 
)
double oofem::ConcreteDPM2::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 1690 of file concretedpm2.C.

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

Referenced by performVertexReturn().

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

Compute tempKappa.

Definition at line 1649 of file concretedpm2.C.

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

Referenced by performVertexReturn().

void oofem::ConcreteDPM2::computeTrialCoordinates ( const FloatArray stress,
double &  sig,
double &  rho,
double &  theta 
)
double oofem::ConcreteDPM2::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 1999 of file concretedpm2.C.

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

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

MaterialStatus * oofem::ConcreteDPM2::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 2915 of file concretedpm2.C.

References oofem::FEMComponent::giveDomain().

void oofem::ConcreteDPM2::give1dStressStiffMtrx ( FloatMatrix answer,
MatResponseMode  mmode,
GaussPoint gp,
TimeStep tStep 
)
virtual

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 2645 of file concretedpm2.C.

References oofem::FloatArray::at(), oofem::FloatMatrix::at(), eM, oofem::Material::giveStatus(), oofem::ConcreteDPM2Status::giveTempDamageCompression(), oofem::ConcreteDPM2Status::giveTempDamageTension(), oofem::StructuralMaterialStatus::giveTempStrainVector(), isotropicFlag, OOFEM_WARNING, oofem::FloatMatrix::resize(), and oofem::FloatMatrix::times().

void oofem::ConcreteDPM2::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 2675 of file concretedpm2.C.

References compute3dSecantStiffness(), oofem::StructuralMaterial::give3dMaterialStiffnessMatrix(), giveLinearElasticMaterial(), and OOFEM_WARNING.

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

Reimplemented from oofem::StructuralMaterial.

Definition at line 704 of file concretedpm2.h.

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

Implements oofem::FEMComponent.

Definition at line 705 of file concretedpm2.h.

References _IFT_ConcreteDPM2_Name.

int oofem::ConcreteDPM2::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 2872 of file concretedpm2.C.

References oofem::FloatArray::at(), oofem::ConcreteDPM2Status::giveDamageCompression(), oofem::ConcreteDPM2Status::giveDamageTension(), oofem::ConcreteDPM2Status::giveDissWork(), oofem::StructuralMaterial::giveIPValue(), oofem::ConcreteDPM2Status::givePlasticStrain(), oofem::Material::giveStatus(), oofem::ConcreteDPM2Status::giveStressWork(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().

IRResultType oofem::ConcreteDPM2::initializeFrom ( InputRecord ir)
virtual

Initializes receiver according to object description stored in input record.

This function is called immediately after creating object using constructor. Input record can be imagined as data record in component database belonging to receiver. Receiver may use value-name extracting functions to extract particular field from record.

See also
IR_GIVE_FIELD
IR_GIVE_OPTIONAL_FIELD
Parameters
irInput record to initialize from.
Returns
IRResultType

Reimplemented from oofem::StructuralMaterial.

Definition at line 519 of file concretedpm2.C.

References _IFT_ConcreteDPM2_ahard, _IFT_ConcreteDPM2_asoft, _IFT_ConcreteDPM2_bhard, _IFT_ConcreteDPM2_chard, _IFT_ConcreteDPM2_dhard, _IFT_ConcreteDPM2_dilation, _IFT_ConcreteDPM2_ecc, _IFT_ConcreteDPM2_efc, _IFT_ConcreteDPM2_fc, _IFT_ConcreteDPM2_fcZero, _IFT_ConcreteDPM2_ft, _IFT_ConcreteDPM2_ftOne, _IFT_ConcreteDPM2_helem, _IFT_ConcreteDPM2_hp, _IFT_ConcreteDPM2_isoflag, _IFT_ConcreteDPM2_kinit, _IFT_ConcreteDPM2_newtoniter, _IFT_ConcreteDPM2_rateFlag, _IFT_ConcreteDPM2_softeningType, _IFT_ConcreteDPM2_timeFactor, _IFT_ConcreteDPM2_wf, _IFT_ConcreteDPM2_wfOne, _IFT_ConcreteDPM2_yieldtol, _IFT_IsotropicLinearElasticMaterial_e, _IFT_IsotropicLinearElasticMaterial_n, _IFT_IsotropicLinearElasticMaterial_talpha, oofem::Dictionary::add(), AHard, ASoft, BHard, CHard, DHard, dilationConst, e0, ecc, efCompression, eM, fc, fcZero, ft, ftOne, gM, helem, oofem::LinearElasticMaterial::initializeFrom(), oofem::StructuralMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, oofem::IRRT_BAD_FORMAT, oofem::IRRT_OK, isotropicFlag, kM, linearElasticMaterial, m, newtonIter, nu, OOFEM_WARNING, oofem::Material::propertyDictionary, softeningType, strainRateFlag, tAlpha, timeFactor, wf, wfOne, yieldHardInitial, yieldHardPrimePeak, yieldTol, and yieldTolDamage.

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

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

Reimplemented from oofem::Material.

Definition at line 1053 of file concretedpm2.h.

double oofem::ConcreteDPM2::performVertexReturn ( FloatArray stress,
double  apexStress,
double  tempKappaP,
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.
tempKappaPtemporary cummulative plastic strain
gpGauss point.
Returns
updated temporary cummulative plastic strain

Definition at line 1544 of file concretedpm2.C.

References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeRatioPotential(), oofem::StructuralMaterial::computeSecondCoordinate(), computeTempKappa(), computeYieldValue(), oofem::FloatArray::giveSize(), returnResult, returnType, RR_Converged, RR_NotConverged, RT_Compression, RT_Regular, RT_Tension, and yieldTol.

Referenced by performPlasticityReturn().

Member Data Documentation

double oofem::ConcreteDPM2::AHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 607 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::ASoft
protected

Parameter of the ductilityMeasure of the damage model.

Definition at line 619 of file concretedpm2.h.

Referenced by computeDuctilityMeasureDamage(), and initializeFrom().

double oofem::ConcreteDPM2::BHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 609 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::CHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 611 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::DHard
protected

Parameter of the ductilityMeasure of the plasticity model.

Definition at line 613 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::dilationConst
protected

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

Definition at line 628 of file concretedpm2.h.

Referenced by computeDDGDDInv(), computeDDGDDInv1d(), computeDDGDInvDKappa(), computeDDGDInvDKappa1d(), computeDeltaPlasticStrainNormCompression(), computeDGDInv(), computeDGDInv1d(), computeRatioPotential(), and initializeFrom().

double oofem::ConcreteDPM2::ecc
protected
double oofem::ConcreteDPM2::efCompression
protected

Control parameter for the exponential softening law.

Definition at line 664 of file concretedpm2.h.

Referenced by computeDamageParamCompression(), and initializeFrom().

double oofem::ConcreteDPM2::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 600 of file concretedpm2.h.

Referenced by computeDDGDDInv(), computeDDGDDInv1d(), computeDDGDInvDKappa(), computeDDGDInvDKappa1d(), computeDDuctilityMeasureDInv(), computeDDuctilityMeasureDInv1d(), computeDFDInv(), computeDFDInv1d(), computeDFDKappa(), computeDGDInv(), computeDGDInv1d(), computeDuctilityMeasure(), computeEquivalentStrain(), computeRateFactor(), computeRatioPotential(), computeYieldValue(), and initializeFrom().

double oofem::ConcreteDPM2::fcZero
protected

This parameter is needed for the rate dependence. It should be read in if rate dependence is considered.

Definition at line 691 of file concretedpm2.h.

Referenced by computeRateFactor(), and initializeFrom().

double oofem::ConcreteDPM2::ftOne
protected

Control parameter for the bilinear softening law.

Definition at line 673 of file concretedpm2.h.

Referenced by computeDamageParamTension(), and initializeFrom().

double oofem::ConcreteDPM2::gM
protected

Elastic shear modulus.

Definition at line 657 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::hardeningModulus
protected

Hardening modulus.

Definition at line 616 of file concretedpm2.h.

double oofem::ConcreteDPM2::helem
protected

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

Definition at line 649 of file concretedpm2.h.

Referenced by initDamaged(), and initializeFrom().

int oofem::ConcreteDPM2::isotropicFlag
protected
double oofem::ConcreteDPM2::kM
protected

Elastic bulk modulus.

Definition at line 659 of file concretedpm2.h.

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

LinearElasticMaterial* oofem::ConcreteDPM2::linearElasticMaterial
protected

Pointer for linear elastic material.

Definition at line 652 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::mQ
protected

Dilation parameter of the plastic potential.

Definition at line 646 of file concretedpm2.h.

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

int oofem::ConcreteDPM2::newtonIter
protected

Maximum number of iterations for stress return.

Definition at line 682 of file concretedpm2.h.

Referenced by computeDamageParamCompression(), computeDamageParamTension(), ConcreteDPM2(), initializeFrom(), and performRegularReturn().

double oofem::ConcreteDPM2::nu
protected

Elastic poisson's ration.

Definition at line 661 of file concretedpm2.h.

Referenced by computeRatioPotential(), and initializeFrom().

ConcreteDPM2_ReturnResult oofem::ConcreteDPM2::returnResult
protected
ConcreteDPM2_ReturnType oofem::ConcreteDPM2::returnType
protected
double oofem::ConcreteDPM2::sig
protected

Volumetric stress.

Todo:
These values should not be stored by the material model itself. As temp values in the status they would be OK, but this will be thread unsafe (and makes the model into a complete spagetti code) / Mikael

Definition at line 633 of file concretedpm2.h.

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

int oofem::ConcreteDPM2::softeningType
protected

Type of softening function used.

Definition at line 685 of file concretedpm2.h.

Referenced by computeDamageParamTension(), and initializeFrom().

int oofem::ConcreteDPM2::strainRateFlag
protected

Flag which signals if strainRate effects should be considered.

Definition at line 694 of file concretedpm2.h.

Referenced by computeRateFactor(), giveRealStressVector_1d(), giveRealStressVector_3d(), and initializeFrom().

double oofem::ConcreteDPM2::timeFactor
protected

Input parameter which simulates a loading rate. Only for debugging purposes.

Definition at line 688 of file concretedpm2.h.

Referenced by giveRealStressVector_1d(), giveRealStressVector_3d(), and initializeFrom().

double oofem::ConcreteDPM2::wf
protected

Control parameter for the linear/bilinear softening law in tension.

Definition at line 667 of file concretedpm2.h.

Referenced by computeDamageParamTension(), and initializeFrom().

double oofem::ConcreteDPM2::wfOne
protected

Control parameter for the bilinear softening law in tension.

Definition at line 670 of file concretedpm2.h.

Referenced by computeDamageParamTension(), and initializeFrom().

double oofem::ConcreteDPM2::yieldHardInitial
protected

Parameter of the hardening law of the plasticity model.

Definition at line 625 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::yieldHardPrimePeak
protected

Parameter of the hardening law of the plasticity model.

Definition at line 622 of file concretedpm2.h.

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

double oofem::ConcreteDPM2::yieldTol
protected

yield tolerance for the plasticity model.

Definition at line 676 of file concretedpm2.h.

Referenced by ConcreteDPM2(), initializeFrom(), performRegularReturn(), and performVertexReturn().

double oofem::ConcreteDPM2::yieldTolDamage
protected

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