|
OOFEM 3.0
|
#include <concretedpm.h>
Public Member Functions | |
| ConcreteDPM (int n, Domain *d) | |
| Constructor. | |
| void | initializeFrom (InputRecord &ir) override |
| const char * | giveClassName () const override |
| const char * | giveInputRecordName () const override |
| ConcreteDPMStatus * | giveConcreteDPMStatus (GaussPoint *gp) const |
| FloatArrayF< 6 > | giveRealStressVector_3d (const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress. | |
| void | performPlasticityReturn (GaussPoint *gp, FloatArrayF< 6 > &strain) const |
| void | checkForVertexCase (double &answer, double sig, double tempKappa, GaussPoint *gp) const |
| void | performRegularReturn (FloatArrayF< 6 > &stress, GaussPoint *gp, double theta) const |
| void | performVertexReturn (FloatArrayF< 6 > &stress, double apexStress, 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 | 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 |
| FloatMatrixF< 3, 3 > | computeAMatrix (double sig, double rho, double theta, double tempKappa, double deltaLambda) const |
| FloatArrayF< 2 > | computeDGDInv (double sig, double rho, double tempKappa) const |
| double | computeRatioPotential (double sig, double tempKappa) const |
| FloatMatrixF< 2, 2 > | computeDDGDDInv (double sig, double rho, double tempKappa) const |
| FloatArrayF< 2 > | computeDDGDInvDKappa (double sig, double rho, double tempKappa) const |
| FloatArrayF< 2 > | computeDDKappaDDeltaLambdaDInv (double sig, double rho, double theta, 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 |
| double | computeTempKappa (double kappaInitial, double sigTrial, double rhoTrial, double sig) const |
| std::pair< double, double > | computeDamage (const FloatArrayF< 6 > &strain, GaussPoint *gp, TimeStep *tStep) const |
| virtual double | computeDamageParam (double kappa, GaussPoint *gp) const |
| Compute damage parameter. | |
| double | computeInverseDamage (double dam, GaussPoint *gp) const |
| Compute the damage-driving variable from given damage. | |
| virtual double | computeEquivalentStrain (const FloatArrayF< 6 > &elasticStrain, GaussPoint *gp, TimeStep *tStep) const |
| Compute equivalent strain value. | |
| double | computeDuctilityMeasureDamage (const FloatArrayF< 6 > &strain, GaussPoint *gp) const |
| Compute the ductility measure for the damage model. | |
| void | initDamaged (double kappa, const FloatArrayF< 6 > &elasticStrain, GaussPoint *gp) const |
| std::tuple< double, double, double > | computeTrialCoordinates (const FloatArrayF< 6 > &stress, GaussPoint *gp) const |
| Compute the trial 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. | |
| void | computeDSigDStress (FloatArrayF< 6 > &answer) 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. | |
| FloatArrayF< 6 > | computeDCosThetaDStress (const FloatArrayF< 6 > &stress) const |
| Computes the derivative of costheta with respect to the stress. | |
| double | computeDRDCosTheta (double theta, double ecc) const |
| Compute the derivative of R with respect to costheta. | |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const override |
| int | setIPValue (const FloatArray &value, GaussPoint *gp, InternalStateType type) override |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| void | restoreConsistency (GaussPoint *gp) 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 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 | |
| IsotropicLinearElasticMaterial | linearElasticMaterial |
| Linear elastic material. | |
| double | fc = 0. |
| double | ft = 0. |
| double | ecc = 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 | ASoft = 0. |
| Parameter of the ductilityMeasure of the damage 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. |
| Plastic multiplier of the plasticity model. | |
| double | mQ = 0. |
| The dilation parameter of the plastic potential. | |
| double | helem = 0. |
| Element size (to be used in fracture energy approach (crack band). | |
| 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 | ef = 0. |
| Hardening variable of plasticity model. | |
| double | yieldTol = 0. |
| Yield tolerance for the plasticity model. | |
| int | newtonIter = 0 |
| Maximum number of iterations for stress return. | |
| double | href = 0. |
| Stress and its deviatoric part. | |
| 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. The yield surface of the plasticity model is based on the extension of the Menetrey and Willam yield criterion. The flow rule is nonassociated. The evolution laws of the hardening variables depend on the stress state. The plasticity model describes only hardening and perfect plasticity. It is based on h effective stress. The damage parameter of the isotropic damage model is based on the total volumetric strain. An exponential softening law is implemented.
Definition at line 382 of file concretedpm.h.
| oofem::ConcreteDPM::ConcreteDPM | ( | int | n, |
| Domain * | d ) |
Constructor.
Definition at line 304 of file concretedpm.C.
References linearElasticMaterial, and oofem::StructuralMaterial::StructuralMaterial().
| void oofem::ConcreteDPM::assignStateFlag | ( | GaussPoint * | gp | ) | const |
Assign state flag.
Definition at line 1549 of file concretedpm.C.
References oofem::ConcreteDPMStatus::ConcreteDPM_Damage, oofem::ConcreteDPMStatus::ConcreteDPM_Elastic, oofem::ConcreteDPMStatus::ConcreteDPM_Plastic, oofem::ConcreteDPMStatus::ConcreteDPM_PlasticDamage, oofem::ConcreteDPMStatus::ConcreteDPM_Unloading, and giveConcreteDPMStatus().
Referenced by giveRealStressVector_3d().
| void oofem::ConcreteDPM::checkForVertexCase | ( | double & | answer, |
| 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 | The volumetric apex stress. |
| sig | The volumetric stress. |
| tempKappa | The hardening variable. |
Definition at line 746 of file concretedpm.C.
References computeHardeningOne(), fc, oofem::Material::giveStatus(), m, newtonIter, oofem::ConcreteDPMStatus::VT_Compression, oofem::ConcreteDPMStatus::VT_Regular, oofem::ConcreteDPMStatus::VT_Tension, and yieldTol.
Referenced by performPlasticityReturn().
| FloatMatrixF< 3, 3 > oofem::ConcreteDPM::computeAMatrix | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa, | ||
| double | deltaLambda ) const |
This matrix is the core of the closest point projection and collects the derivative of the flow rule and the hardening parameters.
Definition at line 1458 of file concretedpm.C.
References computeDDGDDInv(), computeDDGDInvDKappa(), computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), gM, oofem::inv(), and kM.
Referenced by performRegularReturn().
| std::pair< double, double > oofem::ConcreteDPM::computeDamage | ( | const FloatArrayF< 6 > & | strain, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Perform stress return for the damage model, i.e. if the trial stress state does not violate the plasticity surface.
| strain | Strain. |
| gp | Gauss point. |
| tStep | Time step. |
Definition at line 446 of file concretedpm.C.
References computeDamageParam(), computeEquivalentStrain(), giveConcreteDPMStatus(), href, and initDamaged().
Referenced by giveRealStressVector_3d().
|
virtual |
Compute damage parameter.
Definition at line 547 of file concretedpm.C.
References DPM_DAMAGE_TOLERANCE, ef, eM, ft, giveConcreteDPMStatus(), oofem::ConcreteDPMStatus::giveEpsLoc(), oofem::ConcreteDPMStatus::giveLe(), helem, href, and OOFEM_ERROR.
Referenced by computeDamage().
| FloatArrayF< 6 > oofem::ConcreteDPM::computeDCosThetaDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the derivative of costheta with respect to the stress.
Definition at line 1660 of file concretedpm.C.
References oofem::StructuralMaterial::computeDeviator(), computeDRhoDStress(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::from_voigt_stress_6().
| FloatMatrixF< 2, 2 > oofem::ConcreteDPM::computeDDGDDInv | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
Here, the second derivative of the plastic potential with respect to the invariants sig and rho are computed.
Definition at line 1409 of file concretedpm.C.
References computeHardeningOne(), dilationConst, fc, ft, and m.
Referenced by computeAMatrix(), and computeDDKappaDDeltaLambdaDInv().
| FloatArrayF< 2 > oofem::ConcreteDPM::computeDDGDInvDKappa | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
Here, the mixed derivative of the plastic potential with respect to the invariants and the hardening parameter are determined.
Definition at line 1373 of file concretedpm.C.
References computeHardeningOne(), computeHardeningOnePrime(), dilationConst, fc, ft, m, and mQ.
Referenced by computeAMatrix(), and computeDDKappaDDeltaLambdaDKappa().
| FloatArrayF< 2 > oofem::ConcreteDPM::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 1191 of file concretedpm.C.
References computeDDGDDInv(), computeDDuctilityMeasureDInv(), computeDGDInv(), and computeDuctilityMeasure().
Referenced by computeAMatrix().
| double oofem::ConcreteDPM::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 1224 of file concretedpm.C.
References computeDDGDInvDKappa(), computeDGDInv(), and computeDuctilityMeasure().
Referenced by computeAMatrix().
| FloatMatrixF< 6, 6 > oofem::ConcreteDPM::computeDDRhoDDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the second derivative of rho with the respect to the stress.
Definition at line 1620 of file concretedpm.C.
References oofem::StructuralMaterial::computeDeviator(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::dyad().
| FloatArrayF< 2 > oofem::ConcreteDPM::computeDDuctilityMeasureDInv | ( | double | sig, |
| double | rho, | ||
| double | theta, | ||
| double | tempKappa ) const |
Computes the first derivative of the ductility measure with respect to the invariants sig and rho based on the stress state and the hardening parameter.
Definition at line 1281 of file concretedpm.C.
References AHard, BHard, CHard, DHard, and fc.
Referenced by computeDDKappaDDeltaLambdaDInv().
| FloatArrayF< 2 > oofem::ConcreteDPM::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 1143 of file concretedpm.C.
References oofem::AL, computeHardeningOne(), ecc, fc, and m.
Referenced by performRegularReturn().
| double oofem::ConcreteDPM::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 of the stress state. |
| tempKappa | Hardening variable. |
Definition at line 1096 of file concretedpm.C.
References computeHardeningOne(), computeHardeningOnePrime(), ecc, fc, and m.
Referenced by performRegularReturn().
| FloatArrayF< 2 > oofem::ConcreteDPM::computeDGDInv | ( | double | sig, |
| double | rho, | ||
| double | tempKappa ) const |
Here, the first derivative of the plastic potential with respect to the invariants sig and rho are computed.
Definition at line 1312 of file concretedpm.C.
References computeHardeningOne(), dilationConst, fc, ft, m, and mQ.
Referenced by computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDKappaDDeltaLambda(), and performRegularReturn().
| double oofem::ConcreteDPM::computeDKappaDDeltaLambda | ( | double | sig, |
| double | rho, | ||
| double | 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 of the stress state. |
| tempKappa | Hardening variable. |
Definition at line 1173 of file concretedpm.C.
References computeDGDInv(), and computeDuctilityMeasure().
Referenced by performRegularReturn().
| double oofem::ConcreteDPM::computeDRDCosTheta | ( | double | theta, |
| double | ecc ) const |
Compute the derivative of R with respect to costheta.
Definition at line 1690 of file concretedpm.C.
References ecc.
| FloatArrayF< 6 > oofem::ConcreteDPM::computeDRhoDStress | ( | const FloatArrayF< 6 > & | stress | ) | const |
Computes the derivative of rho with respect to the stress.
Definition at line 1588 of file concretedpm.C.
References oofem::StructuralMaterial::computeDeviator(), and oofem::StructuralMaterial::computeSecondCoordinate().
Referenced by computeDCosThetaDStress().
| void oofem::ConcreteDPM::computeDSigDStress | ( | FloatArrayF< 6 > & | answer | ) | const |
Computes the derivative of sig with respect to the stress.
Definition at line 1606 of file concretedpm.C.
|
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 1253 of file concretedpm.C.
References AHard, BHard, CHard, DHard, and fc.
Referenced by computeDDKappaDDeltaLambdaDInv(), computeDDKappaDDeltaLambdaDKappa(), computeDKappaDDeltaLambda(), and computeTempKappa().
| double oofem::ConcreteDPM::computeDuctilityMeasureDamage | ( | const FloatArrayF< 6 > & | strain, |
| GaussPoint * | gp ) const |
Compute the ductility measure for the damage model.
Definition at line 665 of file concretedpm.C.
References ASoft, oofem::StructuralMaterial::computePrincipalValues(), oofem::from_voigt_strain_6(), and giveConcreteDPMStatus().
Referenced by computeEquivalentStrain().
|
virtual |
Compute equivalent strain value.
FIXME
FIXME
Definition at line 478 of file concretedpm.C.
References computeDuctilityMeasureDamage(), and giveConcreteDPMStatus().
Referenced by computeDamage().
| double oofem::ConcreteDPM::computeHardeningOne | ( | double | tempKappa | ) | const |
Compute the value of the hardening function based on the hardening variable.
| tempKappa | Hardening variable. |
Definition at line 1488 of file concretedpm.C.
References yieldHardInitial.
Referenced by checkForVertexCase(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeRatioPotential(), and computeYieldValue().
| double oofem::ConcreteDPM::computeHardeningOnePrime | ( | double | tempKappa | ) | const |
Compute the derivative of the hardening function based on the hardening parameter.
| tempKappa | Hardening variable. |
Definition at line 1502 of file concretedpm.C.
References yieldHardInitial.
Referenced by computeDDGDInvDKappa(), and computeDFDKappa().
| double oofem::ConcreteDPM::computeInverseDamage | ( | double | dam, |
| GaussPoint * | gp ) const |
Compute the damage-driving variable from given damage.
Definition at line 526 of file concretedpm.C.
References oofem::Element::computeMeanSize(), ef, eM, ft, giveConcreteDPMStatus(), oofem::GaussPoint::giveElement(), and helem.
Referenced by restoreConsistency().
| double oofem::ConcreteDPM::computeRatioPotential | ( | double | sig, |
| 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 1343 of file concretedpm.C.
References computeHardeningOne(), dilationConst, fc, ft, m, mQ, and nu.
Referenced by performVertexReturn().
| double oofem::ConcreteDPM::computeTempKappa | ( | double | kappaInitial, |
| double | sigTrial, | ||
| double | rhoTrial, | ||
| double | sig ) const |
Compute temporary kappa.
FIXME: verify
Definition at line 1051 of file concretedpm.C.
References computeDuctilityMeasure(), gM, and kM.
Referenced by performVertexReturn().
| std::tuple< double, double, double > oofem::ConcreteDPM::computeTrialCoordinates | ( | const FloatArrayF< 6 > & | stress, |
| GaussPoint * | gp ) const |
Compute the trial coordinates.
Definition at line 1534 of file concretedpm.C.
References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), oofem::StructuralMaterial::computeSecondCoordinate(), and oofem::StructuralMaterial::computeThirdCoordinate().
Referenced by performPlasticityReturn().
| double oofem::ConcreteDPM::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 1068 of file concretedpm.C.
References computeHardeningOne(), 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 1788 of file concretedpm.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 1514 of file concretedpm.C.
References giveConcreteDPMStatus(), and linearElasticMaterial.
Referenced by restoreConsistency().
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 479 of file concretedpm.h.
|
inline |
Definition at line 482 of file concretedpm.h.
References oofem::Material::giveStatus().
Referenced by assignStateFlag(), computeDamage(), computeDamageParam(), computeDuctilityMeasureDamage(), computeEquivalentStrain(), computeInverseDamage(), give3dMaterialStiffnessMatrix(), giveIPValue(), giveRealStressVector_3d(), initDamaged(), restoreConsistency(), and setIPValue().
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Definition at line 480 of file concretedpm.h.
References _IFT_ConcreteDPM_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 1743 of file concretedpm.C.
References oofem::FloatArray::at(), giveConcreteDPMStatus(), oofem::StructuralMaterial::giveIPValue(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
overridevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
is this right? or was is supposed to use the old converged values?
Reimplemented from oofem::StructuralMaterial.
Definition at line 382 of file concretedpm.C.
References oofem::StructuralMaterial::applyElasticStiffness(), assignStateFlag(), computeDamage(), oofem::StructuralMaterial::computeStressIndependentStrainVector_3d(), oofem::FloatArray::dotProduct(), E, eM, ft, ft_strength, oofem::Material::give(), giveConcreteDPMStatus(), oofem::FloatArray::giveSize(), nu, and performPlasticityReturn().
| void oofem::ConcreteDPM::initDamaged | ( | double | kappa, |
| const FloatArrayF< 6 > & | elasticStrain, | ||
| GaussPoint * | gp ) const |
Initialize the characteristic length, if damage is not yet activated.
Definition at line 618 of file concretedpm.C.
References oofem::FloatArrayF< N >::at(), oofem::Element::computeMeanSize(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::from_voigt_strain_6(), oofem::Element::giveCharacteristicLength(), giveConcreteDPMStatus(), oofem::GaussPoint::giveElement(), and helem.
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 311 of file concretedpm.C.
References _IFT_ConcreteDPM_ahard, _IFT_ConcreteDPM_asoft, _IFT_ConcreteDPM_bhard, _IFT_ConcreteDPM_chard, _IFT_ConcreteDPM_dhard, _IFT_ConcreteDPM_dilation, _IFT_ConcreteDPM_ecc, _IFT_ConcreteDPM_fc, _IFT_ConcreteDPM_ft, _IFT_ConcreteDPM_gf, _IFT_ConcreteDPM_helem, _IFT_ConcreteDPM_href, _IFT_ConcreteDPM_kinit, _IFT_ConcreteDPM_newtoniter, _IFT_ConcreteDPM_wf, _IFT_ConcreteDPM_yieldtol, _IFT_IsotropicLinearElasticMaterial_e, _IFT_IsotropicLinearElasticMaterial_n, AHard, ASoft, BHard, CHard, DHard, dilationConst, ecc, ef, eM, fc, fc_strength, ft, ft_strength, gM, oofem::InputRecord::hasField(), helem, href, oofem::StructuralMaterial::initializeFrom(), IR_GIVE_FIELD, IR_GIVE_OPTIONAL_FIELD, kM, linearElasticMaterial, m, newtonIter, nu, oofem::Material::propertyDictionary, yieldHardInitial, and yieldTol.
|
inlineoverridevirtual |
Returns true if stiffness matrix of receiver is symmetric Default implementation returns true.
Reimplemented from oofem::Material.
Definition at line 738 of file concretedpm.h.
| void oofem::ConcreteDPM::performPlasticityReturn | ( | GaussPoint * | gp, |
| FloatArrayF< 6 > & | strain ) const |
| gp | Gauss point. |
| strain | Strain vector of this Gauss point. |
FIXME
FIXME
Definition at line 692 of file concretedpm.C.
References oofem::StructuralMaterial::applyElasticCompliance(), oofem::StructuralMaterial::applyElasticStiffness(), checkForVertexCase(), computeTrialCoordinates(), computeYieldValue(), eM, oofem::Material::giveStatus(), oofem::ConcreteDPMStatus::initTempStatus(), nu, performRegularReturn(), performVertexReturn(), oofem::ConcreteDPMStatus::VT_Compression, oofem::ConcreteDPMStatus::VT_Regular, oofem::ConcreteDPMStatus::VT_Tension, and yieldTol.
Referenced by giveRealStressVector_3d().
| void oofem::ConcreteDPM::performRegularReturn | ( | FloatArrayF< 6 > & | stress, |
| 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. |
| gp | Gauss point. |
| theta | Lode angle of the stress state. |
FIXME: verify
FIXME: verify
FIXME: verify
FIXME
Definition at line 792 of file concretedpm.C.
References computeAMatrix(), oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeDKappaDDeltaLambda(), oofem::StructuralMaterial::computePrincipalValDir(), oofem::StructuralMaterial::computeSecondCoordinate(), computeYieldValue(), oofem::dot(), oofem::from_voigt_stress_6(), oofem::Material::giveStatus(), gM, kM, M_PI, newtonIter, OOFEM_ERROR, oofem::StructuralMaterial::transformStressVectorTo(), and yieldTol.
Referenced by performPlasticityReturn().
| void oofem::ConcreteDPM::performVertexReturn | ( | FloatArrayF< 6 > & | stress, |
| double | apexStress, | ||
| 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. |
| gp | Gauss point. |
FIXME verify
Definition at line 958 of file concretedpm.C.
References oofem::StructuralMaterial::computeDeviatoricVolumetricSplit(), computeRatioPotential(), oofem::StructuralMaterial::computeSecondCoordinate(), computeTempKappa(), computeYieldValue(), ft, oofem::FloatArrayF< N >::giveSize(), oofem::Material::giveStatus(), oofem::ConcreteDPMStatus::VT_Compression, oofem::ConcreteDPMStatus::VT_Regular, oofem::ConcreteDPMStatus::VT_Tension, and yieldTol.
Referenced by performPlasticityReturn().
|
overridevirtual |
Restores consistency of the status, i.e., computes or corrects the values of certain status variables such that the state is admissible. For instance, if the initial values of some internal variables are read from a file, other internal variables are adjusted accordingly.
Reimplemented from oofem::Material.
Definition at line 1708 of file concretedpm.C.
References computeInverseDamage(), give3dMaterialStiffnessMatrix(), giveConcreteDPMStatus(), and oofem::solve().
|
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 1868 of file concretedpm.C.
References AHard, ASoft, BHard, CHard, oofem::CIO_IOERR, CM_Definition, DHard, dilationConst, ecc, ef, eM, fc, ft, gM, helem, href, kM, linearElasticMaterial, m, mQ, newtonIter, nu, oofem::DataStream::read(), oofem::FEMComponent::restoreContext(), THROW_CIOERR, yieldHardInitial, and yieldTol.
|
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 1794 of file concretedpm.C.
References AHard, ASoft, BHard, CHard, oofem::CIO_IOERR, CM_Definition, DHard, dilationConst, ecc, ef, eM, fc, ft, gM, helem, href, kM, linearElasticMaterial, m, mQ, newtonIter, nu, oofem::FEMComponent::saveContext(), THROW_CIOERR, oofem::DataStream::write(), yieldHardInitial, and yieldTol.
|
overridevirtual |
Sets the value of a certain variable at a given integration point to the given value.
| value | Contains the value(s) to be set (in reduced form). |
| gp | Integration point. |
| type | Determines the type of internal variable. |
| type | Determines the type of internal variable. |
Reimplemented from oofem::Material.
Definition at line 1732 of file concretedpm.C.
References giveConcreteDPMStatus(), and oofem::StructuralMaterial::setIPValue().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 397 of file concretedpm.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the damage model.
Definition at line 406 of file concretedpm.h.
Referenced by computeDuctilityMeasureDamage(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 399 of file concretedpm.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 401 of file concretedpm.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Parameter of the ductilityMeasure of the plasticity model.
Definition at line 403 of file concretedpm.h.
Referenced by computeDDuctilityMeasureDInv(), computeDuctilityMeasure(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Control parameter for te volumetric plastic flow of the plastic potential.
Definition at line 412 of file concretedpm.h.
Referenced by computeDDGDDInv(), computeDDGDInvDKappa(), computeDGDInv(), computeRatioPotential(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Definition at line 394 of file concretedpm.h.
Referenced by computeDFDInv(), computeDFDKappa(), computeDRDCosTheta(), computeYieldValue(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Hardening variable of plasticity model.
Hardening variable of damage model. Damage variable of damage model. Control parameter for the exponential softening law.
Definition at line 456 of file concretedpm.h.
Referenced by computeDamageParam(), computeInverseDamage(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Elastic Young's modulus.
Definition at line 435 of file concretedpm.h.
Referenced by computeDamageParam(), computeInverseDamage(), giveRealStressVector_3d(), initializeFrom(), performPlasticityReturn(), 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 394 of file concretedpm.h.
Referenced by checkForVertexCase(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDDuctilityMeasureDInv(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeDuctilityMeasure(), computeRatioPotential(), computeYieldValue(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Definition at line 394 of file concretedpm.h.
Referenced by computeDamageParam(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDGDInv(), computeInverseDamage(), computeRatioPotential(), giveRealStressVector_3d(), initializeFrom(), performVertexReturn(), restoreContext(), and saveContext().
|
protected |
Elastic shear modulus.
Definition at line 437 of file concretedpm.h.
Referenced by computeAMatrix(), computeTempKappa(), initializeFrom(), performRegularReturn(), restoreContext(), and saveContext().
|
protected |
Element size (to be used in fracture energy approach (crack band).
Definition at line 432 of file concretedpm.h.
Referenced by computeDamageParam(), computeInverseDamage(), initDamaged(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Stress and its deviatoric part.
Material parameter of the size-dependent adjustment (reference element size)
Definition at line 470 of file concretedpm.h.
Referenced by computeDamage(), computeDamageParam(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Elastic bulk modulus.
Definition at line 439 of file concretedpm.h.
Referenced by computeAMatrix(), computeTempKappa(), initializeFrom(), performRegularReturn(), restoreContext(), and saveContext().
|
protected |
Linear elastic material.
Definition at line 389 of file concretedpm.h.
Referenced by ConcreteDPM(), give3dMaterialStiffnessMatrix(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Plastic multiplier of the plasticity model.
the volumetric stress. The length of the deviatoric stress. The lode angle of the trial stress. The friction parameter of the yield surface.
Definition at line 426 of file concretedpm.h.
Referenced by checkForVertexCase(), computeDDGDDInv(), computeDDGDInvDKappa(), computeDFDInv(), computeDFDKappa(), computeDGDInv(), computeRatioPotential(), computeYieldValue(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
The dilation parameter of the plastic potential.
Definition at line 429 of file concretedpm.h.
Referenced by computeDDGDInvDKappa(), computeDGDInv(), computeRatioPotential(), restoreContext(), and saveContext().
|
protected |
Maximum number of iterations for stress return.
Definition at line 462 of file concretedpm.h.
Referenced by checkForVertexCase(), initializeFrom(), performRegularReturn(), restoreContext(), and saveContext().
|
protected |
Elastic Poisson's ration.
Definition at line 441 of file concretedpm.h.
Referenced by computeRatioPotential(), giveRealStressVector_3d(), initializeFrom(), performPlasticityReturn(), restoreContext(), and saveContext().
|
protected |
Parameter of the hardening law of the plasticity model.
Definition at line 409 of file concretedpm.h.
Referenced by computeHardeningOne(), computeHardeningOnePrime(), initializeFrom(), restoreContext(), and saveContext().
|
protected |
Yield tolerance for the plasticity model.
Definition at line 459 of file concretedpm.h.
Referenced by checkForVertexCase(), initializeFrom(), performPlasticityReturn(), performRegularReturn(), performVertexReturn(), restoreContext(), and saveContext().