|
OOFEM 3.0
|
#include <rcm2.h>
Public Member Functions | |
| RCM2Material (int n, Domain *d) | |
| bool | hasMaterialModeCapability (MaterialMode mode) const override |
| const char * | giveClassName () const override |
| void | initializeFrom (InputRecord &ir) override |
| double | give (int aProperty, GaussPoint *gp) const override |
| LinearElasticMaterial * | giveLinearElasticMaterial () const |
| FloatMatrixF< 6, 6 > | give3dMaterialStiffnessMatrix (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const override |
| void | giveRealStressVector (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep) 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. | |
| FloatArrayF< 4 > | giveRealStressVector_PlaneStrain (const FloatArrayF< 4 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_3d. | |
| FloatArrayF< 3 > | giveRealStressVector_PlaneStress (const FloatArrayF< 3 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| FloatArrayF< 1 > | giveRealStressVector_1d (const FloatArrayF< 1 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| FloatArrayF< 2 > | giveRealStressVector_2dBeamLayer (const FloatArrayF< 2 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| FloatArrayF< 5 > | giveRealStressVector_PlateLayer (const FloatArrayF< 5 > &strain, GaussPoint *gp, TimeStep *tStep) const override |
| Default implementation relies on giveRealStressVector_StressControl. | |
| int | giveIPValue (FloatArray &answer, GaussPoint *gp, InternalStateType type, TimeStep *tStep) override |
| std::unique_ptr< MaterialStatus > | CreateStatus (GaussPoint *gp) const override |
| FloatArrayF< 6 > | giveThermalDilatationVector (GaussPoint *gp, TimeStep *tStep) const 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 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< 2 > | giveRealStressVector_Warping (const FloatArrayF< 2 > &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) |
| 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< 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< 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 > | give1dStressStiffnessMatrix_dPdF (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | give1dStressStiffMtrx_dCde (FloatMatrix &answer, MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) |
| 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 bool | isCharacteristicMtrxSymmetric (MatResponseMode rMode) const |
| virtual double | giveCharacteristicValue (MatResponseMode type, GaussPoint *gp, TimeStep *tStep) const |
| Returns characteristic value of the receiver. | |
| 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. | |
| virtual const char * | giveInputRecordName () const =0 |
| 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 | |
| virtual void | checkForNewActiveCracks (IntArray &answer, GaussPoint *gp, const FloatArray &, const FloatArray &, FloatArray &, const FloatArray &) const |
| virtual void | updateCrackStatus (GaussPoint *gp, const FloatArray &crackStrain) const |
| virtual void | checkIfClosedCracks (GaussPoint *gp, FloatArray &crackStrainVector, IntArray &) const |
| virtual int | checkSizeLimit (GaussPoint *gp, double) const |
| virtual double | giveNormalCrackingStress (GaussPoint *gp, double eps_cr, int i) const =0 |
| virtual double | giveMinCrackStrainsForFullyOpenCrack (GaussPoint *gp, int i) const =0 |
| virtual double | computeStrength (GaussPoint *gp, double) const =0 |
| virtual void | updateStatusForNewCrack (GaussPoint *, int, double) const |
| virtual double | giveCharacteristicElementLength (GaussPoint *gp, const FloatArray &crackPlaneNormal) const |
| virtual double | giveCrackingModulus (MatResponseMode rMode, GaussPoint *gp, double effStrain, int i) const |
| virtual void | giveMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode, GaussPoint *gp, TimeStep *tStep) const |
| void | giveCrackedStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const |
| virtual void | giveEffectiveMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const |
| void | giveRealPrincipalStressVector3d (FloatArray &answer, GaussPoint *, FloatArray &, FloatMatrix &, TimeStep *) const |
| void | giveNormalElasticStiffnessMatrix (FloatMatrix &answer, bool reduce, MatResponseMode, GaussPoint *, TimeStep *tStep, const FloatMatrix &) const |
| void | updateActiveCrackMap (GaussPoint *gp, const IntArray *activatedCracks=NULL) const |
| double | giveResidualStrength () |
| FloatMatrixF< 3, 3 > | givePlaneStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 4, 4 > | givePlaneStrainStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 1, 1 > | give1dStressStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 2, 2 > | give2dBeamLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
| FloatMatrixF< 5, 5 > | givePlateLayerStiffMtrx (MatResponseMode mmode, GaussPoint *gp, TimeStep *tStep) const override |
Protected Attributes | |
| LinearElasticMaterial * | linearElasticMaterial = nullptr |
| double | Gf = 0. |
| double | Ft = 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 implements a Rotating Crack Model for fracture in smeared fashion ( only material stiffness modification is required, no changes in mesh topology) coupled with plastic behaviour. In this class we follow a paper written by R. de Borst "Smeared Cracking, plasticity, creep - Unified approach " which allows cracking phenomena to be easily incorporated into existing code for plastic materials.
A model is based on Fracture Energy criterion, with (linear softening phenomena and with shear retention factor), controlled un-re loading.
| oofem::RCM2Material::RCM2Material | ( | int | n, |
| Domain * | d ) |
Definition at line 48 of file rcm2.C.
References oofem::StructuralMaterial::StructuralMaterial().
Referenced by oofem::Concrete3::Concrete3(), give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), give3dMaterialStiffnessMatrix(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), givePlateLayerStiffMtrx(), giveRealStressVector_1d(), giveRealStressVector_2dBeamLayer(), giveRealStressVector_3d(), giveRealStressVector_PlaneStrain(), giveRealStressVector_PlaneStress(), giveRealStressVector_PlateLayer(), oofem::RCSDEMaterial::RCSDEMaterial(), and oofem::RCSDMaterial::RCSDMaterial().
|
protectedvirtual |
Definition at line 332 of file rcm2.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::FloatArray::clear(), computeStrength(), oofem::IntArray::contains(), giveCharacteristicElementLength(), oofem::RCM2MaterialStatus::giveCrackMap(), oofem::Material::giveStatus(), oofem::RCM2MaterialStatus::giveTempCrackDirs(), oofem::RCM2MaterialStatus::giveTempMaxCrackStrain(), oofem::IntArray::resize(), and updateStatusForNewCrack().
Referenced by giveRealPrincipalStressVector3d().
|
protectedvirtual |
Definition at line 500 of file rcm2.C.
References oofem::RCM2MaterialStatus::giveCrackMap(), oofem::Material::giveStatus(), oofem::RCM2MaterialStatus::letCrackMapBe(), and updateActiveCrackMap().
Referenced by giveRealPrincipalStressVector3d().
|
inlineprotectedvirtual |
Reimplemented in oofem::Concrete3, oofem::RCSDEMaterial, and oofem::RCSDMaterial.
|
protectedpure virtual |
Implemented in oofem::Concrete3, oofem::RCSDEMaterial, oofem::RCSDMaterial, and oofem::RCSDNLMaterial.
Referenced by checkForNewActiveCracks().
|
inlineoverridevirtual |
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.
Reimplemented in oofem::RCSDEMaterial, oofem::RCSDMaterial, and oofem::RCSDNLMaterial.
|
overridevirtual |
Returns the value of material property 'aProperty'. Property must be identified by unique int id. Integration point also passed to allow for materials with spatially varying properties
| aProperty | ID of property requested. |
| gp | Integration point, |
Reimplemented from oofem::Material.
Reimplemented in oofem::RCSDEMaterial, and oofem::RCSDMaterial.
Definition at line 792 of file rcm2.C.
References Ft, Gf, linearElasticMaterial, OOFEM_ERROR, oofem::Material::propertyDictionary, pscm_Ee, pscm_Ft, pscm_G, and pscm_Gf.
Referenced by oofem::Concrete3::checkSizeLimit(), oofem::Concrete3::computeStrength(), oofem::Concrete3::giveCrackingModulus(), giveEffectiveMaterialStiffnessMatrix(), oofem::Concrete3::giveMinCrackStrainsForFullyOpenCrack(), oofem::Concrete3::giveNormalCrackingStress(), and giveRealPrincipalStressVector3d().
|
overrideprotectedvirtual |
Method for computing 1d stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 1d stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
FIXME: Temporary const-cast until other routines have been made const.
Reimplemented from oofem::StructuralMaterial.
Definition at line 922 of file rcm2.C.
References giveMaterialStiffnessMatrix(), and RCM2Material().
|
overrideprotectedvirtual |
Method for computing 2d beam layer stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d beam layer stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.
| answer | Stiffness matrix. |
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
FIXME: Temporary const-cast until other routines have been made const.
Reimplemented from oofem::StructuralMaterial.
Definition at line 938 of file rcm2.C.
References giveMaterialStiffnessMatrix(), and RCM2Material().
|
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 875 of file rcm2.C.
References giveMaterialStiffnessMatrix(), and RCM2Material().
|
protectedvirtual |
Reimplemented in oofem::RCSDNLMaterial.
Definition at line 423 of file rcm2.C.
References oofem::Element::giveCharacteristicLength(), and oofem::GaussPoint::giveElement().
Referenced by checkForNewActiveCracks().
|
inlineoverridevirtual |
Implements oofem::FEMComponent.
Reimplemented in oofem::RCSDEMaterial, oofem::RCSDMaterial, and oofem::RCSDNLMaterial.
|
protected |
Definition at line 706 of file rcm2.C.
References oofem::IntArray::at(), giveCrackingModulus(), oofem::RCM2MaterialStatus::giveCrackMap(), oofem::RCM2MaterialStatus::giveCrackStrain(), oofem::RCM2MaterialStatus::giveNumberOfTempActiveCracks(), oofem::Material::giveStatus(), and oofem::IntArray::resize().
Referenced by giveEffectiveMaterialStiffnessMatrix(), and giveRealPrincipalStressVector3d().
|
inlineprotectedvirtual |
Reimplemented in oofem::Concrete3, oofem::RCSDEMaterial, and oofem::RCSDMaterial.
Definition at line 257 of file rcm2.h.
Referenced by giveCrackedStiffnessMatrix().
|
protectedvirtual |
Reimplemented in oofem::RCSDEMaterial, and oofem::RCSDMaterial.
Definition at line 594 of file rcm2.C.
References oofem::FloatMatrix::assemble(), oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::FloatMatrix::beInverseOf(), oofem::IntArray::findFirstIndexOf(), oofem::RCM2MaterialStatus::getPrincipalStrainVector(), oofem::RCM2MaterialStatus::getPrincipalStressVector(), give(), giveCrackedStiffnessMatrix(), giveLinearElasticMaterial(), giveNormalElasticStiffnessMatrix(), oofem::RCM2MaterialStatus::giveNumberOfTempActiveCracks(), oofem::IntArray::giveSize(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStiffnessMatrix(), oofem::StructuralMaterial::giveStressVectorTranformationMtrx(), oofem::RCM2MaterialStatus::giveTempCrackDirs(), oofem::RCM2MaterialStatus::isCrackActive(), pscm_G, rcm2_BIGNUMBER, rcm_SMALL_STRAIN, oofem::FloatMatrix::resize(), oofem::FloatMatrix::rotatedWith(), oofem::StructuralMaterial::StructuralMaterial(), updateActiveCrackMap(), and oofem::FloatMatrix::zero().
Referenced by giveMaterialStiffnessMatrix().
|
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 830 of file rcm2.C.
References oofem::FloatArray::at(), oofem::FloatMatrix::at(), oofem::IntArray::at(), oofem::RCM2MaterialStatus::giveAlreadyCrack(), oofem::RCM2MaterialStatus::giveCharLength(), oofem::RCM2MaterialStatus::giveCrackDirs(), oofem::RCM2MaterialStatus::giveCrackStatus(), oofem::RCM2MaterialStatus::giveCrackStrain(), oofem::RCM2MaterialStatus::giveNumberOfActiveCracks(), oofem::Material::giveStatus(), oofem::max(), oofem::FloatArray::resize(), and oofem::FloatArray::zero().
|
inline |
Definition at line 193 of file rcm2.h.
References linearElasticMaterial.
Referenced by giveEffectiveMaterialStiffnessMatrix(), oofem::RCSDEMaterial::giveEffectiveMaterialStiffnessMatrix(), oofem::RCSDMaterial::giveEffectiveMaterialStiffnessMatrix(), giveNormalElasticStiffnessMatrix(), initializeFrom(), and oofem::RCSDNLMaterial::initializeFrom().
|
protectedvirtual |
Definition at line 83 of file rcm2.C.
References giveEffectiveMaterialStiffnessMatrix().
Referenced by give1dStressStiffMtrx(), give2dBeamLayerStiffMtrx(), give3dMaterialStiffnessMatrix(), givePlaneStrainStiffMtrx(), givePlaneStressStiffMtrx(), givePlateLayerStiffMtrx(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), and oofem::RCSDNLMaterial::giveRealStressVector().
|
protectedpure virtual |
Implemented in oofem::Concrete3, oofem::RCSDEMaterial, oofem::RCSDMaterial, and oofem::RCSDNLMaterial.
Referenced by updateCrackStatus().
|
protectedpure virtual |
Implemented in oofem::Concrete3, oofem::RCSDEMaterial, and oofem::RCSDMaterial.
Referenced by giveRealPrincipalStressVector3d().
|
protected |
Definition at line 532 of file rcm2.C.
References oofem::FloatMatrix::at(), oofem::IntArray::contains(), oofem::IntArray::findFirstIndexOf(), oofem::StructuralMaterial::giveFullSymMatrixForm(), giveLinearElasticMaterial(), oofem::StructuralMaterial::giveStiffnessMatrix(), and oofem::StructuralMaterial::StructuralMaterial().
Referenced by giveEffectiveMaterialStiffnessMatrix(), and giveRealPrincipalStressVector3d().
|
overrideprotectedvirtual |
Method for computing plane strain stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane strain stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix. Note: as already described, if zero strain component is imposed (Plane strain, ..) this condition must be taken into account in geometrical relations, and corresponding component has to be included in reduced vector. (So plane strain conditions are \( \epsilon_z = \gamma_{xz} = \gamma_{yz} = 0 \), but relations for \( \epsilon_z\) and \(\sigma_z\) are included).
| answer | Stiffness matrix. |
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
FIXME: Temporary const-cast until other routines have been made const.
Reimplemented from oofem::StructuralMaterial.
Definition at line 904 of file rcm2.C.
References giveMaterialStiffnessMatrix(), and RCM2Material().
|
overrideprotectedvirtual |
Method for computing plane stress stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to plane stress stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.
| answer | Stiffness matrix. |
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Definition at line 886 of file rcm2.C.
References giveMaterialStiffnessMatrix(), and RCM2Material().
|
overrideprotectedvirtual |
Method for computing 2d plate layer stiffness matrix of receiver. Default implementation computes 3d stiffness matrix using give3dMaterialStiffnessMatrix and reduces it to 2d plate layer stiffness using reduce method described above. However, this reduction is quite time consuming and if it is possible, it is recommended to overload this method and provide direct method for computing particular stiffness matrix.
| answer | Stiffness matrix. |
| mmode | Material response mode. |
| gp | Integration point, which load history is used. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
FIXME: Temporary const-cast until other routines have been made const.
Reimplemented from oofem::StructuralMaterial.
Definition at line 957 of file rcm2.C.
References giveMaterialStiffnessMatrix(), and RCM2Material().
|
protected |
Definition at line 150 of file rcm2.C.
References oofem::FloatMatrix::add(), oofem::FloatArray::assemble(), oofem::FloatMatrix::assemble(), oofem::FloatArray::at(), oofem::IntArray::at(), oofem::FloatArray::beDifferenceOf(), oofem::FloatArray::beProductOf(), checkForNewActiveCracks(), checkIfClosedCracks(), oofem::FloatArray::clear(), give(), oofem::RCM2MaterialStatus::giveCrackDirs(), giveCrackedStiffnessMatrix(), oofem::RCM2MaterialStatus::giveCrackMap(), oofem::RCM2MaterialStatus::giveCrackStrainVector(), giveNormalCrackingStress(), giveNormalElasticStiffnessMatrix(), oofem::RCM2MaterialStatus::giveNumberOfTempActiveCracks(), oofem::RCM2MaterialStatus::givePrevPrincStrainVector(), oofem::FloatArray::giveSize(), oofem::IntArray::giveSize(), oofem::Material::giveStatus(), oofem::RCM2MaterialStatus::letCrackStrainVectorBe(), oofem::RCM2MaterialStatus::letPrincipalStrainVectorBe(), oofem::RCM2MaterialStatus::letPrincipalStressVectorBe(), oofem::RCM2MaterialStatus::letTempCrackDirsBe(), oofem::IntArray::maximum(), OOFEM_ERROR, pscm_Ft, oofem::FloatArray::resize(), oofem::FloatMatrix::resize(), oofem::FloatMatrix::solveForRhs(), oofem::StructuralMaterial::sortPrincDirAndValCloseTo(), oofem::FloatArray::subtract(), updateActiveCrackMap(), oofem::FloatArray::zero(), and oofem::FloatMatrix::zero().
Referenced by giveRealStressVector(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), and oofem::RCSDNLMaterial::giveRealStressVector().
|
overridevirtual |
Computes the real stress vector for given total strain and integration point. The total strain is defined as strain computed directly from displacement field at given time. The stress independent parts (temperature, eigenstrains) are subtracted in constitutive driver. The service should use previously reached equilibrium history variables. Also it should update temporary history variables in status according to newly reached state. The temporary history variables are moved into equilibrium ones after global structure equilibrium has been reached by iteration process.
| answer | Stress vector in reduced form. For large deformations it is treated as the second Piola-Kirchoff stress. |
| gp | Integration point. |
| reducedStrain | Strain vector in reduced form. For large deformations it is treated as the Green-Lagrange strain. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Reimplemented from oofem::StructuralMaterial.
Reimplemented in oofem::RCSDEMaterial, oofem::RCSDMaterial, and oofem::RCSDNLMaterial.
Definition at line 104 of file rcm2.C.
References oofem::StructuralMaterial::computePrincipalValDir(), oofem::RCM2MaterialStatus::giveCrackStrainVector(), giveRealPrincipalStressVector3d(), oofem::Material::giveStatus(), oofem::StructuralMaterial::giveStressDependentPartOfStrainVector(), oofem::RCM2MaterialStatus::giveTempCrackDirs(), oofem::Material::initTempStatus(), oofem::StructuralMaterialStatus::letTempStrainVectorBe(), oofem::StructuralMaterialStatus::letTempStressVectorBe(), oofem::principal_strain, oofem::FloatArray::resizeWithValues(), oofem::StructuralMaterial::transformStressVectorTo(), and updateCrackStatus().
Referenced by giveRealStressVector_1d(), giveRealStressVector_2dBeamLayer(), giveRealStressVector_3d(), giveRealStressVector_PlaneStrain(), giveRealStressVector_PlaneStress(), and giveRealStressVector_PlateLayer().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 218 of file rcm2.h.
References giveRealStressVector(), and RCM2Material().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 224 of file rcm2.h.
References giveRealStressVector(), and RCM2Material().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector for second Piola-Kirchoff stress.
Reimplemented from oofem::StructuralMaterial.
Definition at line 200 of file rcm2.h.
References giveRealStressVector(), and RCM2Material().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_3d.
Reimplemented from oofem::StructuralMaterial.
Definition at line 206 of file rcm2.h.
References giveRealStressVector(), and RCM2Material().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 212 of file rcm2.h.
References giveRealStressVector(), and RCM2Material().
|
inlineoverridevirtual |
Default implementation relies on giveRealStressVector_StressControl.
Reimplemented from oofem::StructuralMaterial.
Definition at line 230 of file rcm2.h.
References giveRealStressVector(), and RCM2Material().
|
inlineprotected |
|
inlineoverridevirtual |
Returns a vector of coefficients of thermal dilatation in direction of each material principal (local) axis.
| answer | Vector of thermal dilatation coefficients. |
| 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 241 of file rcm2.h.
References linearElasticMaterial.
|
overridevirtual |
Tests if material supports material mode.
| mode | Required material mode. |
Reimplemented from oofem::Material.
|
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::Material.
Reimplemented in oofem::RCSDEMaterial, oofem::RCSDMaterial, and oofem::RCSDNLMaterial.
Definition at line 783 of file rcm2.C.
References _IFT_RCM2Material_ft, _IFT_RCM2Material_gf, Ft, Gf, giveLinearElasticMaterial(), and IR_GIVE_FIELD.
|
protected |
Definition at line 753 of file rcm2.C.
References oofem::IntArray::at(), oofem::RCM2MaterialStatus::giveCrackMap(), oofem::Material::giveStatus(), oofem::RCM2MaterialStatus::isCrackActive(), and oofem::RCM2MaterialStatus::letCrackMapBe().
Referenced by checkIfClosedCracks(), giveEffectiveMaterialStiffnessMatrix(), and giveRealPrincipalStressVector3d().
|
protectedvirtual |
Definition at line 429 of file rcm2.C.
References oofem::IntArray::at(), oofem::IntArray::contains(), oofem::RCM2MaterialStatus::giveCrackMap(), giveMinCrackStrainsForFullyOpenCrack(), oofem::Material::giveStatus(), oofem::RCM2MaterialStatus::giveTempCrackStatus(), oofem::RCM2MaterialStatus::giveTempMaxCrackStrain(), pscm_CLOSED, pscm_NONE, pscm_OPEN, pscm_SOFTENING, pscm_UNLOADING, oofem::RCM2MaterialStatus::setTempCrackStatus(), and oofem::RCM2MaterialStatus::setTempMaxCrackStrain().
Referenced by giveRealStressVector(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), and oofem::RCSDNLMaterial::giveRealStressVector().
|
protectedvirtual |
Definition at line 406 of file rcm2.C.
References oofem::Material::giveStatus(), OOFEM_ERROR, and oofem::RCM2MaterialStatus::setCharLength().
Referenced by checkForNewActiveCracks().
|
protected |
Definition at line 179 of file rcm2.h.
Referenced by oofem::Concrete3::checkSizeLimit(), oofem::RCSDEMaterial::checkSizeLimit(), oofem::RCSDMaterial::checkSizeLimit(), oofem::Concrete3::computeStrength(), oofem::RCSDEMaterial::computeStrength(), oofem::RCSDMaterial::computeStrength(), oofem::RCSDNLMaterial::computeStrength(), give(), oofem::Concrete3::giveCrackingModulus(), oofem::RCSDEMaterial::giveCrackingModulus(), oofem::RCSDMaterial::giveCrackingModulus(), oofem::Concrete3::giveMinCrackStrainsForFullyOpenCrack(), oofem::RCSDMaterial::giveMinCrackStrainsForFullyOpenCrack(), oofem::Concrete3::giveNormalCrackingStress(), oofem::RCSDEMaterial::giveNormalCrackingStress(), oofem::RCSDMaterial::giveNormalCrackingStress(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::RCSDNLMaterial::giveRealStressVector(), giveResidualStrength(), initializeFrom(), and oofem::RCSDNLMaterial::initializeFrom().
|
protected |
Definition at line 179 of file rcm2.h.
Referenced by oofem::Concrete3::checkSizeLimit(), oofem::RCSDEMaterial::checkSizeLimit(), oofem::RCSDMaterial::checkSizeLimit(), oofem::Concrete3::computeStrength(), oofem::RCSDEMaterial::computeStrength(), oofem::RCSDMaterial::computeStrength(), give(), oofem::Concrete3::giveCrackingModulus(), oofem::RCSDEMaterial::giveCrackingModulus(), oofem::Concrete3::giveMinCrackStrainsForFullyOpenCrack(), oofem::RCSDMaterial::giveMinCrackStrainsForFullyOpenCrack(), oofem::Concrete3::giveNormalCrackingStress(), oofem::RCSDEMaterial::giveNormalCrackingStress(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDNLMaterial::giveRealStressVector(), initializeFrom(), and oofem::RCSDNLMaterial::initializeFrom().
|
protected |
Definition at line 178 of file rcm2.h.
Referenced by oofem::RCSDEMaterial::computeCurrEquivStrain(), oofem::RCSDMaterial::computeCurrEquivStrain(), oofem::Concrete3::Concrete3(), give(), giveLinearElasticMaterial(), oofem::RCSDEMaterial::giveRealStressVector(), oofem::RCSDMaterial::giveRealStressVector(), oofem::RCSDNLMaterial::giveRealStressVector(), giveThermalDilatationVector(), oofem::RCSDEMaterial::RCSDEMaterial(), oofem::RCSDMaterial::RCSDMaterial(), oofem::Concrete3::~Concrete3(), oofem::RCSDEMaterial::~RCSDEMaterial(), and oofem::RCSDMaterial::~RCSDMaterial().