32 #ifndef _FIB3DLINESECTION_HPP_ 33 #define _FIB3DLINESECTION_HPP_ 38 #include <Eigen/Dense> 64 Fib3DLineSection(
double h,
double b,
const std::vector<std::unique_ptr<Material> > &fibers,
const std::vector<double> &
zi,
const std::vector<double> &
yi,
const std::vector<double> &
Ai,
double k2=1.0,
double k3=1.0,
double theta=0.0,
unsigned int ip=10);
99 Eigen::VectorXd
GetStrainAt(
double x3=0,
double x2=0);
105 Eigen::VectorXd
GetStressAt(
double x3=0,
double x2=0);
123 void UpdateState(
const Eigen::VectorXd strain,
const unsigned int cond=0);
168 std::vector<double>
zi;
171 std::vector<double>
yi;
174 std::vector<double>
Ai;
Eigen::MatrixXd GetTangentStiffness()
Returns the section stiffness matrix.
std::vector< double > Ai
The discretized area value of each fiber.
Definition: Fib3DLineSection.hpp:174
~Fib3DLineSection()
Destroys this Fib3DLineSection object.
unsigned int FiberIndexAt(double x3=0, double x2=0)
Return the fiber strain closer to (x3, x2)
double zcm
Rotation angle.
Definition: Fib3DLineSection.hpp:144
Eigen::VectorXd GetStrainAt(double x3=0, double x2=0)
Returns the section strain at given position.
Fib3DLineSection(double h, double b, const std::vector< std::unique_ptr< Material > > &fibers, const std::vector< double > &zi, const std::vector< double > &yi, const std::vector< double > &Ai, double k2=1.0, double k3=1.0, double theta=0.0, unsigned int ip=10)
Creates a Section to be specified at a Gauss-point in an Element.
double kappa3
Shear area factor (As3/A) along the x3-axis.
Definition: Fib3DLineSection.hpp:153
Eigen::MatrixXd GetDensity()
Access the section density matrix.
Virtual class for creating a section object.
Definition: Section.hpp:45
void ComputeSectionCenter()
Gets the section centroid.
Eigen::VectorXd GetStrain()
Returns the resultant (generalised) strain vector over the section.
This file contains the abstract "Material object" declarations, which computes the strain...
This file contains the abstract "Section object" declarations, which computes the strain...
double Theta
Rotation angle.
Definition: Fib3DLineSection.hpp:156
double ycm
Rotation angle.
Definition: Fib3DLineSection.hpp:147
Eigen::MatrixXd GetInitialTangentStiffness()
Returns the section initial stiffness matrix.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond=0)
Update the section state for this iteration.
unsigned int InsertPoint
Local axis coordinate.
Definition: Fib3DLineSection.hpp:159
std::unique_ptr< Section > CopySection()
Clone the 'Fib3DLineSection' section.
void ReverseState()
Reverse the section states to previous converged state.
double kappa2
Shear area factor (As2/A) along the x2-axis.
Definition: Fib3DLineSection.hpp:150
unsigned int NumberOfFibers
Local axis coordinate.
Definition: Fib3DLineSection.hpp:162
Eigen::VectorXd GetStress()
Returns the resultant (generalised) stress vector over the section.
std::vector< double > zi
The horizontal (x3) local coordinates of each fiber.
Definition: Fib3DLineSection.hpp:168
void InitialState()
Brings the section states to its initial state.
double h
Height.
Definition: Fib3DLineSection.hpp:138
void CommitState()
Perform converged section state update.
std::vector< std::unique_ptr< Material > > theFiber
The section's fiber list.
Definition: Fib3DLineSection.hpp:177
Eigen::VectorXd GetStressAt(double x3=0, double x2=0)
Returns the section stress at given position.
Eigen::VectorXd Strain
Generalized Strain vector.
Definition: Fib3DLineSection.hpp:165
Class for creating a generic fiber section for 3D analysis with linear/nonlinear material for a frame...
Definition: Fib3DLineSection.hpp:50
std::vector< double > yi
The vertical (x2) local coordinates of each fiber.
Definition: Fib3DLineSection.hpp:171
double b
Width.
Definition: Fib3DLineSection.hpp:141