Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Lin3DAngle Class Reference

Class for creating a angle section for 3D analysis with linear elastic material for a frame elements. More...

#include <Lin3DAngle.hpp>

Inheritance diagram for Lin3DAngle:
Collaboration diagram for Lin3DAngle:

Public Member Functions

 Lin3DAngle (double h, double b, double tw, double tf, std::unique_ptr< Material > &material, double theta=0.0, unsigned int ip=10)
 Creates a Section to be specified at a Gauss-point in an Element. More...
 
 ~Lin3DAngle ()
 Destroys this Lin3DAngle object. More...
 
std::unique_ptr< SectionCopySection ()
 Clone the 'Lin3DAngle' section. More...
 
Eigen::VectorXd GetStrain ()
 Returns the resultant (generalised) strain vector over the section. More...
 
Eigen::VectorXd GetStress ()
 Returns the resultant (generalised) stress vector over the section. More...
 
Eigen::MatrixXd GetDensity ()
 Access the section density matrix. More...
 
Eigen::MatrixXd GetTangentStiffness ()
 Returns the section stiffness matrix. More...
 
Eigen::MatrixXd GetInitialTangentStiffness ()
 Returns the section initial stiffness matrix. More...
 
Eigen::VectorXd GetStrainAt (double x3=0, double x2=0)
 Returns the section strain at given position. More...
 
Eigen::VectorXd GetStressAt (double x3=0, double x2=0)
 Returns the section stress at given position. More...
 
void CommitState ()
 Perform converged section state update. More...
 
void ReverseState ()
 Reverse the section states to previous converged state. More...
 
void InitialState ()
 Brings the section states to its initial state. More...
 
void UpdateState (const Eigen::VectorXd strain, const unsigned int cond)
 Update the section state for this iteration. More...
 
- Public Member Functions inherited from Section
 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...
 
std::string GetName ()
 Gets Section name. More...
 

Protected Member Functions

double GetArea ()
 Computes the Angle area. More...
 
double GetShearArea2 ()
 Computes the Angle shear area. More...
 
double GetShearArea3 ()
 Computes the Angle shear area. More...
 
double GetInertiaAxis1 ()
 
double GetInertiaAxis2 ()
 Computes the Angle flexural inertia. More...
 
double GetInertiaAxis3 ()
 Computes the Angle flexural inertia. More...
 
double GetInertiaAxis23 ()
 Computes the Angle product of inertia. More...
 
void ComputeSectionCenter (double &zcm, double &ycm)
 Gets the section centroid. More...
 
- Protected Member Functions inherited from Section
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

double h
 Total height. More...
 
double b
 Total width. More...
 
double tw
 Web thickness. More...
 
double tf
 Flange thickness. More...
 
double Theta
 Rotation angle. More...
 
unsigned int InsertPoint
 Local axis coordinate. More...
 
Eigen::VectorXd Strain
 Generalized Strain vector. More...
 
std::unique_ptr< MaterialtheMaterial
 The section's material. More...
 

Detailed Description

Class for creating a angle section for 3D analysis with linear elastic material for a frame elements.

See also
Section.hpp

Constructor & Destructor Documentation

◆ Lin3DAngle()

Lin3DAngle::Lin3DAngle ( double  h,
double  b,
double  tw,
double  tf,
std::unique_ptr< Material > &  material,
double  theta = 0.0,
unsigned int  ip = 10 
)

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

Parameters
hThe total height of the section.
bThe total width of the section.
twThe thickness of the web.
tfThe thickness of the flange.
materialThe Material that this Section is made out of.
thetaThe angle of rotation about x3-axis.
ipThe insertion point where the section axes are located.
Note
Details about insertion point and local axis rotation can be found at Insertion Point and Local Axes.
See also
Lin3DAngle::h, Lin3DAngle::b, Lin3DAngle::tw, Lin3DAngle::tf, Lin3DAngle::theMaterial, Lin3DAngle::Theta, Lin3DAngle::InsertPoint.

◆ ~Lin3DAngle()

Lin3DAngle::~Lin3DAngle ( )

Destroys this Lin3DAngle object.

Member Function Documentation

◆ CommitState()

void Lin3DAngle::CommitState ( )
virtual

Perform converged section state update.

Note
This function sets the trail stress and strain as converged.

Implements Section.

