|
OOFEM 3.0
|
#include <concretedpm2.h>
Public Types | |
| enum | ConcreteDPM2_ReturnResult { RR_Unknown , RR_NotConverged , RR_Converged } |
| enum | ConcreteDPM2_ReturnType { RT_Unknown , RT_Regular , RT_Tension , RT_Compression , RT_Auxiliary } |
Public Member Functions | |
| ConcreteDPM2 (int n, Domain *d) | |
| Constructor. | |
| void | initializeFrom (InputRecord &ir) override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| FloatArrayF< 6 > | giveRealStressVector_3d (const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| FloatArrayF< 6 > | performPlasticityReturn (GaussPoint *gp, const FloatMatrixF< 6, 6 > &D, const FloatArrayF< 6 > &strain) const |
| void | checkForVertexCase (double &answer, ConcreteDPM2_ReturnType &returnType, double sig, double tempKappa, GaussPoint *gp) const |
| double | performRegularReturn (FloatArrayF< 6 > &stress, ConcreteDPM2_ReturnResult &returnResult, ConcreteDPM2_ReturnType &returnType, double kappaP, GaussPoint *gp, double theta) const |
| FloatMatrixF< 4, 4 > | computeJacobian (double sig, double rho, double theta, double tempKappa, double deltaLambda, GaussPoint *gp) const |
| double | performVertexReturn (FloatArrayF< 6 > &stress, ConcreteDPM2_ReturnResult &returnResult, ConcreteDPM2_ReturnType &returnType, double apexStress, double tempKappaP, GaussPoint *gp) const |
| double | computeYieldValue (double sig, double rho, double theta, double tempKappa) const |
| double | computeHardeningOne (double tempKappa) const |
| double | computeHardeningOnePrime (double tempKappa) const |
| double | computeHardeningTwo (double tempKappa) const |
| double | computeHardeningTwoPrime (double tempKappa) const |
| double | computeDFDKappa (double sig, double rho, double theta, double tempKappa) const |
| double | computeDKappaDDeltaLambda (double sig, double rho, double theta, double tempKappa) const |
| virtual double | computeDuctilityMeasure (double sig, double rho, double theta) const |
| FloatArrayF< 2 > | computeDDuctilityMeasureDInv (double sig, double rho, double theta, double tempKappa) const |
| FloatArrayF< 2 > | computeDGDInv (double sig, double rho, double tempKappa) const |
| double | computeRatioPotential (double sig, double rho, double tempKappa) const |
| double | computeRateFactor (double alpha, double timeFactor, GaussPoint *gp, TimeStep *deltaTime) const |
| FloatMatrixF< 2, 2 > | computeDDGDDInv (double sig, double rho, double tempKappa) const |
| FloatMatrixF< 6, 6 > | computeDDGDDStress (const FloatArrayF< 6 > &stress, const double tempKappa) const |
| FloatArrayF< 6 > | computeDCosThetaDStress (const FloatArrayF< 6 > &stress) const |
| FloatMatrixF< 6, 6 > | computeDDCosThetaDDStress (const FloatArrayF< 6 > &stress) const |
| FloatArrayF< 2 > | computeDDGDInvDKappa (double sig, double rho, double tempKappa) const |
| FloatArrayF< 2 > | computeDDKappaDDeltaLambdaDInv (double sig, double rho, double theta, double tempKappa) const |
| FloatArrayF< 6 > | computeDDKappaDDeltaLambdaDStress (const FloatArrayF< 6 > &stress, double tempKappa) const |
| FloatArrayF< 6 > | computeDDGDStressDKappa (const FloatArrayF< 6 > &stress, double tempKappa) const |
| double | computeDDKappaDDeltaLambdaDKappa (double sig, double rho, double theta, double tempKappa) const |
| FloatArrayF< 2 > | computeDFDInv (double sig, double rho, double theta, double tempKappa) const |
| FloatArrayF< 6 > | computeDFDStress (const FloatArrayF< 6 > &stress, const double tempKappa) const |
| FloatArrayF< 6 > | computeDGDStress (const FloatArrayF< 6 > &stress, const double tempKappa) const |
| FloatMatrixF< 8, 8 > | computeFullJacobian (const FloatArrayF< 6 > &stress, const double deltaLambda, GaussPoint *gp, TimeStep *atTime, const double tempKappa) const |
| double | computeTempKappa (double kappaInitial, double sigTrial, double rhoTrial, double sig) const |
| Compute tempKappa. | |
| FloatArrayF< 2 > | computeDamage (const FloatArrayF< 6 > &strain, const FloatMatrixF< 6, 6 > &D, double timeFactor, GaussPoint *gp, TimeStep *tStep, double alpha, const FloatArrayF< 6 > &effectiveStress) const |
| Compute damage parameters. | |
| int | checkForUnAndReloading (double &tempEquivStrain, double &minEquivStrain, const FloatMatrixF< 6, 6 > &D, GaussPoint *gp) const |
| Check for un- and reloading in the damage part. | |
| double | computeAlpha (FloatArrayF< 6 > &effectiveStressTension, FloatArrayF< 6 > &effectiveStressCompression, const FloatArrayF< 6 > &effectiveStress) const |
| Compute alpha for rate effect. | |
| virtual double | computeDamageParamTension (double equivStrain, double kappaOne, double kappaTwo, double le, double omegaOld, double rateFactor) const |
| Compute damage parameter in tension. | |
| virtual double | computeDamageParamCompression (double equivStrain, double kappaOne, double kappaTwo, double omegaOld, double rateFactor) const |
| Compute damage parameter in compression. | |
| double | computeDeltaPlasticStrainNormTension (double tempKappaD, double kappaD, GaussPoint *gp) const |
| Compute equivalent strain value for tension. | |
| double | computeDeltaPlasticStrainNormCompression (double tempAlpha, double tempKappaD, double kappaD, GaussPoint *gp, const double rho) const |
| Compute equivalent strain value for compression. | |
| virtual double | computeEquivalentStrain (double sig, double rho, double theta) const |
| Compute the base equivalent strain value. | |
| double | computeDuctilityMeasureDamage (GaussPoint *gp, const double sig, const double rho) const |
| Compute the ductility measure for the damage model. | |
| void | initDamaged (double kappa, const FloatArrayF< 6 > &strain, GaussPoint *gp) const |
| void | computeCoordinates (const FloatArrayF< 6 > &stress, double &sig, double &rho, double &theta) const |
| Compute the Haigh-Westergaard coordinates. | |
| void | assignStateFlag (GaussPoint *gp) const |
| Assign state flag. | |
| FloatArrayF< 6 > | computeDRhoDStress (const FloatArrayF< 6 > &stress) const |
| Computes the derivative of rho with respect to the stress. | |
| double | computeDRDCosTheta (const double theta, const double ecc) const |
| Computes the derivative of function r with respect to cos theta. | |
| double | computeDDRDDCosTheta (const double theta, const double ecc) const |
| Computes the second derivative of function r with respect to cos theta. | |
| FloatArrayF< 6 > | computeDSigDStress () const |
| Computes the derivative of sig with respect to the stress. | |
| FloatMatrixF< 6, 6 > | computeDDRhoDDStress (const FloatArrayF< 6 > &stress) const |
| Computes the second derivative of rho with the respect to the stress. | |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 6, 6 > | compute3dSecantStiffness (GaussPoint *gp, TimeStep *tStep) const |
| Compute the 3d secant stiffness matrix. | |
| FloatMatrixF< 6, 6 > | compute3dTangentStiffness (GaussPoint *gp, TimeStep *tStep) const |
| Compute the 3d tangent stiffness matrix. | |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| Public Member Functions inherited from oofem::StructuralMaterial | |
| StructuralMaterial (int n, Domain *d) | |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| const char * | giveClassName () const override |
| void | initializeFrom (InputRecord &ir) override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| void | giveCharacteristicMatrix (FloatMatrix &answer, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override |
| Returns characteristic matrix of the receiver. | |
| void | giveCharacteristicVector (FloatArray &answer, FloatArray &flux, MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const override |
| Returns characteristic vector of the receiver. | |
| virtual void | giveStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) const |
| virtual FloatArrayF< 4 > | giveRealStressVector_PlaneStrain (const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_3d. | |
| virtual FloatArray | giveRealStressVector_StressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| virtual FloatArray | giveRealStressVector_ShellStressControl (const FloatArray &reducedE, const IntArray &strainControl, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 3 > | giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 1 > | giveRealStressVector_1d (const FloatArrayF< 1 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 2 > | giveRealStressVector_Warping (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 2 > | giveRealStressVector_2dBeamLayer (const FloatArrayF< 2 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 5 > | giveRealStressVector_PlateLayer (const FloatArrayF< 5 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 3 > | giveRealStressVector_Fiber (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector_StressControl. | |
| virtual FloatArrayF< 3 > | giveRealStressVector_2dPlateSubSoil (const FloatArrayF< 3 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation is not provided. | |
| virtual FloatArrayF< 6 > | giveRealStressVector_3dBeamSubSoil (const FloatArrayF< 6 > &reducedE, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 9 > | giveFirstPKStressVector_3d (const FloatArrayF< 9 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| virtual FloatArrayF< 5 > | giveFirstPKStressVector_PlaneStrain (const FloatArrayF< 5 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_3d. | |
| virtual FloatArray | giveFirstPKStressVector_StressControl (const FloatArray &reducedvF, const IntArray &FControl, GaussPoint *gp, TimeStep *tStep) const |
| Iteratively calls giveRealStressVector_3d to find the stress controlled equal to zero· | |
| virtual FloatArrayF< 4 > | giveFirstPKStressVector_PlaneStress (const FloatArrayF< 4 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_StressControl. | |
| virtual FloatArrayF< 1 > | giveFirstPKStressVector_1d (const FloatArrayF< 1 > &vF, GaussPoint *gp, TimeStep *tStep) const |
| Default implementation relies on giveFirstPKStressVector_StressControl. | |
| 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 | giveEshelbyStressVector_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedF, TimeStep *tStep) |
| virtual FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const |
| double | giveReferenceTemperature () |
| virtual FloatArray | computeStressIndependentStrainVector (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const |
| FloatArrayF< 6 > | computeStressIndependentStrainVector_3d (GaussPoint *gp, TimeStep *tStep, ValueModeType mode) const |
| virtual void | giveStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 9, 9 > | give3dMaterialStiffnessMatrix_dPdF (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| 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) const |
| int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| virtual FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 4, 4 > | givePlaneStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | givePlaneStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 4, 4 > | givePlaneStrainStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | givePlaneStrainStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | givePlaneStrainStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 1, 1 > | give1dStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| virtual FloatMatrixF< 2, 2 > | give2dBeamLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | givePlateLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 3, 3 > | giveFiberStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 3, 3 > | give2dPlateSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 6, 6 > | give3dBeamSubSoilStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| Public Member Functions inherited from oofem::Material | |
| Material (int n, Domain *d) | |
| virtual | ~Material ()=default |
| Destructor. | |
| virtual double | giveCharacteristicValue (MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const |
| Returns characteristic value of the receiver. | |
| virtual double | give (int aProperty, GaussPoint *gp) const |
| virtual bool | hasProperty (int aProperty, GaussPoint *gp) const |
| virtual void | modifyProperty (int aProperty, double value, GaussPoint *gp) |
| double | giveCastingTime () const |
| virtual bool | isActivated (TimeStep *tStep) const |
| virtual bool | hasCastingTimeSupport () const |
| void | printYourself () override |
| Prints receiver state on stdout. Useful for debugging. | |
| virtual void | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| virtual void | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| int | checkConsistency () override |
| virtual void | restoreConsistency (GaussPoint *gp) |
| virtual int | initMaterial (Element *element) |
| virtual MaterialStatus * | giveStatus (GaussPoint *gp) const |
| virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
| virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip) |
| virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip) |
| virtual double | predictRelativeComputationalCost (GaussPoint *gp) |
| virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
| virtual void | initTempStatus (GaussPoint *gp) const |
| void | saveContext (DataStream &stream, ContextMode mode) override |
| void | restoreContext (DataStream &stream, ContextMode mode) override |
| Public Member Functions inherited from oofem::FEMComponent | |
| FEMComponent (int n, Domain *d) | |
| virtual | ~FEMComponent ()=default |
| Virtual destructor. | |
| Domain * | giveDomain () const |
| virtual void | setDomain (Domain *d) |
| int | giveNumber () const |
| void | setNumber (int num) |
| virtual void | updateLocalNumbering (EntityRenumberingFunctor &f) |
| virtual void | initializeFrom (InputRecord &ir, int priority) |
| virtual void | initializeFinish () |
| virtual void | postInitialize () |
| Performs post initialization steps. Called after all components are created and initialized. | |
| virtual void | printOutputAt (FILE *file, TimeStep *tStep) |
| virtual Interface * | giveInterface (InterfaceType t) |
| std::string | errorInfo (const char *func) const |
| Returns string for prepending output (used by error reporting macros). | |
Protected Member Functions | |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
Protected Attributes | |
| double | fc = 0. |
| 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. | |
| double | ft = 0. |
| double | ecc = 0. |
| int | damageFlag = 0 |
| double | e0 = 0. |
| double | AHard = 0. |
| Parameter of the ductilityMeasure of the plasticity model. | |
| double | BHard = 0. |
| Parameter of the ductilityMeasure of the plasticity model. | |
| double | CHard = 0. |
| Parameter of the ductilityMeasure of the plasticity model. | |
| double | DHard = 0. |
| Parameter of the ductilityMeasure of the plasticity model. | |
| double | hardeningModulus = 0. |
| Hardening modulus. | |
| double | ASoft = 0. |
| Parameter of the ductilityMeasure of the damage model. | |
| double | yieldHardPrimePeak = 0. |
| Parameter of the hardening law of the plasticity model. | |
| double | yieldHardInitial = 0. |
| Parameter of the hardening law of the plasticity model. | |
| double | dilationConst = 0. |
| Control parameter for te volumetric plastic flow of the plastic potential. | |
| double | m = 0. |
| Friction parameter of the yield surface. | |
| double | mQ = 0. |
| Dilation parameter of the plastic potential. | |
| double | helem = 0. |
| Element size (to be used in fracture energy approach (crack band). | |
| IsotropicLinearElasticMaterial | linearElasticMaterial |
| Pointer for linear elastic material. | |
| double | eM = 0. |
| Elastic Young's modulus. | |
| double | gM = 0. |
| Elastic shear modulus. | |
| double | kM = 0. |
| Elastic bulk modulus. | |
| double | nu = 0. |
| Elastic poisson's ration. | |
| double | efCompression = 0. |
| Control parameter for the exponential softening law. | |
| double | wf = 0. |
| Control parameter for the linear/bilinear softening law in tension. | |
| double | wfOne = 0. |
| Control parameter for the bilinear softening law in tension. | |
| double | ftOne = 0. |
| Control parameter for the bilinear softening law. | |
| double | yieldTol = 0. |
| yield tolerance for the plasticity model. | |
| double | yieldTolDamage = 0. |
| yield tolerance for the damage model. | |
| int | newtonIter = 0 |
| Maximum number of iterations for stress return. | |
| int | softeningType = 0 |
| Type of softening function used. | |
| double | deltaTime = 0. |
| Input parameter which simulates a loading rate. Only for debugging purposes. | |
| int | strengthRateType = 0 |
| int | energyRateType = 0 |
| Protected Attributes inherited from oofem::StructuralMaterial | |
| double | referenceTemperature = 0. |
| Reference temperature (temperature, when material has been built into structure). | |
| MatResponseMode | SCStiffMode = TangentStiffness |
| stifness mode used in stress control | |
| double | SCRelTol = 1.e-3 |
| relative tolerance for stress control | |
| double | SCAbsTol = 1.e-12 |
| absolute stress tolerance for stress control | |
| int | SCMaxiter = 100000 |
| maximum iterations for stress-control | |
| Protected Attributes inherited from oofem::Material | |
| Dictionary | propertyDictionary |
| double | castingTime |
| int | preCastingTimeMat |
| Material existing before casting time - optional parameter, zero by default. | |
| Protected Attributes inherited from oofem::FEMComponent | |
| int | number |
| Component number. | |
| Domain * | domain |
| Link to domain object, useful for communicating with other FEM components. | |
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 FloatMatrixF< 9, 9 > | convert_dSdE_2_dPdF_3D (const FloatMatrixF< 6, 6 > &dSdE, const FloatArrayF< 6 > &S, const FloatArrayF< 9 > &F) |
| static FloatMatrixF< 5, 5 > | convert_dSdE_2_dPdF_PlaneStrain (const FloatMatrixF< 4, 4 > &dSdE, const FloatArrayF< 4 > &S, const FloatArrayF< 5 > &F) |
| static FloatMatrixF< 4, 4 > | convert_dSdE_2_dPdF_PlaneStress (const FloatMatrixF< 3, 3 > &dSdE, const FloatArrayF< 3 > &S, const FloatArrayF< 4 > &F) |
| static FloatMatrixF< 1, 1 > | convert_dSdE_2_dPdF_1D (const FloatMatrixF< 1, 1 > &dSdE, const FloatArrayF< 1 > &S, const FloatArrayF< 1 > &F) |
| static void | computePrincipalValues (FloatArray &answer, const FloatArray &s, stressStrainPrincMode mode) |
| Common functions for convenience. | |
| static FloatArrayF< 3 > | computePrincipalValues (const FloatMatrixF< 3, 3 > &s) |
| static FloatArrayF< 3 > | computePrincipalValues (double I1, double I2, double I3) |
| static void | computePrincipalValDir (FloatArray &answer, FloatMatrix &dir, const FloatArray &s, stressStrainPrincMode mode) |
| static std::pair< FloatArrayF< 3 >, FloatMatrixF< 3, 3 > > | computePrincipalValDir (const FloatMatrixF< 3, 3 > &s) |
| static FloatArrayF< 6 > | computeDeviator (const FloatArrayF< 6 > &s) |
| static std::pair< FloatArrayF< 6 >, double > | computeDeviatoricVolumetricSplit (const FloatArrayF< 6 > &s) |
| static FloatArrayF< 6 > | computeDeviatoricVolumetricSum (const FloatArrayF< 6 > &dev, double mean) |
| static FloatArrayF< 6 > | applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyDeviatoricElasticCompliance (const FloatArrayF< 6 > &stress, double GModulus) |
| static FloatArrayF< 6 > | applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyDeviatoricElasticStiffness (const FloatArrayF< 6 > &strain, double GModulus) |
| static FloatArrayF< 6 > | applyElasticStiffness (const FloatArrayF< 6 > &strain, double EModulus, double nu) |
| static FloatArrayF< 6 > | applyElasticCompliance (const FloatArrayF< 6 > &stress, double EModulus, double nu) |
| static double | computeStressNorm (const FloatArrayF< 6 > &stress) |
| static double | computeFirstInvariant (const FloatArrayF< 6 > &s) |
| static double | computeSecondStressInvariant (const FloatArrayF< 6 > &s) |
| static double | computeThirdStressInvariant (const FloatArrayF< 6 > &s) |
| static double | computeFirstCoordinate (const FloatArrayF< 6 > &s) |
| static double | computeSecondCoordinate (const FloatArrayF< 6 > &s) |
| static double | computeThirdCoordinate (const FloatArrayF< 6 > &s) |
| static int | giveVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| static int | giveVoigtSymVectorMask (IntArray &answer, MaterialMode mmode) |
| static void | giveInvertedVoigtVectorMask (IntArray &answer, MaterialMode mmode) |
| static int | giveSizeOfVoigtVector (MaterialMode mmode) |
| static int | giveSizeOfVoigtSymVector (MaterialMode mmode) |
| static void | giveFullVectorForm (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced symmetric Voigt vector (2nd order tensor) to full form. | |
| static void | giveFullVectorFormF (FloatArray &answer, const FloatArray &strainVector, MaterialMode matMode) |
| Converts the reduced deformation gradient Voigt vector (2nd order tensor). | |
| static void | giveFullSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the reduced unsymmetric Voigt vector (2nd order tensor) to full form. | |
| static void | giveReducedVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full symmetric Voigt vector (2nd order tensor) to reduced form. | |
| static void | giveReducedSymVectorForm (FloatArray &answer, const FloatArray &vec, MaterialMode matMode) |
| Converts the full unsymmetric Voigt vector (2nd order tensor) to reduced form. | |
| static void | giveFullSymMatrixForm (FloatMatrix &answer, const FloatMatrix &red, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. | |
| static void | giveReducedMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full symmetric Voigt matrix (4th order tensor) to reduced form. | |
| static void | giveReducedSymMatrixForm (FloatMatrix &answer, const FloatMatrix &full, MaterialMode matMode) |
| Converts the full unsymmetric Voigt matrix (4th order tensor) to reduced form. | |
| static FloatArrayF< 6 > | transformStrainVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &strain, bool transpose=false) |
| static FloatArrayF< 6 > | transformStressVectorTo (const FloatMatrixF< 3, 3 > &base, const FloatArrayF< 6 > &stress, bool transpose=false) |
| static double | computeVonMisesStress (const FloatArray ¤tStress) |
| static double | computeVonMisesStress_3D (const FloatArrayF< 6 > &stress) |
| static double | computeVonMisesStress_PlaneStress (const FloatArrayF< 3 > &stress) |
| static FloatMatrixF< 6, 6 > | giveStrainVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false) |
| static FloatMatrixF< 3, 3 > | give2DStrainVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false) |
| static FloatMatrixF< 6, 6 > | giveStressVectorTranformationMtrx (const FloatMatrixF< 3, 3 > &base, bool transpose=false) |
| static FloatMatrixF< 3, 3 > | givePlaneStressVectorTranformationMtrx (const FloatMatrixF< 2, 2 > &base, bool transpose=false) |
| static void | sortPrincDirAndValCloseTo (FloatArray &pVal, FloatMatrix &pDir, const FloatMatrix &toPDir) |
| Static Public Attributes inherited from oofem::StructuralMaterial | |
| static std::array< std::array< int, 3 >, 3 > | vIindex |
| Voigt index map. | |
| static std::array< std::array< int, 3 >, 3 > | svIndex |
| Symmetric Voigt index map. | |
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 627 of file concretedpm2.h.
| Enumerator | |
|---|---|
| RR_Unknown | |
| RR_NotConverged | |
| RR_Converged | |
Definition at line 631 of file concretedpm2.h.
| Enumerator | |
|---|---|
| RT_Unknown | |
| RT_Regular | |
| RT_Tension | |
| RT_Compression | |
| RT_Auxiliary | |
Definition at line 637 of file concretedpm2.h.
| oofem::ConcreteDPM2::ConcreteDPM2 | ( | int | n, |
| Domain * | d ) |
Constructor.
Definition at line 472 of file concretedpm2.C.
References linearElasticMaterial, and oofem::StructuralMaterial::StructuralMaterial().
| void oofem::ConcreteDPM2::assignStateFlag | ( | GaussPoint * | gp | ) | const |
Assign state flag.
Definition at line 2767 of file concretedpm2.C.
References oofem::ConcreteDPM2Status::ConcreteDPM2_Damage, oofem::ConcreteDPM2Status::ConcreteDPM2_Elastic, oofem::ConcreteDPM2Status::ConcreteDPM2_Plastic, oofem::ConcreteDPM2Status::ConcreteDPM2_PlasticDamage, oofem::ConcreteDPM2Status::ConcreteDPM2_Unloading, oofem::ConcreteDPM2Status::ConcreteDPM2_VertexCompressionDamage, oofem::ConcreteDPM2Status::ConcreteDPM2_VertexTension, oofem::ConcreteDPM2Status::ConcreteDPM2_VertexTensionDamage, oofem::ConcreteDPM2Status::giveDamageTension(), and oofem::Material::giveStatus().
Referenced by giveRealStressVector_3d().
| int oofem::ConcreteDPM2::checkForUnAndReloading | ( | double & | tempEquivStrain, |
| double & | minEquivStrain, | ||
| const FloatMatrixF< 6, 6 > & | D, | ||
| GaussPoint * | gp ) const |
Check for un- and reloading in the damage part.
Definition at line 1046 of file concretedpm2.C.
References computeCoordinates(), computeEquivalentStrain(), oofem::dot(), oofem::ConcreteDPM2Status::giveReducedStrain(), oofem::Material::giveStatus(), and yieldTolDamage.
Referenced by computeDamage().
| void oofem::ConcreteDPM2::checkForVertexCase | ( | double & | answer, |
| ConcreteDPM2_ReturnType & | returnType, | ||
| double | sig, | ||
| double | tempKappa, | ||
| GaussPoint * | gp ) const |
Check if the trial stress state falls within the vertex region of the plasticity model at the apex of triaxial extension or triaxial compression.
| answer | Volumetric apex stress. |
| sig | Volumetric stress. |
| tempKappa | Hardening variable. |
| gp | Gauss point |
Definition at line 1604 of file concretedpm2.C.
References RT_Compression, RT_Regular, and RT_Tension.
Referenced by performPlasticityReturn().
| FloatMatrixF< 6, 6 > oofem::ConcreteDPM2::compute3dSecantStiffness | ( | GaussPoint * | gp, |
| TimeStep * | tStep ) const |
Compute the 3d secant stiffness matrix.
Definition at line 2692 of file concretedpm2.C.
References damageFlag, oofem::Material::giveStatus(), linearElasticMaterial, and oofem::min().
Referenced by give3dMaterialStiffnessMatrix().
| FloatMatrixF< 6, 6 > oofem::ConcreteDPM2::compute3dTangentStiffness | ( | GaussPoint * | gp, |
| TimeStep * | tStep ) const |
Compute the 3d tangent stiffness matrix.
Definition at line 2714 of file concretedpm2.C.
References oofem::FloatMatrixF< N, M >::at(), computeFullJacobian(), damageFlag, oofem::dot(), oofem::Material::giveStatus(), oofem::inv(), linearElasticMaterial, and oofem::min().
Referenced by give3dMaterialStiffnessMatrix().
| double oofem::ConcreteDPM2::computeAlpha | ( | FloatArrayF< 6 > & | effectiveStressTension, |
| FloatArrayF< 6 > & | effectiveStressCompression, | ||
| const FloatArrayF< 6 > & | effectiveStress ) const |
Compute alpha for rate effect.
Definition at line 2565 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::from_voigt_stress_6(), oofem::norm_squared(), and oofem::StructuralMaterial::transformStressVectorTo().
Referenced by giveRealStressVector_3d().
| void oofem::ConcreteDPM2::computeCoordinates | ( | const FloatArrayF< 6 > & | stress, |
| double & | sig, | ||
| double & | rho, | ||
| double & | theta ) const |
Compute the Haigh-Westergaard coordinates.
Definition at line 2757 of file concretedpm2.C.
References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::StructuralMaterial::computeThirdCoordinate().
Referenced by checkForUnAndReloading(), computeDamage(), computeDFDStress(), computeFullJacobian(), and performPlasticityReturn().
| FloatArrayF< 2 > oofem::ConcreteDPM2::computeDamage | ( | const FloatArrayF< 6 > & | strain, |
| const FloatMatrixF< 6, 6 > & | D, | ||
| double | timeFactor, | ||
| GaussPoint * | gp, | ||
| TimeStep * | tStep, | ||
| double | alpha, | ||
| const FloatArrayF< 6 > & | effectiveStress ) const |
Compute damage parameters.
Definition at line 876 of file concretedpm2.C.
References checkForUnAndReloading(), computeCoordinates(), computeDamageParamCompression(), computeDamageParamTension(), computeDeltaPlasticStrainNormCompression(), computeDeltaPlasticStrainNormTension(), computeDuctilityMeasureDamage(), computeRateFactor(), deltaTime, e0, oofem::Material::giveStatus(), initDamaged(), oofem::TimeStep::isTheFirstStep(), and yieldTolDamage.
Referenced by giveRealStressVector_3d().
|
virtual |
Compute damage parameter in compression.
Definition at line 1361 of file concretedpm2.C.
References damageFlag, e0, efCompression, eM, energyRateType, ft, newtonIter, OOFEM_ERROR, strengthRateType, and yieldTolDamage.
Referenced by computeDamage().
|
virtual |
Compute damage parameter in tension.
Definition at line 1283 of file concretedpm2.C.
References e0, eM, energyRateType, ft, ftOne, newtonIter, OOFEM_ERROR, softeningType, strengthRateType, wf, wfOne, and yieldTolDamage.
Referenced by computeDamage().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDCosThetaDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the derivative of cos theta with respect to the stress.
Definition at line 2863 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), oofem::StructuralMaterial::computeDeviator(), computeDRhoDStress(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::from_voigt_stress_6().
Referenced by computeDDKappaDDeltaLambdaDStress(), and computeDFDStress().
| FloatMatrixF< 6, 6 > oofem::ConcreteDPM2::computeDDCosThetaDDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the second derivative of cos theta with respect to the stress.
Definition at line 2894 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), oofem::FloatMatrixF< N, M >::at(), computeDDRhoDDStress(), oofem::StructuralMaterial::computeDeviator(), computeDRhoDStress(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::from_voigt_stress_6().
| FloatMatrixF< 2, 2 > oofem::ConcreteDPM2::computeDDGDDInv | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
Computes the second derivative of the plastic potential with respect to the invariants sig and rho are computed.
Definition at line 2509 of file concretedpm2.C.
References oofem::FloatMatrixF< N, M >::at(), computeHardeningOne(), computeHardeningTwo(), dilationConst, fc, ft, and m.
Referenced by computeDDGDDStress(), computeDDKappaDDeltaLambdaDInv(), and computeJacobian().
| FloatMatrixF< 6, 6 > oofem::ConcreteDPM2::computeDDGDDStress | ( | const FloatArrayF< 6 > & | stress, |
| const double | tempKappa ) const |
Computes the second derivative of the plastic potential with respect to the stress.
Definition at line 2379 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), oofem::FloatMatrixF< N, M >::at(), computeDDGDDInv(), computeDDRhoDDStress(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeDGDInv(), computeDRhoDStress(), computeDSigDStress(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::dyad().
Referenced by computeFullJacobian().
| FloatArrayF< 2 > oofem::ConcreteDPM2::computeDDGDInvDKappa | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
The mixed derivative of the plastic potential with respect to the invariants and the hardening parameter are determined.
Definition at line 2443 of file concretedpm2.C.
References computeHardeningOne(), computeHardeningOnePrime(), computeHardeningTwo(), computeHardeningTwoPrime(), dilationConst, fc, ft, m, and mQ.
Referenced by computeDDGDStressDKappa(), computeDDKappaDDeltaLambdaDKappa(), and computeJacobian().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDDGDStressDKappa | ( | const FloatArrayF< 6 > & | stress, |
| double | tempKappa ) const |
Computes the second mixed derivative of plastic potential with respect to the stress and kappa.
Definition at line 2150 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), computeDDGDInvDKappa(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeDRhoDStress(), computeDSigDStress(), and oofem::StructuralMaterial::computeSecondCoordinate().
Referenced by computeFullJacobian().
| FloatArrayF< 2 > oofem::ConcreteDPM2::computeDDKappaDDeltaLambdaDInv | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
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 2078 of file concretedpm2.C.
References computeDDGDDInv(), computeDDuctilityMeasureDInv(), computeDGDInv(), and computeDuctilityMeasure().
Referenced by computeDDKappaDDeltaLambdaDStress(), and computeJacobian().
| double oofem::ConcreteDPM2::computeDDKappaDDeltaLambdaDKappa | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Computes the derivative of the evolution law of the hardening parameter kappa with respect to the hardening variable kappa.
Definition at line 2180 of file concretedpm2.C.
References computeDDGDInvDKappa(), computeDGDInv(), and computeDuctilityMeasure().
Referenced by computeFullJacobian(), and computeJacobian().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDDKappaDDeltaLambdaDStress | ( | const FloatArrayF< 6 > & | stress, |
| double | tempKappa ) const |
Computes the derivative of evolution law of hardening variable with respect to the stress.
Definition at line 2111 of file concretedpm2.C.
References computeDCosThetaDStress(), computeDDKappaDDeltaLambdaDInv(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeDGDInv(), computeDRhoDStress(), computeDSigDStress(), computeDuctilityMeasure(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::StructuralMaterial::computeThirdCoordinate().
Referenced by computeFullJacobian().
| double oofem::ConcreteDPM2::computeDDRDDCosTheta | ( | const double | theta, |
| const double | ecc ) const |
Computes the second derivative of function r with respect to cos theta.
Definition at line 2827 of file concretedpm2.C.
References ecc.
| FloatMatrixF< 6, 6 > oofem::ConcreteDPM2::computeDDRhoDDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the second derivative of rho with the respect to the stress.
Definition at line 2992 of file concretedpm2.C.
References oofem::FloatMatrixF< N, M >::at(), oofem::StructuralMaterial::computeDeviator(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::dyad().
Referenced by computeDDCosThetaDDStress(), and computeDDGDDStress().
| FloatArrayF< 2 > oofem::ConcreteDPM2::computeDDuctilityMeasureDInv | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Compute derivative the ductility measure with respect to the stress state.
| answer | array of the derivative of the ductility measure with respect to volumetric and deviatoric stress |
| sig | Volumetric stress. |
| rho | Length of the deviatoric strength. |
| theta | Lode angle. |
| tempKappa | plastic strain |
Definition at line 2197 of file concretedpm2.C.
References AHard, BHard, CHard, DHard, and fc.
Referenced by computeDDKappaDDeltaLambdaDInv().
| double oofem::ConcreteDPM2::computeDeltaPlasticStrainNormCompression | ( | double | tempAlpha, |
| double | tempKappaD, | ||
| double | kappaD, | ||
| GaussPoint * | gp, | ||
| const double | rho ) const |
Compute equivalent strain value for compression.
Definition at line 1226 of file concretedpm2.C.
References computeHardeningTwo(), dilationConst, e0, ft, oofem::Material::giveStatus(), oofem::ConcreteDPM2Status::giveTempPlasticStrain(), oofem::norm(), and yieldTolDamage.
Referenced by computeDamage().
| double oofem::ConcreteDPM2::computeDeltaPlasticStrainNormTension | ( | double | tempKappaD, |
| double | kappaD, | ||
| GaussPoint * | gp ) const |
Compute equivalent strain value for tension.
Definition at line 1197 of file concretedpm2.C.
References e0, oofem::Material::giveStatus(), oofem::ConcreteDPM2Status::giveTempPlasticStrain(), oofem::norm(), and yieldTolDamage.
Referenced by computeDamage().
| FloatArrayF< 2 > oofem::ConcreteDPM2::computeDFDInv | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Computes the derivative of the yield surface with respect to the invariants sig and rho.
Definition at line 2035 of file concretedpm2.C.
References oofem::AL, computeHardeningOne(), computeHardeningTwo(), ecc, fc, and m.
Referenced by computeDFDStress(), and computeJacobian().
| double oofem::ConcreteDPM2::computeDFDKappa | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Compute the derivative of the yield surface with respect to the hardening variable based on the stress state and the hardening variable
| sig | Volumetric stress. |
| rho | Deviatoric length. |
| theta | Lode angle. |
| tempKappa | Hardening variable. |
Definition at line 1990 of file concretedpm2.C.
References computeHardeningOne(), computeHardeningOnePrime(), computeHardeningTwo(), computeHardeningTwoPrime(), ecc, fc, and m.
Referenced by computeFullJacobian(), and computeJacobian().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDFDStress | ( | const FloatArrayF< 6 > & | stress, |
| const double | tempKappa ) const |
Computes the derivative of the yield surface with respect to the stress.
Definition at line 2258 of file concretedpm2.C.
References computeCoordinates(), computeDCosThetaDStress(), computeDFDInv(), computeDRDCosTheta(), computeDRhoDStress(), computeDSigDStress(), computeHardeningOne(), computeHardeningTwo(), ecc, fc, and m.
Referenced by computeFullJacobian().
| FloatArrayF< 2 > oofem::ConcreteDPM2::computeDGDInv | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
Compute derivative the palstic potential function with respect to the stress state.
| answer | array of the derivative of the plastic potential with respect to volumetric and deviatoric stress |
| sig | volumetric stress. |
| rho | deviatoric stress. |
| tempKappa | plastic strain |
Definition at line 2226 of file concretedpm2.C.
References computeHardeningOne(), computeHardeningTwo(), dilationConst, fc, ft, m, and mQ.
Referenced by computeDDGDDStress(), computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDDKappaDDeltaLambdaDStress(), computeDGDStress(), computeDKappaDDeltaLambda(), computeJacobian(), and performRegularReturn().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDGDStress | ( | const FloatArrayF< 6 > & | stress, |
| const double | tempKappa ) const |
Computes the derivative of the plastic potential with respect to the stress.
Definition at line 2287 of file concretedpm2.C.
References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeDGDInv(), computeDRhoDStress(), computeDSigDStress(), and oofem::StructuralMaterial::computeSecondCoordinate().
Referenced by computeFullJacobian().
| double oofem::ConcreteDPM2::computeDKappaDDeltaLambda | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Compute the derivative of kappa with respect of delta lambda based on the stress state and the hardening variable.
| sig | Volumetric stress. |
| rho | Length of the deviatoric stress. |
| theta | Lode angle |
| tempKappa | Hardening variable. |
Definition at line 2066 of file concretedpm2.C.
References computeDGDInv(), and computeDuctilityMeasure().
Referenced by computeFullJacobian(), computeJacobian(), and performRegularReturn().
| double oofem::ConcreteDPM2::computeDRDCosTheta | ( | const double | theta, |
| const double | ecc ) const |
Computes the derivative of function r with respect to cos theta.
Definition at line 2809 of file concretedpm2.C.
References ecc.
Referenced by computeDFDStress().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDRhoDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the derivative of rho with respect to the stress.
Definition at line 2966 of file concretedpm2.C.
References oofem::StructuralMaterial::computeDeviator(), and oofem::StructuralMaterial::computeSecondCoordinate().
Referenced by computeDCosThetaDStress(), computeDDCosThetaDDStress(), computeDDGDDStress(), computeDDGDStressDKappa(), computeDDKappaDDeltaLambdaDStress(), computeDFDStress(), and computeDGDStress().
| FloatArrayF< 6 > oofem::ConcreteDPM2::computeDSigDStress | ( | ) | const |
Computes the derivative of sig with respect to the stress.
Definition at line 2983 of file concretedpm2.C.
Referenced by computeDDGDDStress(), computeDDGDStressDKappa(), computeDDKappaDDeltaLambdaDStress(), computeDFDStress(), and computeDGDStress().
|
virtual |
Compute the ductility measure based on the stress state.
| sig | Volumetric stress. |
| rho | Length of the deviatoric strength. |
| theta | Lode angle of stress state. |
Definition at line 1739 of file concretedpm2.C.
References AHard, BHard, CHard, DHard, and fc.
Referenced by computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDDKappaDDeltaLambdaDStress(), computeDKappaDDeltaLambda(), and computeTempKappa().
| double oofem::ConcreteDPM2::computeDuctilityMeasureDamage | ( | GaussPoint * | gp, |
| const double | sig, | ||
| const double | rho ) const |
Compute the ductility measure for the damage model.
Definition at line 1462 of file concretedpm2.C.
References ASoft.
Referenced by computeDamage().
|
virtual |
Compute the base equivalent strain value.
Definition at line 1262 of file concretedpm2.C.
References e0, ecc, fc, and m.
Referenced by checkForUnAndReloading().
| FloatMatrixF< 8, 8 > oofem::ConcreteDPM2::computeFullJacobian | ( | const FloatArrayF< 6 > & | stress, |
| const double | deltaLambda, | ||
| GaussPoint * | gp, | ||
| TimeStep * | atTime, | ||
| const double | tempKappa ) const |
Computes full Jacobian used for the algorithmic tangent stiffness
Definition at line 2310 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), oofem::FloatMatrixF< N, M >::at(), computeCoordinates(), computeDDGDDStress(), computeDDGDStressDKappa(), computeDDKappaDDeltaLambdaDKappa(), computeDDKappaDDeltaLambdaDStress(), computeDFDKappa(), computeDFDStress(), computeDGDStress(), computeDKappaDDeltaLambda(), oofem::dot(), and linearElasticMaterial.
Referenced by compute3dTangentStiffness().
| double oofem::ConcreteDPM2::computeHardeningOne | ( | double | tempKappa | ) | const |
Compute the value of the hardening function based on the hardening variable.
| tempKappa | Hardening variable. |
Definition at line 2613 of file concretedpm2.C.
References yieldHardInitial, and yieldHardPrimePeak.
Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDFDInv(), computeDFDKappa(), computeDFDStress(), computeDGDInv(), computeRatioPotential(), and computeYieldValue().
| double oofem::ConcreteDPM2::computeHardeningOnePrime | ( | double | tempKappa | ) | const |
Compute the derivative of the hardening function based on the hardening parameter.
| tempKappa | Hardening variable. |
Definition at line 2630 of file concretedpm2.C.
References yieldHardInitial, and yieldHardPrimePeak.
Referenced by computeDDGDInvDKappa(), and computeDFDKappa().
| double oofem::ConcreteDPM2::computeHardeningTwo | ( | double | tempKappa | ) | const |
Compute the value of the hardening function based on the hardening variable.
| tempKappa | Hardening variable. |
Definition at line 2646 of file concretedpm2.C.
References yieldHardPrimePeak.
Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDeltaPlasticStrainNormCompression(), computeDFDInv(), computeDFDKappa(), computeDFDStress(), computeDGDInv(), computeRatioPotential(), and computeYieldValue().
| double oofem::ConcreteDPM2::computeHardeningTwoPrime | ( | double | tempKappa | ) | const |
Compute the derivative of the hardening function based on the hardening parameter.
| tempKappa | Hardening variable. |
Definition at line 2659 of file concretedpm2.C.
References yieldHardPrimePeak.
Referenced by computeDDGDInvDKappa(), and computeDFDKappa().
| FloatMatrixF< 4, 4 > oofem::ConcreteDPM2::computeJacobian | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa, | ||
| double | deltaLambda, | ||
| GaussPoint * | gp ) const |
Compute jacobian for 2D(plane strain) and 3d cases
| sig | volumetric strain |
| rho | deviatoric |
| tempKappa | plastic strain |
| deltaLambda | plastic multiplier |
| gp | Gauss point |
Definition at line 1914 of file concretedpm2.C.
References oofem::FloatMatrixF< N, M >::at(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeDKappaDDeltaLambda(), gM, and kM.
Referenced by performRegularReturn().
| double oofem::ConcreteDPM2::computeRateFactor | ( | double | alpha, |
| double | timeFactor, | ||
| GaussPoint * | gp, | ||
| TimeStep * | deltaTime ) const |
This function computes the rate factor which is used to take into account the strain rate dependence of the material.
Definition at line 1109 of file concretedpm2.C.
References CDPM2_TOL, oofem::StructuralMaterial::computePrincipalValues(), deltaTime, oofem::from_voigt_strain_6(), oofem::Material::giveStatus(), oofem::ConcreteDPM2Status::giveTempReducedStrain(), and strengthRateType.
Referenced by computeDamage().
| double oofem::ConcreteDPM2::computeRatioPotential | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
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 1761 of file concretedpm2.C.
References computeHardeningOne(), computeHardeningTwo(), dilationConst, fc, ft, m, mQ, and nu.
Referenced by performVertexReturn().
| double oofem::ConcreteDPM2::computeTempKappa | ( | double | kappaInitial, |
| double | sigTrial, | ||
| double | rhoTrial, | ||
| double | sig ) const |
Compute tempKappa.
Definition at line 1722 of file concretedpm2.C.
References computeDuctilityMeasure(), gM, kM, and M_PI.
Referenced by performVertexReturn().
| double oofem::ConcreteDPM2::computeYieldValue | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Compute the yield value based on stress and hardening variable.
| sig | Volumetric stress. |
| rho | Length of the deviatoric stress. |
| theta | Lode angle of the stress state. |
| tempKappa | Hardening variable. |
Definition at line 1961 of file concretedpm2.C.
References computeHardeningOne(), computeHardeningTwo(), ecc, fc, and m.
Referenced by performPlasticityReturn(), performRegularReturn(), and performVertexReturn().
|
overrideprotectedvirtual |
Creates new copy of associated status and inserts it into given integration point.
| gp | Integration point where newly created status will be stored. |
Reimplemented from oofem::Material.
Definition at line 3070 of file concretedpm2.C.
|
overridevirtual |
Computes full 3d material stiffness matrix at given integration point, time, respecting load history in integration point.
| answer | Computed results. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 2671 of file concretedpm2.C.
References compute3dSecantStiffness(), compute3dTangentStiffness(), oofem::ConcreteDPM2Status::ConcreteDPM2_Plastic, oofem::ConcreteDPM2Status::ConcreteDPM2_PlasticDamage, oofem::Material::giveStatus(), and linearElasticMaterial.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 750 of file concretedpm2.h.
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 751 of file concretedpm2.h.
References _IFT_ConcreteDPM2_Name.
|
overridevirtual |
Returns the integration point corresponding value in Reduced form.
| answer | Contain corresponding ip value, zero sized if not available. |
| gp | Integration point to which the value refers. |
| type | Determines the type of internal variable. |
| tStep | Determines the time step. |
Reimplemented from oofem::Material.
Definition at line 3027 of file concretedpm2.C.
References oofem::FloatArray::at(), oofem::StructuralMaterial::giveIPValue(), oofem::Material::giveStatus(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 811 of file concretedpm2.C.
References assignStateFlag(), computeAlpha(), computeDamage(), oofem::StructuralMaterial::computeStressIndependentStrainVector_3d(), damageFlag, deltaTime, eM, ft, oofem::Material::giveStatus(), oofem::TimeStep::giveTimeIncrement(), oofem::ConcreteDPM2Status::initTempStatus(), linearElasticMaterial, OOFEM_ERROR, and performPlasticityReturn().
|
overridevirtual |
Tests if material supports material mode.
| mode | Required material mode. |
Reimplemented from oofem::Material.
Definition at line 479 of file concretedpm2.C.
| void oofem::ConcreteDPM2::initDamaged | ( | double | kappa, |
| const FloatArrayF< 6 > & | strain, | ||
| GaussPoint * | gp ) const |
Initialize the characteristic length, if damage is not yet activated Set the increase factor for the strain rate dependence
Definition at line 1412 of file concretedpm2.C.
References oofem::FloatArray::at(), oofem::Element::computeMeanSize(), oofem::StructuralMaterial::computePrincipalValDir(), e0, oofem::from_voigt_strain_6(), oofem::Element::giveCharacteristicLength(), oofem::GaussPoint::giveElement(), oofem::Material::giveStatus(), helem, and yieldTolDamage.
Referenced by computeDamage().
|
overridevirtual |
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. Note that initializeFrom may be called mutiple times.
| ir | Input record to initialize from. |
| priority | Priority of the input record. This is used to determine the order of initialization |
Reimplemented from oofem::FEMComponent.
Definition at line 489 of file concretedpm2.C.
References _IFT_ConcreteDPM2_ahard, _IFT_ConcreteDPM2_asoft, _IFT_ConcreteDPM2_bhard, _IFT_ConcreteDPM2_chard, _IFT_ConcreteDPM2_damflag, _IFT_ConcreteDPM2_deltatime, _IFT_ConcreteDPM2_dhard, _IFT_ConcreteDPM2_dilation, _IFT_ConcreteDPM2_ecc, _IFT_ConcreteDPM2_efc, _IFT_ConcreteDPM2_energyratetype, _IFT_ConcreteDPM2_fc, _IFT_ConcreteDPM2_ft, _IFT_ConcreteDPM2_ftOne, _IFT_ConcreteDPM2_helem, _IFT_ConcreteDPM2_hp, _IFT_ConcreteDPM2_kinit, _IFT_ConcreteDPM2_newtoniter, _IFT_ConcreteDPM2_softeningType, _IFT_ConcreteDPM2_strengthratetype, _IFT_ConcreteDPM2_wf, _IFT_ConcreteDPM2_wfOne, _IFT_ConcreteDPM2_yieldtol, _IFT_IsotropicLinearElasticMaterial_e, _IFT_IsotropicLinearElasticMaterial_n, AHard, ASoft, BHard, CHard, damageFlag, deltaTime, DHard, dilationConst, e0, ecc, efCompression, eM, energyRateType, fc, ft, ftOne, gM, helem, oofem::StructuralMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, kM, linearElasticMaterial, m, newtonIter, nu, OOFEM_ERROR, OOFEM_WARNING, oofem::Material::propertyDictionary, softeningType, strengthRateType, wf, wfOne, yieldHardInitial, yieldHardPrimePeak, yieldTol, and yieldTolDamage.
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 1119 of file concretedpm2.h.
| FloatArrayF< 6 > oofem::ConcreteDPM2::performPlasticityReturn | ( | GaussPoint * | gp, |
| const FloatMatrixF< 6, 6 > & | D, | ||
| const FloatArrayF< 6 > & | strain ) const |
Perform stress return of the plasticity model and compute history variables.
| gp | Gauss point. |
| D | stiffness matrix |
| strain | Strain vector of this Gauss point. |
Definition at line 1483 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), checkForVertexCase(), computeCoordinates(), computeYieldValue(), oofem::ConcreteDPM2Status::ConcreteDPM2_VertexCompression, oofem::ConcreteDPM2Status::ConcreteDPM2_VertexTension, oofem::dot(), oofem::GaussPoint::giveElement(), oofem::Element::giveGlobalNumber(), oofem::Material::giveStatus(), oofem::inv(), OOFEM_ERROR, OOFEM_LOG_INFO, performRegularReturn(), performVertexReturn(), RR_Converged, RR_NotConverged, RR_Unknown, RT_Compression, RT_Regular, RT_Tension, and RT_Unknown.
Referenced by giveRealStressVector_3d().
| double oofem::ConcreteDPM2::performRegularReturn | ( | FloatArrayF< 6 > & | stress, |
| ConcreteDPM2_ReturnResult & | returnResult, | ||
| ConcreteDPM2_ReturnType & | returnType, | ||
| double | kappaP, | ||
| GaussPoint * | gp, | ||
| double | theta ) const |
Perform regular stress return for the plasticity model, i.e. if the trial stress state does not lie in the vertex region.
| stress | Stress vector which is computed. |
| kappaP | Initial guess for kappa P (i.e. previous kappa) |
| gp | Gauss point. |
| theta | Load angle of trial stress (remains constant throughout return). |
Definition at line 1793 of file concretedpm2.C.
References oofem::FloatArray::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeDGDInv(), computeDKappaDDeltaLambda(), computeJacobian(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondCoordinate(), computeYieldValue(), oofem::from_voigt_stress_6(), oofem::Material::giveStatus(), gM, kM, M_PI, oofem::max(), newtonIter, oofem::norm(), oofem::FloatArray::resize(), RR_Converged, RR_NotConverged, oofem::solve(), oofem::StructuralMaterial::transformStressVectorTo(), and yieldTol.
Referenced by performPlasticityReturn().
| double oofem::ConcreteDPM2::performVertexReturn | ( | FloatArrayF< 6 > & | stress, |
| ConcreteDPM2_ReturnResult & | returnResult, | ||
| ConcreteDPM2_ReturnType & | returnType, | ||
| double | apexStress, | ||
| double | tempKappaP, | ||
| GaussPoint * | gp ) const |
Perform stress return for vertex case of the plasticity model, i.e. if the trial stress state lies within the vertex region.
| stress | Stress vector of this Gauss point. |
| apexStress | Volumetric stress at the apex of the yield surface. |
| theta | Lode angle. |
| tempKappaP | temporary cummulative plastic strain |
| gp | Gauss point. |
Definition at line 1624 of file concretedpm2.C.
References oofem::FloatArrayF< N >::at(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeRatioPotential(), oofem::StructuralMaterial::computeSecondCoordinate(), computeTempKappa(), computeYieldValue(), OOFEM_WARNING, RR_Converged, RR_NotConverged, RT_Compression, RT_Regular, RT_Tension, and yieldTol.
Referenced by performPlasticityReturn().
|
overridevirtual |
Restores the receiver state previously written in stream.
| stream | Input stream. |
| mode | Determines amount of info available in stream (state, definition, ...). |
| throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::FEMComponent.
Definition at line 703 of file concretedpm2.C.
References AHard, ASoft, BHard, CHard, oofem::CIO_IOERR, CM_Definition, damageFlag, deltaTime, DHard, dilationConst, e0, ecc, efCompression, eM, energyRateType, fc, ft, ftOne, gM, hardeningModulus, helem, kM, linearElasticMaterial, m, mQ, newtonIter, nu, oofem::DataStream::read(), oofem::FEMComponent::restoreContext(), softeningType, strengthRateType, THROW_CIOERR, wf, wfOne, yieldHardInitial, yieldHardPrimePeak, yieldTol, and yieldTolDamage.
|
overridevirtual |
Stores receiver state to output stream.
| stream | Output stream. |
| mode | Determines amount of info required in stream (state, definition, ...). |
| throws | an ContextIOERR exception if error encountered. |
Reimplemented from oofem::FEMComponent.
Definition at line 596 of file concretedpm2.C.
References AHard, ASoft, BHard, CHard, oofem::CIO_IOERR, CM_Definition, damageFlag, deltaTime, DHard, dilationConst, e0, ecc, efCompression, eM, energyRateType, fc, ft, ftOne, gM, hardeningModulus, helem, kM, linearElasticMaterial, m, mQ, newtonIter, nu, oofem::FEMComponent::saveContext(), softeningType, strengthRateType, THROW_CIOERR, wf, wfOne, oofem::DataStream::write(), yieldHardInitial, yieldHardPrimePeak, yieldTol, and yieldTolDamage.
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 658 of file concretedpm2.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the damage model.
Definition at line 670 of file concretedpm2.h.
Referenced by computeDuctilityMeasureDamage(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 660 of file concretedpm2.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 662 of file concretedpm2.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter which controls the type of damage that is used together with plasticity 0: tensile and compression damage based on split of stress state
Definition at line 653 of file concretedpm2.h.
Referenced by compute3dSecantStiffness(), compute3dTangentStiffness(), computeDamageParamCompression(), giveRealStressVector_3d(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Input parameter which simulates a loading rate. Only for debugging purposes.
Definition at line 727 of file concretedpm2.h.
Referenced by computeDamage(), computeRateFactor(), giveRealStressVector_3d(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 664 of file concretedpm2.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Control parameter for te volumetric plastic flow of the plastic potential.
Definition at line 679 of file concretedpm2.h.
Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDeltaPlasticStrainNormCompression(), computeDGDInv(), computeRatioPotential(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Definition at line 655 of file concretedpm2.h.
Referenced by computeDamage(), computeDamageParamCompression(), computeDamageParamTension(), computeDeltaPlasticStrainNormCompression(), computeDeltaPlasticStrainNormTension(), computeEquivalentStrain(), initDamaged(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Definition at line 648 of file concretedpm2.h.
Referenced by computeDDRDDCosTheta(), computeDFDInv(), computeDFDKappa(), computeDFDStress(), computeDRDCosTheta(), computeEquivalentStrain(), computeYieldValue(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Control parameter for the exponential softening law.
Definition at line 703 of file concretedpm2.h.
Referenced by computeDamageParamCompression(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Elastic Young's modulus.
Definition at line 694 of file concretedpm2.h.
Referenced by computeDamageParamCompression(), computeDamageParamTension(), giveRealStressVector_3d(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Type of energy strain rate dependence used if strengthRateType >0 . 0 = mod. CEB strain rate effect for strength with constant fracture energy 1 = CEB strain rate effect for strength and linear for fracture energy 2 = mod. CEB strain rate effect for strength and squared for fracture energy
Definition at line 742 of file concretedpm2.h.
Referenced by computeDamageParamCompression(), computeDamageParamTension(), initializeFrom(), restoreContext(), and saveContext().
|
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 648 of file concretedpm2.h.
Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDDuctilityMeasureDInv(), computeDFDInv(), computeDFDKappa(), computeDFDStress(), computeDGDInv(), computeDuctilityMeasure(), computeEquivalentStrain(), computeRatioPotential(), computeYieldValue(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Definition at line 648 of file concretedpm2.h.
Referenced by computeDamageParamCompression(), computeDamageParamTension(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDeltaPlasticStrainNormCompression(), computeDGDInv(), computeRatioPotential(), giveRealStressVector_3d(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Control parameter for the bilinear softening law.
Definition at line 712 of file concretedpm2.h.
Referenced by computeDamageParamTension(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Elastic shear modulus.
Definition at line 696 of file concretedpm2.h.
Referenced by computeJacobian(), computeTempKappa(), initializeFrom(), performRegularReturn(), restoreContext(), and saveContext().
|
protected |
Hardening modulus.
Definition at line 667 of file concretedpm2.h.
Referenced by restoreContext(), and saveContext().
|
protected |
Element size (to be used in fracture energy approach (crack band).
Definition at line 688 of file concretedpm2.h.
Referenced by initDamaged(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Elastic bulk modulus.
Definition at line 698 of file concretedpm2.h.
Referenced by computeJacobian(), computeTempKappa(), initializeFrom(), performRegularReturn(), restoreContext(), and saveContext().
|
protected |
Pointer for linear elastic material.
Definition at line 691 of file concretedpm2.h.
Referenced by compute3dSecantStiffness(), compute3dTangentStiffness(), computeFullJacobian(), ConcreteDPM2(), give3dMaterialStiffnessMatrix(), giveRealStressVector_3d(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Friction parameter of the yield surface.
Definition at line 682 of file concretedpm2.h.
Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDFDInv(), computeDFDKappa(), computeDFDStress(), computeDGDInv(), computeEquivalentStrain(), computeRatioPotential(), computeYieldValue(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Dilation parameter of the plastic potential.
Definition at line 685 of file concretedpm2.h.
Referenced by computeDDGDInvDKappa(), computeDGDInv(), computeRatioPotential(), restoreContext(), and saveContext().
|
protected |
Maximum number of iterations for stress return.
Definition at line 721 of file concretedpm2.h.
Referenced by computeDamageParamCompression(), computeDamageParamTension(), initializeFrom(), performRegularReturn(), restoreContext(), and saveContext().
|
protected |
Elastic poisson's ration.
Definition at line 700 of file concretedpm2.h.
Referenced by computeRatioPotential(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Type of softening function used.
Definition at line 724 of file concretedpm2.h.
Referenced by computeDamageParamTension(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Type of strength strain rate dependence used. 0 = no strain rate (default) 1 = Model Code 2010 initial branch of strain rate effect for strength 2 = Model Code 2010 initial and second branch of strain rate effect for strength
Definition at line 734 of file concretedpm2.h.
Referenced by computeDamageParamCompression(), computeDamageParamTension(), computeRateFactor(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Control parameter for the linear/bilinear softening law in tension.
Definition at line 706 of file concretedpm2.h.
Referenced by computeDamageParamTension(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Control parameter for the bilinear softening law in tension.
Definition at line 709 of file concretedpm2.h.
Referenced by computeDamageParamTension(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the hardening law of the plasticity model.
Definition at line 676 of file concretedpm2.h.
Referenced by computeHardeningOne(), computeHardeningOnePrime(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the hardening law of the plasticity model.
Definition at line 673 of file concretedpm2.h.
Referenced by computeHardeningOne(), computeHardeningOnePrime(), computeHardeningTwo(), computeHardeningTwoPrime(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
yield tolerance for the plasticity model.
Definition at line 715 of file concretedpm2.h.
Referenced by initializeFrom(), performRegularReturn(), performVertexReturn(), restoreContext(), and saveContext().
|
protected |
yield tolerance for the damage model.
Definition at line 718 of file concretedpm2.h.
Referenced by checkForUnAndReloading(), computeDamage(), computeDamageParamCompression(), computeDamageParamTension(), computeDeltaPlasticStrainNormCompression(), computeDeltaPlasticStrainNormTension(), initDamaged(), initializeFrom(), restoreContext(), and saveContext().