|
OOFEM 3.0
|
#include <structuralcrosssection.h>
Public Member Functions | |
| StructuralCrossSection (int n, Domain *d) | |
| FloatArray | giveRealStresses (const FloatArray &reducedStrain, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 6 > | giveRealStress_3d (const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 6 > | giveRealStress_3dDegeneratedShell (const FloatArrayF< 6 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 4 > | giveRealStress_PlaneStrain (const FloatArrayF< 4 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 3 > | giveRealStress_PlaneStress (const FloatArrayF< 3 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 1 > | giveRealStress_1d (const FloatArrayF< 1 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 2 > | giveRealStress_Warping (const FloatArrayF< 2 > &reducedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 6, 6 > | giveStiffnessMatrix_3d (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 3, 3 > | giveStiffnessMatrix_PlaneStress (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 4, 4 > | giveStiffnessMatrix_PlaneStrain (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 1, 1 > | giveStiffnessMatrix_1d (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 3 > | giveGeneralizedStress_Beam2d (const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 6 > | giveGeneralizedStress_Beam3d (const FloatArrayF< 6 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 5 > | giveGeneralizedStress_Plate (const FloatArrayF< 5 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 8 > | giveGeneralizedStress_Shell (const FloatArrayF< 8 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 9 > | giveGeneralizedStress_ShellRot (const FloatArrayF< 9 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 4 > | giveGeneralizedStress_MembraneRot (const FloatArrayF< 4 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 3 > | giveGeneralizedStress_PlateSubSoil (const FloatArrayF< 3 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 6 > | giveGeneralizedStress_3dBeamSubSoil (const FloatArrayF< 6 > &generalizedStrain, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArray | giveFirstPKStresses (const FloatArray &reducedF, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatArrayF< 9 > | giveFirstPKStress_3d (const FloatArrayF< 9 > &reducedF, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 5 > | giveFirstPKStress_PlaneStrain (const FloatArrayF< 5 > &reducedF, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 4 > | giveFirstPKStress_PlaneStress (const FloatArrayF< 4 > &reducedF, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArrayF< 1 > | giveFirstPKStress_1d (const FloatArrayF< 1 > &reducedF, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual void | giveCauchyStresses (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)=0 |
| virtual void | giveEshelbyStresses (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedvF, TimeStep *tStep) |
| virtual void | giveStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual void | giveCharMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual FloatMatrixF< 3, 3 > | give2dBeamStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 6, 6 > | give3dBeamStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 6, 6 > | give3dBeamSubSoilStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 5, 5 > | give2dPlateStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 8, 8 > | give3dShellStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 9, 9 > | give3dShellRotStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 6, 6 > | give3dDegeneratedShellStiffMtrx (MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep) const |
| virtual FloatMatrixF< 4, 4 > | giveMembraneRotStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 3, 3 > | give2dPlateSubSoilStiffMtrx (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual void | giveCharMaterialStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0 |
| virtual FloatMatrixF< 9, 9 > | giveStiffnessMatrix_dPdF_3d (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 5, 5 > | giveStiffnessMatrix_dPdF_PlaneStrain (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 4, 4 > | giveStiffnessMatrix_dPdF_PlaneStress (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatMatrixF< 1, 1 > | giveStiffnessMatrix_dPdF_1d (MatResponseMode mode, GaussPoint *gp, TimeStep *tStep) const =0 |
| virtual FloatArray * | imposeStressConstrainsOnGradient (GaussPoint *gp, FloatArray *gradientStressVector3d) |
| virtual FloatArray * | imposeStrainConstrainsOnGradient (GaussPoint *gp, FloatArray *gradientStressVector3d) |
| int | testCrossSectionExtension (CrossSectExtension ext) override |
| Material * | giveMaterial (IntegrationPoint *ip) const override |
| virtual Interface * | giveMaterialInterface (InterfaceType t, IntegrationPoint *ip) |
| virtual void | createMaterialStatus (GaussPoint &iGP)=0 |
| int | checkConsistency () override=0 |
| bool | isCharacteristicMtrxSymmetric (MatResponseMode mode) const override=0 |
| Public Member Functions inherited from oofem::CrossSection | |
| CrossSection (int n, Domain *d) | |
| virtual | ~CrossSection () |
| Destructor. | |
| int | giveSetNumber () const |
| virtual bool | hasProperty (CrossSectionProperty a) |
| virtual double | give (CrossSectionProperty a, GaussPoint *gp) const |
| virtual double | give (CrossSectionProperty a, const FloatArray &coords, Element *elem, bool local=true) const |
| virtual double | give (int aProperty, GaussPoint *gp) const |
| void | printYourself () override |
| Prints receiver state on stdout. Useful for debugging. | |
| virtual int | setupIntegrationPoints (IntegrationRule &irule, int npoints, Element *element) |
| virtual int | setupIntegrationPoints (IntegrationRule &irule, int npointsXY, int npointsZ, Element *element) |
| virtual int | giveIPValue (FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep) |
| virtual int | packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)=0 |
| virtual int | unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)=0 |
| virtual int | estimatePackSize (DataStream &buff, GaussPoint *ip)=0 |
| virtual double | predictRelativeComputationalCost (GaussPoint *ip) |
| virtual double | giveRelativeSelfComputationalCost () |
| virtual double | predictRelativeRedistributionCost (GaussPoint *gp) |
| void | initializeFrom (InputRecord &ir) override |
| void | giveInputRecord (DynamicInputRecord &input) override |
| virtual void | saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| virtual void | restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp) |
| 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 * | giveClassName () const =0 |
| 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). | |
Additional Inherited Members | |
| Protected Attributes inherited from oofem::CrossSection | |
| Dictionary | propertyDictionary |
| int | setNumber |
| Protected Attributes inherited from oofem::FEMComponent | |
| int | number |
| Component number. | |
| Domain * | domain |
| Link to domain object, useful for communicating with other FEM components. | |
Abstract base class for all structural cross section models. It declares commons services provided by all structural cross section models. The implementation of this services is left on derived classes, which will implement cross section model dependent part. However, some general services are implemented here. For information, how to introduce integration points in cross section volume for macro integration point, see CrossSection reference manual.
At structural level of cross section or constitutive models are introduced several stress/strain modes. Full and reduced formats of stress/strain vectors are also introduced for convenience. The full format includes all components, even if they are zero due to stress/strain mode nature, but in the reduced format, only generally nonzero components are stored. (full format must used only if absolutely necessary, to avoid wasting of space. It is used by output routines to print results in general form). Methods for converting vectors between full and reduced format are provided. General full strain vector has one of the following forms:
Definition at line 75 of file structuralcrosssection.h.
|
inline |
Constructor. Creates cross section with given number, belonging to given domain.
| n | Cross section number. |
| d | Domain to which new cross section will belong. |
Definition at line 83 of file structuralcrosssection.h.
References oofem::CrossSection::CrossSection().
Referenced by oofem::FiberedCrossSection::FiberedCrossSection(), oofem::LayeredCrossSection::LayeredCrossSection(), and oofem::SimpleCrossSection::SimpleCrossSection().
|
overridepure virtual |
Allows programmer to test some internal data, before computation begins. For example, one may use this function, to ensure that element has material with required capabilities is assigned to element. This must be done after all mesh components are instanciated.
Reimplemented from oofem::FEMComponent.
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
|
pure virtual |
|
pure virtual |
Computes the stiffness matrix for 2d beams.
| answer | The requested matrix. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step. |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::LIBeam2dNL::computeConstitutiveMatrixAt().
|
pure virtual |
Method for computing 2d plate stiffness matrix.
| answer | Stiffness matrix. |
| mode | 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). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::CCTPlate::computeConstitutiveMatrixAt(), oofem::DKTPlate::computeConstitutiveMatrixAt(), oofem::QDKTPlate::computeConstitutiveMatrixAt(), and oofem::Quad1Mindlin::computeConstitutiveMatrixAt().
|
pure virtual |
Method for computing subsoil stiffness matrix for plates.
| answer | Stiffness matrix. |
| mode | 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). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
|
pure virtual |
Computes the stiffness matrix for 2d beams.
| answer | The requested matrix. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step. |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
References give3dBeamSubSoilStiffMtrx().
Referenced by oofem::LIBeam3d2::computeConstitutiveMatrixAt(), oofem::LIBeam3dNL2::computeConstitutiveMatrixAt(), and oofem::LIBeam3dNL::computeConstitutiveMatrixAt().
|
virtual |
Method for computing subsoil stiffness matrix for 2d beams.
| answer | Stiffness matrix. |
| mode | 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). |
Definition at line 162 of file structuralcrosssection.C.
References oofem::StructuralMaterial::give3dBeamSubSoilStiffMtrx(), and giveMaterial().
Referenced by give3dBeamStiffMtrx().
|
inlinevirtual |
Method for computing 3d shell stiffness matrix on degenerated shell elements.
| answer | Stiffness matrix. |
| mode | 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 in oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Definition at line 274 of file structuralcrosssection.h.
References OOFEM_ERROR.
|
pure virtual |
Method for computing 3d shell (with normal rotation) stiffness matrix.
| answer | Stiffness matrix. |
| mode | 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). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
|
pure virtual |
Method for computing 3d shell stiffness matrix.
| answer | Stiffness matrix. |
| mode | 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). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::Quad1MindlinShell3D::computeConstitutiveMatrixAt().
|
pure virtual |
Computes the Cauchy stress vector for a given increment of deformation gradient and given integration point. 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. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains the Cauchy stress. |
| gp | Integration point. |
| reducedFIncrement | Increment of the deformation gradient vector in reduced form. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::NLStructuralElement::computeCauchyStressVector().
|
pure virtual |
Computes the stiffness matrix of receiver in given integration point, respecting its history. The algorithm should use temporary or equilibrium history variables stored in integration point status to compute and return required result. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains result. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::PhaseFieldElement::computeStiffnessMatrix_uu().
|
pure virtual |
Computes the large strain stiffness matrix of receiver in given integration point, respecting its history. The algorithm should use temporary or equilibrium history variables stored in integration point status to compute and return required result. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains result. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
|
inlinevirtual |
Computes the Eshelby stress vector. Does not update history variables, this is a postprocesing computation.
| answer | Contains the Eshelby stress. |
| gp | Integration point. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Reimplemented in oofem::SimpleCrossSection.
Definition at line 187 of file structuralcrosssection.h.
Referenced by oofem::MaterialForceEvaluator::computeMaterialForce().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by giveFirstPKStresses().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by giveFirstPKStresses().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by giveFirstPKStresses().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by giveFirstPKStresses().
|
virtual |
Computes the First Piola-Kirchoff stress vector for a given deformation gradient and integration point. 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. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains the First Piola-Kirchoff stresses. |
| gp | Integration point. |
| reducedFIncrement | Increment of the deformation gradient vector in reduced form. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Definition at line 83 of file structuralcrosssection.C.
References giveFirstPKStress_1d(), giveFirstPKStress_3d(), giveFirstPKStress_PlaneStrain(), giveFirstPKStress_PlaneStress(), oofem::GaussPoint::giveMaterialMode(), and OOFEM_ERROR.
Referenced by oofem::NLStructuralElement::computeFirstPKStressVector(), giveGeneralizedStress_PlateSubSoil(), and oofem::SolidShell::giveInternalForcesVector().
|
virtual |
Definition at line 155 of file structuralcrosssection.C.
References giveMaterial(), and oofem::StructuralMaterial::giveRealStressVector_3dBeamSubSoil().
Referenced by giveGeneralizedStress_PlateSubSoil().
|
pure virtual |
Computes the generalized stress vector for given strain and integration point.
| answer | Contains result. |
| gp | Integration point. |
| generalizedStrain | Strain vector in reduced generalized form. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::LIBeam2dNL::computeStressVector(), and giveRealStresses().
|
pure virtual |
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::CCTPlate::computeStressVector(), oofem::DKTPlate::computeStressVector(), oofem::QDKTPlate::computeStressVector(), oofem::Quad1Mindlin::computeStressVector(), oofem::DKTPlate::computeVertexBendingMoments(), and giveRealStresses().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
References giveFirstPKStresses(), and giveGeneralizedStress_3dBeamSubSoil().
|
pure virtual |
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by giveRealStresses().
|
inlineoverridevirtual |
Implements oofem::CrossSection.
Definition at line 356 of file structuralcrosssection.h.
References OOFEM_ERROR.
Referenced by give3dBeamSubSoilStiffMtrx(), oofem::MITC4Shell::giveCharacteristicTensor(), giveGeneralizedStress_3dBeamSubSoil(), oofem::QPlaneStress2dSlip::giveHomogenizedFields(), oofem::QTrPlaneStress2dSlip::giveHomogenizedFields(), giveRealStresses(), oofem::QPlaneStress2dSlip::giveSensitivities(), and oofem::QTrPlaneStress2dSlip::giveSensitivities().
|
inlinevirtual |
Reimplemented in oofem::SimpleCrossSection.
Definition at line 357 of file structuralcrosssection.h.
Referenced by oofem::GradientDamageElement::computeInternalForces_dB(), oofem::GradientDamageElement::computeInternalForces_dN(), oofem::GradientDamageElement::computeStiffnessMatrix_dd(), oofem::GradientDamageElement::computeStiffnessMatrix_du(), oofem::GradDpElement::computeStiffnessMatrix_kk(), oofem::GradDpElement::computeStiffnessMatrix_ku(), oofem::BaseMixedPressureElement::computeStiffnessMatrix_pp(), oofem::GradientDamageElement::computeStiffnessMatrix_ud(), oofem::GradDpElement::computeStiffnessMatrix_uk(), oofem::BaseMixedPressureElement::computeStiffnessMatrix_uu(), oofem::GradDpElement::computeStiffnessMatrix_uu(), oofem::GradientDamageElement::computeStiffnessMatrix_uu(), oofem::BaseMixedPressureElement::computeStressVector(), oofem::GradientDamageElement::computeStressVector_and_localDamageDrivingVariable(), oofem::GradDpElement::computeStressVectorAndLocalCumulatedStrain(), oofem::GradientDamageElement::giveInternalForcesVector_d(), oofem::BaseMixedPressureElement::giveInternalForcesVector_p(), and oofem::BaseMixedPressureElement::giveInternalForcesVector_u().
|
pure virtual |
Method for computing membrane stiffness matrix with added drilling stiffness.
| answer | Stiffness matrix. |
| mode | 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). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
|
pure virtual |
|
pure virtual |
|
inlinevirtual |
Reimplemented in oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Definition at line 102 of file structuralcrosssection.h.
References OOFEM_ERROR.
Referenced by oofem::MITC4Shell::computeStressVector().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::PlaneStrainElement::computeStressVector(), and giveRealStresses().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::PlaneStressElement::computeStressVector(), and giveRealStresses().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by giveRealStresses().
| FloatArray oofem::StructuralCrossSection::giveRealStresses | ( | const FloatArray & | reducedStrain, |
| GaussPoint * | gp, | ||
| TimeStep * | tStep ) const |
Computes the real stress vector for given strain and integration point. 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. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains result. |
| gp | Integration point. |
| reducedStrain | Strain vector in reduced form. |
| tStep | Current time step (most models are able to respond only when tStep is current time step). |
Definition at line 44 of file structuralcrosssection.C.
References giveGeneralizedStress_Beam2d(), giveGeneralizedStress_Beam3d(), giveGeneralizedStress_Plate(), giveGeneralizedStress_Shell(), giveGeneralizedStress_ShellRot(), giveMaterial(), oofem::GaussPoint::giveMaterialMode(), giveRealStress_1d(), giveRealStress_3d(), giveRealStress_PlaneStrain(), giveRealStress_PlaneStress(), giveRealStress_Warping(), and OOFEM_ERROR.
Referenced by oofem::SolidShell::giveInternalForcesVector().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::QTruss1d::computeConstitutiveMatrixAt(), oofem::Truss1d::computeConstitutiveMatrixAt(), oofem::Truss2d::computeConstitutiveMatrixAt(), and oofem::Truss3d::computeConstitutiveMatrixAt().
|
pure virtual |
Method for computing the stiffness matrix.
| answer | Stiffness matrix. |
| mode | 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). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::AxisymElement::computeConstitutiveMatrixAt(), oofem::Structural3DElement::computeConstitutiveMatrixAt(), and oofem::tet21ghostsolid::computeConstitutiveMatrixAt().
|
pure virtual |
Computes the material stiffness matrix dCde of receiver in a given integration point, respecting its history. The algorithm should use temporary or equilibrium history variables stored in integration point status to compute and return required result. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains result. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::NLStructuralElement::computeStiffnessMatrix(), and oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
References imposeStrainConstrainsOnGradient(), and imposeStressConstrainsOnGradient().
Referenced by oofem::QTruss1d::computeConstitutiveMatrix_dPdF_At(), oofem::Truss1d::computeConstitutiveMatrix_dPdF_At(), oofem::Truss2d::computeConstitutiveMatrix_dPdF_At(), and oofem::Truss3d::computeConstitutiveMatrix_dPdF_At().
|
pure virtual |
Computes the material stiffness matrix dPdF of receiver in a given integration point, respecting its history. The algorithm should use temporary or equilibrium history variables stored in integration point status to compute and return required result. Elements should always pass their requests to their cross section model, which performs necessary integration over its volume and invokes necessary material services for corresponding material model defined for given integration point.
| answer | Contains result. |
| mode | Material response mode. |
| gp | Integration point. |
| tStep | Time step (most models are able to respond only when tStep is current time step). |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::AxisymElement::computeConstitutiveMatrix_dPdF_At(), oofem::Structural3DElement::computeConstitutiveMatrix_dPdF_At(), and oofem::tet21ghostsolid::computeConstitutiveMatrix_dPdF_At().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::PlaneStrainElement::computeConstitutiveMatrix_dPdF_At().
|
pure virtual |
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::PlaneStressElement::computeConstitutiveMatrix_dPdF_At().
|
pure virtual |
|
pure virtual |
|
virtual |
Returns modified gradient of strain vector, which is used to compute plastic strain increment. Imposes zeros on places, where zero strain occurs or energetically connected stress is prescribed to be zero.
| gp | Integration point. |
| gradientStressVector3d | General 3d stress gradient. |
Reimplemented in oofem::FiberedCrossSection, and oofem::LayeredCrossSection.
Definition at line 171 of file structuralcrosssection.C.
References oofem::FloatArray::at(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatArray::giveSize(), and OOFEM_ERROR.
Referenced by oofem::PerfectlyPlasticMaterial::computePlasticStiffnessAt(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), giveStiffnessMatrix_dPdF_1d(), oofem::FiberedCrossSection::imposeStrainConstrainsOnGradient(), and oofem::LayeredCrossSection::imposeStrainConstrainsOnGradient().
|
virtual |
Returns modified gradient of stress vector, which is used to bring stresses back to yield surface. Method imposes zeros on places, where zero stress occurs. if energetically connected strain is zero, we do not impose zero there, because stress exist and must be taken into account when computing yield function. In such case a problem is assumed to be full 3d with some explicit strain equal to 0. On the other hand, if some stress is imposed to be zero, we understand such case as subspace of 3d case (like a classical plane stress problem, with no tracing of e_z, sigma_z)
| gp | Integration point. |
| gradientStressVector3d | General 3d stress gradient. |
Reimplemented in oofem::FiberedCrossSection, and oofem::LayeredCrossSection.
Definition at line 103 of file structuralcrosssection.C.
References oofem::FloatArray::at(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatArray::giveSize(), and OOFEM_ERROR.
Referenced by oofem::PerfectlyPlasticMaterial::computePlasticStiffnessAt(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), giveStiffnessMatrix_dPdF_1d(), oofem::PerfectlyPlasticMaterial::GiveStressCorrectionBackToYieldSurface(), oofem::FiberedCrossSection::imposeStressConstrainsOnGradient(), and oofem::LayeredCrossSection::imposeStressConstrainsOnGradient().
|
overridepure virtual |
Check for symmetry of stiffness matrix. Default implementation returns true. It can be moved to base Cross section class in the future.
| rMode | Response mode of material. |
Reimplemented from oofem::CrossSection.
Implemented in oofem::FiberedCrossSection, oofem::LayeredCrossSection, and oofem::SimpleCrossSection.
Referenced by oofem::NLStructuralElement::computeStiffnessMatrix(), oofem::StructuralElementEvaluator::computeStiffnessMatrix(), oofem::Truss3dnl2::computeStiffnessMatrix(), oofem::Truss3dnl::computeStiffnessMatrix(), oofem::PhaseFieldElement::computeStiffnessMatrix_uu(), and oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells().
|
inlineoverridevirtual |
Returns nonzero, if receiver implements required extension.
| ext | Required extension. |
Reimplemented from oofem::CrossSection.
Definition at line 353 of file structuralcrosssection.h.
References oofem::CS_StructuralCapability.