Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Lin2DAngle.hpp
Go to the documentation of this file.
1 //==============================================================================
2 //
3 // Seismo Virtual Laboratory
4 // Module for Serial and Parallel Analysis of seismic
5 // wave propagation and soil-structure interaction simulation
6 // Copyright (C) 2018-2021, The California Institute of Technology
7 // All Rights Reserved.
8 //
9 // Commercial use of this program without express permission of the California
10 // Institute of Technology, is strictly prohibited. See file "COPYRIGHT" in
11 // main directory for information on usage and redistribution, and for a
12 // DISCLAIMER OF ALL WARRANTIES.
13 //
14 //==============================================================================
15 //
16 // Written by:
17 // Danilo S. Kusanovic (dkusanov@caltech.edu)
18 // Elnaz E. Seylabi (elnaze@unr.edu)
19 //
20 // Supervised by:
21 // Domniki M. Asimaki (domniki@caltech.edu)
22 //
23 // References :
24 // [1]
25 //
26 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _LIN2DANGLE_HPP_
33 #define _LIN2DANGLE_HPP_
34 
35 #include <string>
36 #include <memory>
37 #include <Eigen/Dense>
38 
39 #include "Section.hpp"
40 #include "Material.hpp"
41 
49 class Lin2DAngle : public Section{
50 
51  public:
62  Lin2DAngle(double h, double b, double tw, double tf, std::unique_ptr<Material> &material, double theta=0.0, unsigned int ip=10);
63 
65  ~Lin2DAngle();
66 
70  std::unique_ptr<Section> CopySection();
71 
74  Eigen::VectorXd GetStrain();
75 
78  Eigen::VectorXd GetStress();
79 
82  Eigen::MatrixXd GetDensity();
83 
86  Eigen::MatrixXd GetTangentStiffness();
87 
91  Eigen::MatrixXd GetInitialTangentStiffness();
92 
97  Eigen::VectorXd GetStrainAt(double x3=0, double x2=0);
98 
103  Eigen::VectorXd GetStressAt(double x3=0, double x2=0);
104 
107  void CommitState();
108 
111  void ReverseState();
112 
115  void InitialState();
116 
121  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
122 
123  protected:
127  double GetArea();
128 
132  double GetShearArea2();
133 
137  double GetShearArea3();
138 
142  double GetInertiaAxis2();
143 
147  double GetInertiaAxis3();
148 
152  double GetInertiaAxis23();
153 
158  void ComputeSectionCenter(double &zcm, double &ycm);
159 
160  private:
162  double h;
163 
165  double b;
166 
168  double tw;
169 
171  double tf;
172 
174  double Theta;
175 
177  unsigned int InsertPoint;
178 
180  Eigen::VectorXd Strain;
181 
183  std::unique_ptr<Material> theMaterial;
184 };
185 
186 #endif
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 &#39;Lin2DAngle&#39; 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&#39;s material.
Definition: Lin2DAngle.hpp:183