Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Section.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:
29 //------------------------------------------------------------------------------
30 
31 #ifndef _SECTION_HPP_
32 #define _SECTION_HPP_
33 
34 #include <string>
35 #include <memory>
36 #include <Eigen/Dense>
37 
45 class Section {
46 
47  public:
51  Section(std::string name);
52 
54  virtual ~Section() = 0;
55 
59  virtual std::unique_ptr<Section> CopySection() = 0;
60 
63  virtual Eigen::VectorXd GetStrain() = 0;
64 
67  virtual Eigen::VectorXd GetStress() = 0;
68 
71  virtual Eigen::MatrixXd GetDensity() = 0;
72 
75  virtual Eigen::MatrixXd GetTangentStiffness() = 0;
76 
80  virtual Eigen::MatrixXd GetInitialTangentStiffness() = 0;
81 
86  virtual Eigen::VectorXd GetStrainAt(double x3=0, double x2=0) = 0;
87 
92  virtual Eigen::VectorXd GetStressAt(double x3=0, double x2=0) = 0;
93 
96  virtual void CommitState() = 0;
97 
100  virtual void ReverseState() = 0;
101 
104  virtual void InitialState() = 0;
105 
110  virtual void UpdateState(Eigen::VectorXd strain, unsigned int cond=0) = 0;
111 
114  std::string GetName();
115 
116  protected:
121  Eigen::MatrixXd GetLineRotationMatrix(double theta);
122 
127  Eigen::MatrixXd GetAreaRotationMatrix(double Theta);
128 
137  Eigen::MatrixXd GetLineTranslationMatrix(double h, double b, double zc, double yc, unsigned int ip);
138 
149  void InsertionPointCoordinates(double &x3, double &x2, double h, double b, double zc, double yc, unsigned int ip);
150 
159  Eigen::MatrixXd ComputeLineLocalAxes(double h, double b, double zcm, double ycm, double angle, unsigned int ip);
160 
161  private:
163  std::string Name;
164 };
165 
166 #endif
virtual Eigen::MatrixXd GetDensity()=0
Access the section density matrix.
std::string Name
Section Name.
Definition: Section.hpp:163
Eigen::MatrixXd ComputeLineLocalAxes(double h, double b, double zcm, double ycm, double angle, unsigned int ip)
Transforms generalised strain/stresses from Element to Section local coordinate.
virtual void InitialState()=0
Brings the section states to its initial state.
Virtual class for creating a section object.
Definition: Section.hpp:45
virtual void ReverseState()=0
Reverse the section states to previous converged state.
virtual Eigen::VectorXd GetStrain()=0
Returns the resultant (generalised) strain vector over the section.
virtual void CommitState()=0
Perform converged section state update.
virtual Eigen::VectorXd GetStrainAt(double x3=0, double x2=0)=0
Returns the section strain at given position.
virtual std::unique_ptr< Section > CopySection()=0
Clone the selected Section.
Eigen::MatrixXd GetAreaRotationMatrix(double Theta)
Gets Local Axis Rotation for Area Section according to provided angle.
Eigen::MatrixXd GetLineTranslationMatrix(double h, double b, double zc, double yc, unsigned int ip)
Gets centroid translation axis for Line Section according to insertion point.
Eigen::MatrixXd GetLineRotationMatrix(double theta)
Gets Local Axis Rotation for Line Section according to provided angle.
virtual ~Section()=0
Destroys this Section object.
std::string GetName()
Gets Section name.
virtual void UpdateState(Eigen::VectorXd strain, unsigned int cond=0)=0
Update the section state for this iteration.
Section(std::string name)
Creates a Section to be specified at a Gauss-point in an Element.
void InsertionPointCoordinates(double &x3, double &x2, double h, double b, double zc, double yc, unsigned int ip)
Gets the coordinate according to insertion point.
virtual Eigen::VectorXd GetStressAt(double x3=0, double x2=0)=0
Returns the section stress at given position.
virtual Eigen::VectorXd GetStress()=0
Returns the resultant (generalised) stress vector over the section.
virtual Eigen::MatrixXd GetTangentStiffness()=0
Returns the section stiffness matrix.
virtual Eigen::MatrixXd GetInitialTangentStiffness()=0
Returns the section initial stiffness matrix.