Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Fib3DLineSection.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 _FIB3DLINESECTION_HPP_
33 #define _FIB3DLINESECTION_HPP_
34 
35 #include <string>
36 #include <memory>
37 #include <vector>
38 #include <Eigen/Dense>
39 
40 #include "Material.hpp"
41 #include "Section.hpp"
42 
50 class Fib3DLineSection : public Section{
51 
52  public:
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);
65 
68 
72  std::unique_ptr<Section> CopySection();
73 
76  Eigen::VectorXd GetStrain();
77 
80  Eigen::VectorXd GetStress();
81 
84  Eigen::MatrixXd GetDensity();
85 
88  Eigen::MatrixXd GetTangentStiffness();
89 
93  Eigen::MatrixXd GetInitialTangentStiffness();
94 
99  Eigen::VectorXd GetStrainAt(double x3=0, double x2=0);
100 
105  Eigen::VectorXd GetStressAt(double x3=0, double x2=0);
106 
109  void CommitState();
110 
113  void ReverseState();
114 
117  void InitialState();
118 
123  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond=0);
124 
125  protected:
128  void ComputeSectionCenter();
129 
134  unsigned int FiberIndexAt(double x3=0, double x2=0);
135 
136  private:
138  double h;
139 
141  double b;
142 
144  double zcm;
145 
147  double ycm;
148 
150  double kappa2;
151 
153  double kappa3;
154 
156  double Theta;
157 
159  unsigned int InsertPoint;
160 
162  unsigned int NumberOfFibers;
163 
165  Eigen::VectorXd Strain;
166 
168  std::vector<double> zi;
169 
171  std::vector<double> yi;
172 
174  std::vector<double> Ai;
175 
177  std::vector<std::unique_ptr<Material> > theFiber;
178 };
179 
180 #endif
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 &#39;Fib3DLineSection&#39; 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&#39;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