◆ ComputeSectionCenter()

void Lin3DAngle::ComputeSectionCenter ( double &  zcm,
double &  ycm 
)
protected

Gets the section centroid.

Parameters
zcmThe centroid position along x3 axis.
ycmThe centroid position along x2 axis.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ CopySection()

std::unique_ptr<Section> Lin3DAngle::CopySection ( )
virtual

Clone the 'Lin3DAngle' section.

Returns
A Section pointer to the derived class.
See also
Section.

Implements Section.

◆ GetArea()

double Lin3DAngle::GetArea ( )
protected

Computes the Angle area.

Returns
A constant with the cross-section area.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetDensity()

Eigen::MatrixXd Lin3DAngle::GetDensity ( )
virtual

Access the section density matrix.

Returns
The section density matrix.

Implements Section.

◆ GetInertiaAxis1()

double Lin3DAngle::GetInertiaAxis1 ( )
protected
Returns
A constant with polar moment of inertia about axis x1.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetInertiaAxis2()

double Lin3DAngle::GetInertiaAxis2 ( )
protected

Computes the Angle flexural inertia.

Returns
A constant with the moment of inertia about axis x2.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetInertiaAxis23()

double Lin3DAngle::GetInertiaAxis23 ( )
protected

Computes the Angle product of inertia.

Returns
A constant with the product of inertia.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetInertiaAxis3()

double Lin3DAngle::GetInertiaAxis3 ( )
protected

Computes the Angle flexural inertia.

Returns
A constant with the moment of inertia about axis x3.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetInitialTangentStiffness()

Eigen::MatrixXd Lin3DAngle::GetInitialTangentStiffness ( )
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.

Implements Section.

◆ GetShearArea2()

double Lin3DAngle::GetShearArea2 ( )
protected

Computes the Angle shear area.

Returns
A constant with the shear area along axis x2.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetShearArea3()

double Lin3DAngle::GetShearArea3 ( )
protected

Computes the Angle shear area.

Returns
A constant with the shear area along axis x3.
See also
The geometric properties can be revisited at Lin3DAngle.

◆ GetStrain()

Eigen::VectorXd Lin3DAngle::GetStrain ( )
virtual

Returns the resultant (generalised) strain vector over the section.

Returns
Vector with the resultant strain components.

Implements Section.

◆ GetStrainAt()

Eigen::VectorXd Lin3DAngle::GetStrainAt ( double  x3 = 0,
double  x2 = 0 
)
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).

Implements Section.

◆ GetStress()

Eigen::VectorXd Lin3DAngle::GetStress ( )
virtual

Returns the resultant (generalised) stress vector over the section.

Returns
Vector with the resultant stress components.

Implements Section.

◆ GetStressAt()

Eigen::VectorXd Lin3DAngle::GetStressAt ( double  x3 = 0,
double  x2 = 0 
)
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).

Implements Section.

◆ GetTangentStiffness()

Eigen::MatrixXd Lin3DAngle::GetTangentStiffness ( )
virtual

Returns the section stiffness matrix.

Returns
Matrix with the section consistent tangent stiffness matrix.

Implements Section.

◆ InitialState()

void Lin3DAngle::InitialState ( )
virtual

Brings the section states to its initial state.

Note
This function returns the material states to the beginning.

Implements Section.

◆ ReverseState()

void Lin3DAngle::ReverseState ( )
virtual

Reverse the section states to previous converged state.

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

Implements Section.

◆ UpdateState()

void Lin3DAngle::UpdateState ( const Eigen::VectorXd  strain,
const unsigned int  cond 
)
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.

Implements Section.

Member Data Documentation

◆ b

double Lin3DAngle::b
private

Total width.

◆ h

double Lin3DAngle::h
private

Total height.

◆ InsertPoint

unsigned int Lin3DAngle::InsertPoint
private

Local axis coordinate.

◆ Strain

Eigen::VectorXd Lin3DAngle::Strain
private

Generalized Strain vector.

◆ tf

double Lin3DAngle::tf
private

Flange thickness.

◆ theMaterial

std::unique_ptr<Material> Lin3DAngle::theMaterial
private

The section's material.

◆ Theta

double Lin3DAngle::Theta
private

Rotation angle.

◆ tw

double Lin3DAngle::tw
private

Web thickness.