Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Elastic3DLinear.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] Finite Element Procedures, Bathe, K.J., Chapter 4: pages 194, Table 4.3,
25 // Prentice-Hall, 1996.
26 //
27 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _ELASTIC3DLINEAR_HPP_
33 #define _ELASTIC3DLINEAR_HPP_
34 
35 #include <string>
36 #include <memory>
37 #include <Eigen/Dense>
38 
39 #include "Material.hpp"
40 
48 class Elastic3DLinear : public Material{
49 
50  public:
56  Elastic3DLinear(const double E, const double nu, const double rho = 0.0);
57 
60 
61  //Clone the 'Elastic3DLinear' material.
64  std::unique_ptr<Material> CopyMaterial();
65 
69  double GetDensity() const;
70 
74  double GetPoissonRatio() const;
75 
78  double GetBulkModulus() const;
79 
82  double GetShearModulus() const;
83 
87  double GetElasticityModulus() const;
88 
91  double GetEnergy() const;
92 
95  Eigen::MatrixXd GetDamping() const;
96 
100  Eigen::VectorXd GetStrain() const;
101 
105  Eigen::VectorXd GetStress() const;
106 
110  Eigen::VectorXd GetStrainRate() const;
111 
114  Eigen::VectorXd GetTotalStress() const;
115 
119  Eigen::MatrixXd GetTangentStiffness() const;
120 
124  Eigen::MatrixXd GetInitialTangentStiffness() const;
125 
128  void CommitState();
129 
132  void ReverseState();
133 
136  void InitialState();
137 
142  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
143 
144  private:
146  double E;
147 
149  double nu;
150 
152  double Rho;
153 
155  Eigen::VectorXd Strain;
156 
158  Eigen::VectorXd newStrain;
159 
161  Eigen::MatrixXd TangentStiffness;
162 };
163 
164 #endif
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
double GetElasticityModulus() const
Access modulus of elasticity.
void CommitState()
Perform converged material state update.
Eigen::VectorXd GetStress() const
Returns the material stress.
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
double GetDensity() const
Access material density.
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
This file contains the abstract "Material object" declarations, which computes the strain...
std::unique_ptr< Material > CopyMaterial()
Class for creating a triaxial isotropic linear elastic material in for tri-dimensional elements...
Definition: Elastic3DLinear.hpp:48
double E
Modulus of elasticity.
Definition: Elastic3DLinear.hpp:146
Eigen::VectorXd newStrain
Commited Strain vector.
Definition: Elastic3DLinear.hpp:158
~Elastic3DLinear()
Destroys this Elastic3DLinear material.
double Rho
Material density.
Definition: Elastic3DLinear.hpp:152
double GetShearModulus() const
Access shear modulus.
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
Elastic3DLinear(const double E, const double nu, const double rho=0.0)
Creates a Elastic3DLinear material to be specified at a Gauss-point in an Element.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the material state for this iteration.
Eigen::VectorXd Strain
Strain vector.
Definition: Elastic3DLinear.hpp:155
void ReverseState()
Reverse the material states to previous converged state.
double GetEnergy() const
Access the material&#39;s energy at current strain.
Virtual class for creating a material object.
Definition: Material.hpp:45
double GetBulkModulus() const
Access bulk modulus.
Eigen::VectorXd GetStrain() const
Returns the material strain.
double nu
Poisson&#39;s ratio.
Definition: Elastic3DLinear.hpp:149
Eigen::MatrixXd TangentStiffness
Tangent stiffness matrix.
Definition: Elastic3DLinear.hpp:161
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.
void InitialState()
Brings the material states to its initial state in the element.