Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Plastic1DJ2.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] J.C. Simo, T.J.R. Hughes, Computational Inelasticity, Springer, 1997,
25 // pp.45
26 //
27 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _PLASTIC1DJ2_HPP_
33 #define _PLASTIC1DJ2_HPP_
34 
35 #include <string>
36 #include <memory>
37 #include <Eigen/Dense>
38 
39 #include "Material.hpp"
40 
48 class Plastic1DJ2 : public Material{
49 
50  public:
59  Plastic1DJ2(const double E, const double nu, const double rho, const double K, const double H, const double SigmaY);
60 
62  ~Plastic1DJ2();
63 
64  //Clone the 'Plastic1DJ2' material.
67  std::unique_ptr<Material> CopyMaterial();
68 
72  double GetDensity() const;
73 
77  double GetPoissonRatio() const;
78 
81  double GetBulkModulus() const;
82 
85  double GetShearModulus() const;
86 
90  double GetElasticityModulus() const;
91 
94  double GetEnergy() const;
95 
99  Eigen::MatrixXd GetDamping() const;
100 
105  Eigen::VectorXd GetStrain() const;
106 
111  Eigen::VectorXd GetStress() const;
112 
117  Eigen::VectorXd GetStrainRate() const;
118 
122  Eigen::VectorXd GetTotalStress() const;
123 
128  Eigen::MatrixXd GetTangentStiffness() const;
129 
133  Eigen::MatrixXd GetInitialTangentStiffness() const;
134 
137  void CommitState();
138 
141  void ReverseState();
142 
145  void InitialState();
146 
151  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
152 
153  private:
155  double E;
156 
158  double nu;
159 
161  double Rho;
162 
164  double K;
165 
167  double H;
168 
170  double SigmaY;
171 
173  double alpha;
174 
176  double alpha_n;
177 
179  Eigen::VectorXd Strain;
180 
182  Eigen::VectorXd Stress;
183 
185  Eigen::VectorXd BackStress;
186 
188  Eigen::VectorXd PlasticStrain;
189 
191  Eigen::MatrixXd TangentStiffness;
192 
194  Eigen::VectorXd BackStress_n;
195 
197  Eigen::VectorXd PlasticStrain_n;
198 
200  Eigen::VectorXd Strain_n;
201 
203  Eigen::VectorXd Stress_n;
204 
206  Eigen::MatrixXd TangentStiffness_n;
207 };
208 
209 #endif
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
Eigen::VectorXd PlasticStrain_n
Plastic strain at time n.
Definition: Plastic1DJ2.hpp:197
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
Eigen::VectorXd PlasticStrain
Plastic strain.
Definition: Plastic1DJ2.hpp:188
Class for creating an uniaxial isotropic plastic material for one-dimensional elements.
Definition: Plastic1DJ2.hpp:48
void InitialState()
Brings the material states to its initial state in the element.
std::unique_ptr< Material > CopyMaterial()
Eigen::VectorXd BackStress
Back stress.
Definition: Plastic1DJ2.hpp:185
Eigen::VectorXd GetStrain() const
Returns the material strain.
double nu
Poisson&#39;s ratio.
Definition: Plastic1DJ2.hpp:158
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::VectorXd Strain
Strain vector.
Definition: Plastic1DJ2.hpp:179
double GetDensity() const
Access material density.
Eigen::VectorXd Stress
Stress vector.
Definition: Plastic1DJ2.hpp:182
double E
Modulus of elasticity.
Definition: Plastic1DJ2.hpp:155
Eigen::VectorXd BackStress_n
Back stress at time n.
Definition: Plastic1DJ2.hpp:194
double GetEnergy() const
Access the material&#39;s energy at current strain.
Eigen::VectorXd Strain_n
Strain at time n.
Definition: Plastic1DJ2.hpp:200
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
double Rho
Material density.
Definition: Plastic1DJ2.hpp:161
double GetBulkModulus() const
Access bulk modulus.
void CommitState()
Perform converged material state update.
Virtual class for creating a material object.
Definition: Material.hpp:45
double H
Kinematic hardening modulus.
Definition: Plastic1DJ2.hpp:167
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
Eigen::MatrixXd TangentStiffness_n
Consistent tangent stiffness at time n.
Definition: Plastic1DJ2.hpp:206
~Plastic1DJ2()
Destroys this Plastic1DJ2 material.
double alpha_n
Internal hardening variable at time n.
Definition: Plastic1DJ2.hpp:176
double alpha
Internal hardening variable.
Definition: Plastic1DJ2.hpp:173
double GetElasticityModulus() const
Access modulus of elasticity.
double SigmaY
Yield stress.
Definition: Plastic1DJ2.hpp:170
Eigen::VectorXd Stress_n
Stress vector at time n.
Definition: Plastic1DJ2.hpp:203
double K
Plastic modulus.
Definition: Plastic1DJ2.hpp:164
Eigen::MatrixXd TangentStiffness
Consistent tangent stiffness.
Definition: Plastic1DJ2.hpp:191
Plastic1DJ2(const double E, const double nu, const double rho, const double K, const double H, const double SigmaY)
Creates a Plastic1DJ2 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.
double GetShearModulus() const
Access shear modulus.
void ReverseState()
Reverse the material states to previous converged state.
Eigen::VectorXd GetStress() const
Returns the material stress.
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.