Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Section Class Referenceabstract

Virtual class for creating a section object. More...

#include <Section.hpp>

Inheritance diagram for Section:
Collaboration diagram for Section:

Public Member Functions

 Section (std::string name)
 Creates a Section to be specified at a Gauss-point in an Element. More...
 
virtual ~Section ()=0
 Destroys this Section object. More...
 
virtual std::unique_ptr< SectionCopySection ()=0
 Clone the selected Section. More...
 
virtual Eigen::VectorXd GetStrain ()=0
 Returns the resultant (generalised) strain vector over the section. More...
 
virtual Eigen::VectorXd GetStress ()=0
 Returns the resultant (generalised) stress vector over the section. More...
 
virtual Eigen::MatrixXd GetDensity ()=0
 Access the section density matrix. More...
 
virtual Eigen::MatrixXd GetTangentStiffness ()=0
 Returns the section stiffness matrix. More...
 
virtual Eigen::MatrixXd GetInitialTangentStiffness ()=0
 Returns the section initial stiffness matrix. More...
 
virtual Eigen::VectorXd GetStrainAt (double x3=0, double x2=0)=0
 Returns the section strain at given position. More...
 
virtual Eigen::VectorXd GetStressAt (double x3=0, double x2=0)=0
 Returns the section stress at given position. More...
 
virtual void CommitState ()=0
 Perform converged section state update. More...
 
virtual void ReverseState ()=0
 Reverse the section states to previous converged state. More...
 
virtual void InitialState ()=0
 Brings the section states to its initial state. More...
 
virtual void UpdateState (Eigen::VectorXd strain, unsigned int cond=0)=0
 Update the section state for this iteration. More...
 
std::string GetName ()
 Gets Section name. More...
 

Protected Member Functions

Eigen::MatrixXd GetLineRotationMatrix (double theta)
 Gets Local Axis Rotation for Line Section according to provided angle. More...
 
Eigen::MatrixXd GetAreaRotationMatrix (double Theta)
 Gets Local Axis Rotation for Area Section according to provided angle. More...
 
Eigen::MatrixXd GetLineTranslationMatrix (double h, double b, double zc, double yc, unsigned int ip)
 Gets centroid translation axis for Line Section according to insertion point. More...
 
void InsertionPointCoordinates (double &x3, double &x2, double h, double b, double zc, double yc, unsigned int ip)
 Gets the coordinate according to insertion point. More...
 
Eigen::MatrixXd ComputeLineLocalAxes (double h, double b, double zcm, double ycm, double angle, unsigned int ip)
 Transforms generalised strain/stresses from Element to Section local coordinate. More...
 

Private Attributes

std::string Name
 Section Name. More...
 

Detailed Description

Virtual class for creating a section object.

See also
Mesh.hpp

Constructor & Destructor Documentation

◆ Section()

Section::Section ( std::string  name)

Creates a Section to be specified at a Gauss-point in an Element.

Parameters
nameName of the derived class.
See also
Section::Name.

◆ ~Section()

virtual Section::~Section ( )
pure virtual

Destroys this Section object.

Member Function Documentation

◆ CommitState()

virtual void Section::CommitState ( )
pure virtual

◆ ComputeLineLocalAxes()

Eigen::MatrixXd Section::ComputeLineLocalAxes ( double  h,
double  b,
double  zcm,
double  ycm,
double  angle,
unsigned int  ip 
)
protected

Transforms generalised strain/stresses from Element to Section local coordinate.

Parameters
hThe total height of the section.
bThe total width of the section.
zcmPosition of the section axis along x3-axis.
ycmPosition of the section axis along x2-axis.
ipThe insertion point where the section axes are located.
Returns
The section translation matrix according to ip.
Note
More details can be found at Insertion Point.

◆ CopySection()

◆ GetAreaRotationMatrix()

Eigen::MatrixXd Section::GetAreaRotationMatrix ( double  Theta)
protected

Gets Local Axis Rotation for Area Section according to provided angle.

Parameters
thetaThe angle of rotation about x3-axis.
Returns
The section rotation matrix in theta degrees.
Note
More details can be found at Local Axes.

◆ GetDensity()

virtual Eigen::MatrixXd Section::GetDensity ( )
pure virtual

◆ GetInitialTangentStiffness()

virtual Eigen::MatrixXd Section::GetInitialTangentStiffness ( )
pure virtual

Returns the section initial stiffness matrix.

Returns
Matrix with the initial section tangent stiffness matrix.
Note
The initial tangent stiffness matrix is computed when the generalized strain vector is zero.

Implemented in Fib3DLineSection, Lin2DAngle, Lin2DChannel, Lin2DRectangularTube, Lin2DTee, Lin2DWideFlange, Lin3DAngle, Lin3DChannel, Lin3DRectangularTube, Lin3DTee, Lin3DWideFlange, Lin2DRectangular, Lin3DCircularTube, Lin3DRectangular, Lin2DCircular, Lin2DCircularTube, Lin2DUserDefined, Lin3DCircular, Lin3DUserDefined, and Lin3DThinArea.

◆ GetLineRotationMatrix()

Eigen::MatrixXd Section::GetLineRotationMatrix ( double  theta)
protected

Gets Local Axis Rotation for Line Section according to provided angle.

