Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Steel1DFiber.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] Filippou, F. C., Popov, E. P., Bertero, V. V. (1983). "Effects of Bond
25 // Deterioration on Hysteretic Behavior of Reinforced Concrete Joints".
26 // Report EERC 83-19, Earthquake Engineering Research Center, University
27 // of California, Berkeley.
28 //
29 // Description:
33 //------------------------------------------------------------------------------
34 
35 #ifndef _STEEL1DFIBER_HPP_
36 #define _STEEL1DFIBER_HPP_
37 
38 #include "Material.hpp"
39 
47 class Steel1DFiber : public Material {
48 
49  public:
64  Steel1DFiber(double fy, double E, double b, double R0=15.00, double cR1=0.925, double cR2=0.150, double a1=0.0, double a2=1.0, double a3=0.0, double a4=1.0, double nu=0.33, double rho=0.0);
65 
67  ~Steel1DFiber();
68 
69  //Clone the 'Steel1DFiber' material.
72  std::unique_ptr<Material> CopyMaterial();
73 
77  double GetDensity() const;
78 
82  double GetPoissonRatio() const;
83 
86  double GetBulkModulus() const;
87 
90  double GetShearModulus() const;
91 
95  double GetElasticityModulus() const;
96 
99  double GetEnergy() const;
100 
103  Eigen::MatrixXd GetDamping() const;
104 
109  Eigen::VectorXd GetStrain() const;
110 
115  Eigen::VectorXd GetStress() const;
116 
120  Eigen::VectorXd GetStrainRate() const;
121 
126  Eigen::VectorXd GetTotalStress() const;
127 
132  Eigen::MatrixXd GetTangentStiffness() const;
133 
137  Eigen::MatrixXd GetInitialTangentStiffness() const;
138 
141  void CommitState();
142 
145  void ReverseState();
146 
149  void InitialState();
150 
155  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond=0);
156 
157  private:
158  //Steel properties:
159  double nu;
160  double fy;
161  double E;
162  double b;
163  double R0;
164  double cR1;
165  double cR2;
166  double a1;
167  double a2;
168  double a3;
169  double a4;
170  double Rho;
171 
172  //Steel Commited history variables:
173  double newStrain;
174  double newStress;
176  double newMinStrain;
177  double newMaxStrain;
184 
185  //Steel Updated history variables:
186  double oldStrain;
187  double oldStress;
189  double oldMinStrain;
190  double oldMaxStrain;
197 };
198 
199 #endif
void CommitState()
Perform converged material state update.
double nu
Definition: Steel1DFiber.hpp:159
double oldStrain
Definition: Steel1DFiber.hpp:186
double newStressInterception
Definition: Steel1DFiber.hpp:180
double newStrain
Definition: Steel1DFiber.hpp:173
double newTangentStiffness
Definition: Steel1DFiber.hpp:175
double GetEnergy() const
Access the material&#39;s energy at current strain.
double newStress
Definition: Steel1DFiber.hpp:174
Eigen::VectorXd GetStress() const
Returns the material stress.
Eigen::VectorXd GetStrain() const
Returns the material strain.
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
double newPlasticStrain
Definition: Steel1DFiber.hpp:178
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond=0)
Update the material state for this iteration.
This file contains the abstract "Material object" declarations, which computes the strain...
int oldLoadUnloadIndex
Definition: Steel1DFiber.hpp:196
double oldStressInterception
Definition: Steel1DFiber.hpp:193
double a4
Definition: Steel1DFiber.hpp:169
double cR1
Definition: Steel1DFiber.hpp:164
double oldMinStrain
Definition: Steel1DFiber.hpp:189
double GetBulkModulus() const
Access bulk modulus.
double R0
Definition: Steel1DFiber.hpp:163
void InitialState()
Brings the material states to its initial state in the element.
double b
Definition: Steel1DFiber.hpp:162
double oldStrainLastInversion
Definition: Steel1DFiber.hpp:194
double oldStrainInterception
Definition: Steel1DFiber.hpp:192
double newMinStrain
Definition: Steel1DFiber.hpp:176
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
double GetElasticityModulus() const
Access modulus of elasticity.
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
double newMaxStrain
Definition: Steel1DFiber.hpp:177
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
double GetDensity() const
Access material density.
Steel1DFiber(double fy, double E, double b, double R0=15.00, double cR1=0.925, double cR2=0.150, double a1=0.0, double a2=1.0, double a3=0.0, double a4=1.0, double nu=0.33, double rho=0.0)
Creates a Steel1DFiber material to be specified at a Gauss-point in an Element.
double fy
Definition: Steel1DFiber.hpp:160
double cR2
Definition: Steel1DFiber.hpp:165
double oldStressLastInversion
Definition: Steel1DFiber.hpp:195
Virtual class for creating a material object.
Definition: Material.hpp:45
double a2
Definition: Steel1DFiber.hpp:167
Class for creating an uniaxial steel (Menegotto and Pinto 1973) material compatible with fiber sectio...
Definition: Steel1DFiber.hpp:47
double newStressLastInversion
Definition: Steel1DFiber.hpp:182
double a3
Definition: Steel1DFiber.hpp:168
double Rho
Definition: Steel1DFiber.hpp:170
std::unique_ptr< Material > CopyMaterial()
int newLoadUnloadIndex
Definition: Steel1DFiber.hpp:183
double oldMaxStrain
Definition: Steel1DFiber.hpp:190
double E
Definition: Steel1DFiber.hpp:161
double newStrainLastInversion
Definition: Steel1DFiber.hpp:181
double oldTangentStiffness
Definition: Steel1DFiber.hpp:188
double GetShearModulus() const
Returns linear shear modulus.
void ReverseState()
Reverse the material states to previous converged state.
double a1
Definition: Steel1DFiber.hpp:166
double newStrainInterception
Definition: Steel1DFiber.hpp:179
double oldPlasticStrain
Definition: Steel1DFiber.hpp:191
~Steel1DFiber()
Destroys this Steel1DFiber material.
double oldStress
Definition: Steel1DFiber.hpp:187
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.