Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Lin2DRectangular.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 _LIN2DRECTANGULAR_HPP_
33 #define _LIN2DRECTANGULAR_HPP_
34 
35 #include <string>
36 #include <memory>
37 #include <Eigen/Dense>
38 
39 #include "Section.hpp"
40 #include "Material.hpp"
41 
49 class Lin2DRectangular : public Section{
50 
51  public:
60  Lin2DRectangular(double h, double b, std::unique_ptr<Material> &material, double theta=0.0, unsigned int ip=10);
61 
64 
68  std::unique_ptr<Section> CopySection();
69 
72  Eigen::VectorXd GetStrain();
73 
76  Eigen::VectorXd GetStress();
77 
80  Eigen::MatrixXd GetDensity();
81 
84  Eigen::MatrixXd GetTangentStiffness();
85 
89  Eigen::MatrixXd GetInitialTangentStiffness();
90 
95  Eigen::VectorXd GetStrainAt(double x3=0, double x2=0);
96 
101  Eigen::VectorXd GetStressAt(double x3=0, double x2=0);
102 
105  void CommitState();
106 
109  void ReverseState();
110 
113  void InitialState();
114 
119  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
120 
121  protected:
125  double GetArea();
126 
130  double GetShearArea2();
131 
135  double GetShearArea3();
136 
140  double GetInertiaAxis2();
141 
145  double GetInertiaAxis3();
146 
151  void ComputeSectionCenter(double &zcm, double &ycm);
152 
153  private:
155  double h;
156 
158  double b;
159 
161  double Theta;
162 
164  unsigned int InsertPoint;
165 
167  Eigen::VectorXd Strain;
168 
170  std::unique_ptr<Material> theMaterial;
171 };
172 
173 #endif
Eigen::MatrixXd GetDensity()
Access the section density matrix.
void InitialState()
Brings the section states to its initial state.
void ComputeSectionCenter(double &zcm, double &ycm)
Gets the section centroid.
double GetInertiaAxis2()
Computes the Angle flexural inertia.
Eigen::VectorXd GetStrain()
Returns the resultant (generalised) strain vector over the section.
Virtual class for creating a section object.
Definition: Section.hpp:45
double GetArea()
Computes the Angle area.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the section state for this iteration.
void ReverseState()
Reverse the section states to previous converged state.
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::MatrixXd GetTangentStiffness()
Returns the section stiffness matrix.
This file contains the abstract "Section object" declarations, which computes the strain...
unsigned int InsertPoint
Local axis coordinate.
Definition: Lin2DRectangular.hpp:164
Eigen::VectorXd Strain
Generalized Strain vector.
Definition: Lin2DRectangular.hpp:167
std::unique_ptr< Material > theMaterial
The section&#39;s material.
Definition: Lin2DRectangular.hpp:170
Class for creating a solid rectangular section for 2D analysis with linear elastic material for a fra...
Definition: Lin2DRectangular.hpp:49
double b
Width.
Definition: Lin2DRectangular.hpp:158
double h
Height.
Definition: Lin2DRectangular.hpp:155
double GetShearArea3()
Computes the Angle shear area.
Eigen::MatrixXd GetInitialTangentStiffness()
Returns the section initial stiffness matrix.
void CommitState()
Perform converged section state update.
double Theta
Rotation angle.
Definition: Lin2DRectangular.hpp:161
double GetInertiaAxis3()
Computes the Angle flexural inertia.
Eigen::VectorXd GetStress()
Returns the resultant (generalised) stress vector over the section.
Eigen::VectorXd GetStressAt(double x3=0, double x2=0)
Returns the section stress at given position.
Eigen::VectorXd GetStrainAt(double x3=0, double x2=0)
Returns the section strain at given position.
Lin2DRectangular(double h, double b, 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.
std::unique_ptr< Section > CopySection()
Clone the &#39;Lin2DRectangular&#39; section.
double GetShearArea2()
Computes the Angle shear area.
~Lin2DRectangular()
Destroys this Lin2DCircular object.