OOFEM  2.4
OOFEM.org - Object Oriented Finite Element Solver
oofem::StructuralCrossSection Class Referenceabstract

Abstract base class for all structural cross section models. More...

#include <structuralcrosssection.h>

+ Inheritance diagram for oofem::StructuralCrossSection:
+ Collaboration diagram for oofem::StructuralCrossSection:

Public Member Functions

 StructuralCrossSection (int n, Domain *d)
 Constructor. More...
 
virtual ~StructuralCrossSection ()
 Destructor. More...
 
virtual void giveFirstPKStresses (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)=0
 Computes the First Piola-Kirchoff stress vector for a given deformation gradient and integration point. More...
 
virtual void giveCauchyStresses (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedFIncrement, TimeStep *tStep)=0
 Computes the Cauchy stress vector for a given increment of deformation gradient and given integration point. More...
 
virtual void giveEshelbyStresses (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedvF, TimeStep *tStep)
 Computes the Eshelby stress vector. More...
 
virtual void giveStiffnessMatrix_dPdF (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Computes the material stiffness matrix dPdF of receiver in a given integration point, respecting its history. More...
 
virtual void giveStiffnessMatrix_dCde (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Computes the material stiffness matrix dCde of receiver in a given integration point, respecting its history. More...
 
virtual void giveCharMaterialStiffnessMatrix (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Computes the stiffness matrix of receiver in given integration point, respecting its history. More...
 
virtual void give2dBeamStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Computes the stiffness matrix for 2d beams. More...
 
virtual void give3dBeamStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Computes the stiffness matrix for 2d beams. More...
 
virtual void give3dBeamSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)
 Method for computing subsoil stiffness matrix for 2d beams. More...
 
virtual void give2dPlateStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Method for computing 2d plate stiffness matrix. More...
 
virtual void give3dShellStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Method for computing 3d shell stiffness matrix. More...
 
virtual void give3dDegeneratedShellStiffMtrx (FloatMatrix &answer, MatResponseMode rMode, GaussPoint *gp, TimeStep *tStep)
 Method for computing 3d shell stiffness matrix on degenerated shell elements. More...
 
virtual void giveMembraneRotStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Method for computing membrane stiffness matrix with added drilling stiffness. More...
 
virtual void give2dPlateSubSoilStiffMtrx (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Method for computing subsoil stiffness matrix for plates. More...
 
virtual FloatArrayimposeStressConstrainsOnGradient (GaussPoint *gp, FloatArray *gradientStressVector3d)
 Returns modified gradient of stress vector, which is used to bring stresses back to yield surface. More...
 
virtual FloatArrayimposeStrainConstrainsOnGradient (GaussPoint *gp, FloatArray *gradientStressVector3d)
 Returns modified gradient of strain vector, which is used to compute plastic strain increment. More...
 
virtual int testCrossSectionExtension (CrossSectExtension ext)
 Returns nonzero, if receiver implements required extension. More...
 
virtual MaterialgiveMaterial (IntegrationPoint *ip)
 Returns the material associated with the GP. More...
 
virtual void createMaterialStatus (GaussPoint &iGP)=0
 
virtual int checkConsistency ()=0
 Allows programmer to test some internal data, before computation begins. More...
 
virtual InterfacegiveMaterialInterface (InterfaceType t, IntegrationPoint *ip)
 
virtual bool isCharacteristicMtrxSymmetric (MatResponseMode mode)=0
 Check for symmetry of stiffness matrix. More...
 
void giveRealStresses (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
 Computes the real stress vector for given strain and integration point. More...
 
virtual void giveRealStress_3d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)=0
 
virtual void giveRealStress_3dDegeneratedShell (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)
 
virtual void giveRealStress_PlaneStrain (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)=0
 
virtual void giveRealStress_PlaneStress (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)=0
 
virtual void giveRealStress_1d (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)=0
 
virtual void giveRealStress_Warping (FloatArray &answer, GaussPoint *gp, const FloatArray &reducedStrain, TimeStep *tStep)=0
 
virtual void giveStiffnessMatrix_3d (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 Method for computing the stiffness matrix. More...
 
virtual void giveStiffnessMatrix_PlaneStress (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 
virtual void giveStiffnessMatrix_PlaneStrain (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 
virtual void giveStiffnessMatrix_1d (FloatMatrix &answer, MatResponseMode mode, GaussPoint *gp, TimeStep *tStep)=0
 
virtual void giveGeneralizedStress_Beam2d (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
 Computes the generalized stress vector for given strain and integration point. More...
 
virtual void giveGeneralizedStress_Beam3d (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
 
virtual void giveGeneralizedStress_Plate (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
 
virtual void giveGeneralizedStress_Shell (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
 
virtual void giveGeneralizedStress_MembraneRot (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
 
virtual void giveGeneralizedStress_PlateSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)=0
 
virtual void giveGeneralizedStress_3dBeamSubSoil (FloatArray &answer, GaussPoint *gp, const FloatArray &generalizedStrain, TimeStep *tStep)
 
- Public Member Functions inherited from oofem::CrossSection
 CrossSection (int n, Domain *d)
 Constructor. More...
 
virtual ~CrossSection ()
 Destructor. More...
 
int giveSetNumber () const
 
virtual bool hasProperty (CrossSectionProperty a)
 Returns true if the dictionary contains the requested property. More...
 
virtual double give (CrossSectionProperty a, GaussPoint *gp)
 Returns the value of cross section property at given point. More...
 
virtual double give (CrossSectionProperty a, const FloatArray &coords, Element *elem, bool local=true)
 Returns the value of cross section property at given point (belonging to given element). More...
 
virtual double give (int aProperty, GaussPoint *gp)
 Returns the value of cross section property. More...
 
virtual void printYourself ()
 Prints receiver state on stdout. Useful for debugging. More...
 
virtual int setupIntegrationPoints (IntegrationRule &irule, int npoints, Element *element)
 Sets up integration rule for the given element. More...
 
virtual int setupIntegrationPoints (IntegrationRule &irule, int npointsXY, int npointsZ, Element *element)
 Sets up integration rule for the given element. More...
 
virtual int giveIPValue (FloatArray &answer, GaussPoint *ip, InternalStateType type, TimeStep *tStep)
 Returns the integration point corresponding value in Reduced form. More...
 
virtual int packUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)=0
 Pack all necessary data of integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int unpackAndUpdateUnknowns (DataStream &buff, TimeStep *tStep, GaussPoint *ip)=0
 Unpack and updates all necessary data of given integration point (according to element parallel_mode) into given communication buffer. More...
 
virtual int estimatePackSize (DataStream &buff, GaussPoint *ip)=0
 Estimates the necessary pack size to hold all packed data of receiver. More...
 
virtual double predictRelativeComputationalCost (GaussPoint *ip)
 Returns the weight representing relative computational cost of receiver The reference cross section is integral model in plane stress. More...
 
virtual double giveRelativeSelfComputationalCost ()
 Returns the weight representing relative computational cost of receiver The reference element is integral model in plane stress. More...
 
virtual double predictRelativeRedistributionCost (GaussPoint *gp)
 
virtual IRResultType initializeFrom (InputRecord *ir)
 Initializes receiver according to object description stored in input record. More...
 
virtual void giveInputRecord (DynamicInputRecord &input)
 Setups the input record string of receiver. More...
 
virtual contextIOResultType saveIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Stores integration point state to output stream. More...
 
virtual contextIOResultType restoreIPContext (DataStream &stream, ContextMode mode, GaussPoint *gp)
 Reads integration point state to output stream. More...
 
- Public Member Functions inherited from oofem::FEMComponent
 FEMComponent (int n, Domain *d)
 Regular constructor, creates component with given number and belonging to given domain. More...
 
virtual ~FEMComponent ()
 Virtual destructor. More...
 
virtual const char * giveClassName () const =0
 
virtual const char * giveInputRecordName () const =0
 
DomaingiveDomain () const
 
virtual void setDomain (Domain *d)
 Sets associated Domain. More...
 
int giveNumber () const
 
void setNumber (int num)
 Sets number of receiver. More...
 
virtual void updateLocalNumbering (EntityRenumberingFunctor &f)
 Local renumbering support. More...
 
virtual contextIOResultType saveContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Stores receiver state to output stream. More...
 
virtual contextIOResultType restoreContext (DataStream &stream, ContextMode mode, void *obj=NULL)
 Restores the receiver state previously written in stream. More...
 
virtual void printOutputAt (FILE *file, TimeStep *tStep)
 Prints output of receiver to stream, for given time step. More...
 
virtual InterfacegiveInterface (InterfaceType t)
 Interface requesting service. More...
 
std::string errorInfo (const char *func) const
 Returns string for prepending output (used by error reporting macros). More...
 

Additional Inherited Members

- Protected Attributes inherited from oofem::CrossSection
Dictionary propertyDictionary
 Dictionary for storing cross section parameters (like dimensions). More...
 
int setNumber
 
- Protected Attributes inherited from oofem::FEMComponent
int number
 Component number. More...
 
Domaindomain
 Link to domain object, useful for communicating with other FEM components. More...
 

Detailed Description

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:

  1. strainVector3d {eps_x,eps_y,eps_z,gamma_yz,gamma_zx,gamma_xy}
  2. For integrated cross section models (2d and 3d beams, plates and general shells) strainVectorShell {eps_x,eps_y,gamma_xy, kappa_x, kappa_y, kappa_xy, gamma_zx, gamma_zy}

Definition at line 75 of file structuralcrosssection.h.

Constructor & Destructor Documentation

oofem::StructuralCrossSection::StructuralCrossSection ( int  n,
Domain d 
)
inline

Constructor.

Creates cross section with given number, belonging to given domain.

Parameters
nCross section number.
dDomain to which new cross section will belong.

Definition at line 83 of file structuralcrosssection.h.

virtual oofem::StructuralCrossSection::~StructuralCrossSection ( )
inlinevirtual

Destructor.

Definition at line 85 of file structuralcrosssection.h.

Member Function Documentation

virtual int oofem::StructuralCrossSection::checkConsistency ( )
pure 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.

Returns
Nonzero if receiver is consistent.

Reimplemented from oofem::FEMComponent.

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

virtual void oofem::StructuralCrossSection::give2dBeamStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual

Computes the stiffness matrix for 2d beams.

Parameters
answerThe requested matrix.
modeMaterial response mode.
gpIntegration point.
tStepTime step.

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::LIBeam2dNL::computeConstitutiveMatrixAt(), oofem::LIBeam2d::computeConstitutiveMatrixAt(), and oofem::Beam2d::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::give2dPlateStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual

Method for computing 2d plate stiffness matrix.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::CCTPlate::computeConstitutiveMatrixAt(), oofem::QDKTPlate::computeConstitutiveMatrixAt(), oofem::DKTPlate::computeConstitutiveMatrixAt(), and oofem::Quad1Mindlin::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::give2dPlateSubSoilStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual

Method for computing subsoil stiffness matrix for plates.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::Quad1PlateSubSoil::computeConstitutiveMatrixAt(), and oofem::Tria1PlateSubSoil::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::give3dBeamStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual
void oofem::StructuralCrossSection::give3dBeamSubSoilStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
virtual

Method for computing subsoil stiffness matrix for 2d beams.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Definition at line 138 of file structuralcrosssection.C.

References oofem::StructuralMaterial::give3dBeamSubSoilStiffMtrx(), and giveMaterial().

virtual void oofem::StructuralCrossSection::give3dDegeneratedShellStiffMtrx ( FloatMatrix answer,
MatResponseMode  rMode,
GaussPoint gp,
TimeStep tStep 
)
inlinevirtual

Method for computing 3d shell stiffness matrix on degenerated shell elements.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::LayeredCrossSection, and oofem::SimpleCrossSection.

Definition at line 274 of file structuralcrosssection.h.

References OOFEM_ERROR.

Referenced by oofem::MITC4Shell::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::give3dShellStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual

Method for computing 3d shell stiffness matrix.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::RerShell::computeConstitutiveMatrixAt(), and oofem::Quad1MindlinShell3D::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::giveCauchyStresses ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedFIncrement,
TimeStep tStep 
)
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.

Parameters
answerContains the Cauchy stress.
gpIntegration point.
reducedFIncrementIncrement of the deformation gradient vector in reduced form.
Todo:
should this then be in a multiplicative way? /JB
Parameters
tStepCurrent time step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::NLStructuralElement::computeCauchyStressVector().

virtual void oofem::StructuralCrossSection::giveCharMaterialStiffnessMatrix ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
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.

Parameters
answerContains result.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::Tr_Warp::computeConstitutiveMatrixAt(), oofem::Lattice2d::computeConstitutiveMatrixAt(), oofem::PhaseFieldElement::computeStiffnessMatrix_uu(), and oofem::StructuralElement::~StructuralElement().

virtual void oofem::StructuralCrossSection::giveEshelbyStresses ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedvF,
TimeStep tStep 
)
inlinevirtual

Computes the Eshelby stress vector.

Does not update history variables, this is a postprocesing computation.

Parameters
answerContains the Eshelby stress.
gpIntegration point.
tStepCurrent time step (most models are able to respond only when tStep is current time step).

Reimplemented in oofem::SimpleCrossSection.

Definition at line 181 of file structuralcrosssection.h.

Referenced by oofem::MaterialForceEvaluator::computeMaterialForce().

virtual void oofem::StructuralCrossSection::giveFirstPKStresses ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedFIncrement,
TimeStep tStep 
)
pure 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.

Parameters
answerContains the First Piola-Kirchoff stresses.
gpIntegration point.
reducedFIncrementIncrement of the deformation gradient vector in reduced form.
Todo:
should this then be in a multiplicative way? /JB
Parameters
tStepCurrent time step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::NLStructuralElement::computeFirstPKStressVector(), and oofem::SolidShell::giveInternalForcesVector().

void oofem::StructuralCrossSection::giveGeneralizedStress_3dBeamSubSoil ( FloatArray answer,
GaussPoint gp,
const FloatArray generalizedStrain,
TimeStep tStep 
)
virtual
virtual void oofem::StructuralCrossSection::giveGeneralizedStress_Beam2d ( FloatArray answer,
GaussPoint gp,
const FloatArray generalizedStrain,
TimeStep tStep 
)
pure virtual

Computes the generalized stress vector for given strain and integration point.

Parameters
answerContains result.
gpIntegration point.
generalizedStrainStrain vector in reduced generalized form.
tStepCurrent time step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::LIBeam2dNL::computeStressVector(), oofem::LIBeam2d::computeStressVector(), oofem::Beam2d::computeStressVector(), and giveRealStresses().

virtual void oofem::StructuralCrossSection::giveGeneralizedStress_MembraneRot ( FloatArray answer,
GaussPoint gp,
const FloatArray generalizedStrain,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveGeneralizedStress_PlateSubSoil ( FloatArray answer,
GaussPoint gp,
const FloatArray generalizedStrain,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveGeneralizedStress_Shell ( FloatArray answer,
GaussPoint gp,
const FloatArray generalizedStrain,
TimeStep tStep 
)
pure virtual
virtual Material* oofem::StructuralCrossSection::giveMaterial ( IntegrationPoint ip)
inlinevirtual

Returns the material associated with the GP.

Default implementation uses gp->giveMaterial() for backwards compatibility, but it should be overloaded in each specialized cross-section.

Implements oofem::CrossSection.

Reimplemented in oofem::LayeredCrossSection, and oofem::SimpleCrossSection.

Definition at line 317 of file structuralcrosssection.h.

References OOFEM_ERROR.

Referenced by give3dBeamSubSoilStiffMtrx(), oofem::MITC4Shell::giveCharacteristicTensor(), giveGeneralizedStress_3dBeamSubSoil(), oofem::MITC4Shell::giveMidplaneIPValue(), giveRealStresses(), and oofem::LinearizedDilationForceAssembler::vectorFromElement().

virtual void oofem::StructuralCrossSection::giveMembraneRotStiffMtrx ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual

Method for computing membrane stiffness matrix with added drilling stiffness.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::TrPlaneStrRot::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::giveRealStress_1d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveRealStress_3d ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveRealStress_3dDegeneratedShell ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
inlinevirtual
virtual void oofem::StructuralCrossSection::giveRealStress_PlaneStrain ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveRealStress_PlaneStress ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveRealStress_Warping ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)
pure virtual
void oofem::StructuralCrossSection::giveRealStresses ( FloatArray answer,
GaussPoint gp,
const FloatArray reducedStrain,
TimeStep tStep 
)

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.

Parameters
answerContains result.
gpIntegration point.
reducedStrainStrain vector in reduced form.
tStepCurrent time step (most models are able to respond only when tStep is current time step).
Todo:
this part only works for simple cross section and will be removed soon when new interface elements are done /JB

Definition at line 44 of file structuralcrosssection.C.

References giveGeneralizedStress_Beam2d(), giveGeneralizedStress_Beam3d(), giveGeneralizedStress_Plate(), giveGeneralizedStress_Shell(), giveMaterial(), oofem::GaussPoint::giveMaterialMode(), giveRealStress_1d(), giveRealStress_3d(), giveRealStress_PlaneStrain(), giveRealStress_PlaneStress(), giveRealStress_Warping(), oofem::StructuralMaterial::giveRealStressVector(), oofem::StructuralMaterial::hasMaterialModeCapability(), and OOFEM_ERROR.

Referenced by oofem::StructuralElement::computeStrainVector(), oofem::Lattice2d::computeStressVector(), and oofem::SolidShell::giveInternalForcesVector().

virtual void oofem::StructuralCrossSection::giveStiffnessMatrix_3d ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual

Method for computing the stiffness matrix.

Parameters
answerStiffness matrix.
modeMaterial response mode.
gpIntegration point, which load history is used.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::tet21ghostsolid::computeConstitutiveMatrixAt(), oofem::Structural3DElement::computeConstitutiveMatrixAt(), and oofem::AxisymElement::computeConstitutiveMatrixAt().

virtual void oofem::StructuralCrossSection::giveStiffnessMatrix_dCde ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
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.

Parameters
answerContains result.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::NLStructuralElement::computeStiffnessMatrix(), and oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells().

virtual void oofem::StructuralCrossSection::giveStiffnessMatrix_dPdF ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
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.

Parameters
answerContains result.
modeMaterial response mode.
gpIntegration point.
tStepTime step (most models are able to respond only when tStep is current time step).

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::SolidShell::computeStiffnessMatrix(), oofem::NLStructuralElement::computeStiffnessMatrix(), and oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells().

virtual void oofem::StructuralCrossSection::giveStiffnessMatrix_PlaneStrain ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual
virtual void oofem::StructuralCrossSection::giveStiffnessMatrix_PlaneStress ( FloatMatrix answer,
MatResponseMode  mode,
GaussPoint gp,
TimeStep tStep 
)
pure virtual
FloatArray * oofem::StructuralCrossSection::imposeStrainConstrainsOnGradient ( GaussPoint gp,
FloatArray gradientStressVector3d 
)
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.

See also
imposeStressConstrainsOnGradient
Parameters
gpIntegration point.
gradientStressVector3dGeneral 3d stress gradient.

Reimplemented in oofem::LayeredCrossSection, and oofem::FiberedCrossSection.

Definition at line 148 of file structuralcrosssection.C.

References oofem::__MaterialModeToString(), oofem::FloatArray::at(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatArray::giveSize(), and OOFEM_ERROR.

Referenced by oofem::PerfectlyPlasticMaterial::computePlasticStiffnessAt(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), oofem::FiberedCrossSection::imposeStrainConstrainsOnGradient(), and oofem::LayeredCrossSection::imposeStrainConstrainsOnGradient().

FloatArray * oofem::StructuralCrossSection::imposeStressConstrainsOnGradient ( GaussPoint gp,
FloatArray gradientStressVector3d 
)
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)

Parameters
gpIntegration point.
gradientStressVector3dGeneral 3d stress gradient.

Reimplemented in oofem::LayeredCrossSection, and oofem::FiberedCrossSection.

Definition at line 79 of file structuralcrosssection.C.

References oofem::__MaterialModeToString(), oofem::FloatArray::at(), oofem::GaussPoint::giveMaterialMode(), oofem::FloatArray::giveSize(), and OOFEM_ERROR.

Referenced by oofem::PerfectlyPlasticMaterial::computePlasticStiffnessAt(), oofem::PerfectlyPlasticMaterial::giveRealStressVector(), oofem::PerfectlyPlasticMaterial::GiveStressCorrectionBackToYieldSurface(), oofem::FiberedCrossSection::imposeStressConstrainsOnGradient(), and oofem::LayeredCrossSection::imposeStressConstrainsOnGradient().

virtual bool oofem::StructuralCrossSection::isCharacteristicMtrxSymmetric ( MatResponseMode  rMode)
pure virtual

Check for symmetry of stiffness matrix.

Default implementation returns true. It can be moved to base Cross section class in the future.

Parameters
rModeResponse mode of material.
Returns
True if stiffness matrix of receiver is symmetric.

Reimplemented from oofem::CrossSection.

Implemented in oofem::LayeredCrossSection, oofem::FiberedCrossSection, and oofem::SimpleCrossSection.

Referenced by oofem::StructuralElementEvaluator::computeStiffnessMatrix(), oofem::NLStructuralElement::computeStiffnessMatrix(), oofem::PhaseFieldElement::computeStiffnessMatrix_uu(), and oofem::NLStructuralElement::computeStiffnessMatrix_withIRulesAsSubcells().

virtual int oofem::StructuralCrossSection::testCrossSectionExtension ( CrossSectExtension  ext)
inlinevirtual

Returns nonzero, if receiver implements required extension.

Parameters
extRequired extension.
Returns
Nonzero, if supported, zero otherwise.

Reimplemented from oofem::CrossSection.

Definition at line 315 of file structuralcrosssection.h.

References oofem::CS_StructuralCapability.


The documentation for this class was generated from the following files:

This page is part of the OOFEM documentation. Copyright (c) 2011 Borek Patzak
Project e-mail: info@oofem.org
Generated at Tue Jan 2 2018 20:07:41 for OOFEM by doxygen 1.8.11 written by Dimitri van Heesch, © 1997-2011