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