32 #ifndef _LIN2DANGLE_HPP_ 33 #define _LIN2DANGLE_HPP_ 37 #include <Eigen/Dense> 62 Lin2DAngle(
double h,
double b,
double tw,
double tf, std::unique_ptr<Material> &material,
double theta=0.0,
unsigned int ip=10);
97 Eigen::VectorXd
GetStrainAt(
double x3=0,
double x2=0);
103 Eigen::VectorXd
GetStressAt(
double x3=0,
double x2=0);
121 void UpdateState(
const Eigen::VectorXd strain,
const unsigned int cond);
double b
Total width.
Definition: Lin2DAngle.hpp:165
void ComputeSectionCenter(double &zcm, double &ycm)
Gets the section centroid.
Virtual class for creating a section object.
Definition: Section.hpp:45
double GetInertiaAxis23()
Computes the Angle product of inertia.
Class for creating an angle section for 2D analysis with linear elastic material for a frame elements...
Definition: Lin2DAngle.hpp:49
This file contains the abstract "Material object" declarations, which computes the strain...
~Lin2DAngle()
Destroys this Lin2DAngle object.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the section state for this iteration.
Eigen::VectorXd GetStress()
Returns the resultant (generalised) stress vector over the section.
Eigen::VectorXd GetStrainAt(double x3=0, double x2=0)
Returns the section strain at given position.
This file contains the abstract "Section object" declarations, which computes the strain...
double GetInertiaAxis3()
Computes the Angle flexural inertia.
void InitialState()
Brings the section states to its initial state.
Eigen::MatrixXd GetDensity()
Access the section density matrix.
Eigen::MatrixXd GetTangentStiffness()
Returns the section stiffness matrix.
std::unique_ptr< Section > CopySection()
Clone the 'Lin2DAngle' section.
Lin2DAngle(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.
double GetArea()
Computes the Angle area.
void ReverseState()
Reverse the section states to previous converged state.
void CommitState()
Perform converged section state update.
double tw
Web thickness.
Definition: Lin2DAngle.hpp:168
Eigen::VectorXd GetStressAt(double x3=0, double x2=0)
Returns the section stress at given position.
double GetInertiaAxis2()
Computes the Angle flexural inertia.
Eigen::VectorXd GetStrain()
Returns the resultant (generalised) strain vector over the section.
double h
Total height.
Definition: Lin2DAngle.hpp:162
double Theta
Rotation angle.
Definition: Lin2DAngle.hpp:174
double GetShearArea2()
Computes the Angle shear area.
unsigned int InsertPoint
Local axis coordinate.
Definition: Lin2DAngle.hpp:177
Eigen::VectorXd Strain
Generalized Strain vector.
Definition: Lin2DAngle.hpp:180
Eigen::MatrixXd GetInitialTangentStiffness()
Returns the section initial stiffness matrix.
double GetShearArea3()
Computes the Angle shear area.
double tf
Flange thickness.
Definition: Lin2DAngle.hpp:171
std::unique_ptr< Material > theMaterial
The section's material.
Definition: Lin2DAngle.hpp:183