Parameters
thetaThe angle of rotation about x1-axis.
Returns
The section rotation matrix in theta degrees.
Note
More details can be found at Local Axes.

◆ GetLineTranslationMatrix()

Eigen::MatrixXd Section::GetLineTranslationMatrix ( double  h,
double  b,
double  zc,
double  yc,
unsigned int  ip 
)
protected

Gets centroid translation axis for Line Section according to insertion point.

Parameters
hThe total height of the section.
bThe total width of the section.
zcPosition of the section axis along x3-axis.
ycPosition of the section axis along x2-axis.
ipThe insertion point where the section axes are located.
Returns
The section translation matrix according to ip.
Note
More details can be found at Insertion Point.

◆ GetName()

std::string Section::GetName ( )

Gets Section name.

See also
Section::Name.

◆ GetStrain()

virtual Eigen::VectorXd Section::GetStrain ( )
pure virtual

◆ GetStrainAt()

virtual Eigen::VectorXd Section::GetStrainAt ( double  x3 = 0,
double  x2 = 0 
)
pure virtual

Returns the section strain at given position.

Parameters
x3Position on the x3-axis where the strain is evaluated.
x2Position on the x2-axis where the strain is evaluated.
Returns
Vector with the strain components at coordinate (x3,x2).

Implemented in Fib3DLineSection, Lin2DAngle, Lin2DChannel, Lin2DRectangularTube, Lin2DTee, Lin2DWideFlange, Lin3DAngle, Lin3DChannel, Lin3DRectangularTube, Lin3DTee, Lin3DWideFlange, Lin2DRectangular, Lin3DCircularTube, Lin3DRectangular, Lin2DCircular, Lin2DCircularTube, Lin2DUserDefined, Lin3DCircular, Lin3DUserDefined, and Lin3DThinArea.

◆ GetStress()

virtual Eigen::VectorXd Section::GetStress ( )
pure virtual

◆ GetStressAt()

virtual Eigen::VectorXd Section::GetStressAt ( double  x3 = 0,
double  x2 = 0 
)
pure virtual

Returns the section stress at given position.

Parameters
x3Position on the x3-axis where the stress is evaluated.
x2Position on the x2-axis where the stress is evaluated.
Returns
Vector with the stress components at coordinate (x3,x2).

Implemented in Fib3DLineSection, Lin2DAngle, Lin2DChannel, Lin2DRectangularTube, Lin2DTee, Lin2DWideFlange, Lin3DAngle, Lin3DChannel, Lin3DRectangularTube, Lin3DTee, Lin3DWideFlange, Lin2DRectangular, Lin3DCircularTube, Lin3DRectangular, Lin2DCircular, Lin2DCircularTube, Lin2DUserDefined, Lin3DCircular, Lin3DUserDefined, and Lin3DThinArea.

◆ GetTangentStiffness()

virtual Eigen::MatrixXd Section::GetTangentStiffness ( )
pure virtual

◆ InitialState()

virtual void Section::InitialState ( )
pure virtual

◆ InsertionPointCoordinates()

void Section::InsertionPointCoordinates ( double &  x3,
double &  x2,
double  h,
double  b,
double  zc,
double  yc,
unsigned int  ip 
)
protected

Gets the coordinate according to insertion point.

Parameters
x3The coordinate along the 3-axis.
x2The coordinate along the 2-axis.
hThe total height of the section.
bThe total width of the section.
zcPosition of the section axis along x3-axis.
ycPosition of the section axis along x2-axis.
ipThe insertion point where the section axes are located.
Returns
The section translation matrix according to ip.
Note
More details can be found at Insertion Point.

◆ ReverseState()

virtual void Section::ReverseState ( )
pure virtual

Reverse the section states to previous converged state.

Note
This function returns the material states to previous converged states.

Implemented in Fib3DLineSection, Lin2DAngle, Lin2DChannel, Lin2DRectangularTube, Lin2DTee, Lin2DWideFlange, Lin3DAngle, Lin3DChannel, Lin3DRectangularTube, Lin3DTee, Lin3DWideFlange, Lin2DRectangular, Lin3DCircularTube, Lin3DRectangular, Lin2DCircular, Lin2DCircularTube, Lin2DUserDefined, Lin3DCircular, Lin3DUserDefined, and Lin3DThinArea.

◆ UpdateState()

virtual void Section::UpdateState ( Eigen::VectorXd  strain,
unsigned int  cond = 0 
)
pure virtual

Update the section state for this iteration.

Parameters
strainVector with the strain components at this Gauss-point.
condIf the the elastic/plastic material components will be updated.
Note
This function computes the strain and tanget stiffness matrix once the trial strain converged.

Implemented in Fib3DLineSection, Lin2DAngle, Lin2DChannel, Lin2DRectangularTube, Lin2DTee, Lin2DWideFlange, Lin3DAngle, Lin3DChannel, Lin3DRectangularTube, Lin3DTee, Lin3DWideFlange, Lin2DRectangular, Lin3DCircularTube, Lin3DRectangular, Lin2DCircular, Lin2DCircularTube, Lin2DUserDefined, Lin3DCircular, Lin3DUserDefined, and Lin3DThinArea.

Member Data Documentation

◆ Name

std::string Section::Name
private

Section Name.