Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Plastic3DJ2.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.124
26 // [2] R.I. Borja, Plasticity Modeling & Computation, Springer, 2013, pp. 47
27 //
28 // Description:
31 //------------------------------------------------------------------------------
32 
33 #ifndef _PLASTIC3DJ2_HPP_
34 #define _PLASTIC3DJ2_HPP_
35 
36 #include <string>
37 #include <memory>
38 #include <Eigen/Dense>
39 
40 #include "Material.hpp"
41 
49 class Plastic3DJ2 : public Material{
50 
51  public:
60  Plastic3DJ2(const double K, const double G, const double rho, const double Hbar, const double beta, const double SigmaY);
61 
63  ~Plastic3DJ2();
64 
65  //Clone the 'Plastic3DJ2' material.
68  std::unique_ptr<Material> CopyMaterial();
69 
73  double GetDensity() const;
74 
77  double GetPoissonRatio() const;
78 
82  double GetBulkModulus() const;
83 
87  double GetShearModulus() const;
88 
91  double GetElasticityModulus() const;
92 
95  double GetEnergy() const;
96 
99  Eigen::MatrixXd GetDamping() const;
100 
105  Eigen::VectorXd GetStrain() const;
106 
111  Eigen::VectorXd GetStress() const;
112 
116  Eigen::VectorXd GetStrainRate() const;
117 
121  Eigen::VectorXd GetTotalStress() const;
122 
127  Eigen::MatrixXd GetTangentStiffness() const;
128 
132  Eigen::MatrixXd GetInitialTangentStiffness() const;
133 
136  void CommitState();
137 
140  void ReverseState();
141 
144  void InitialState();
145 
150  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
151 
152  protected:
155  Eigen::VectorXd ComputeIdentityVector() const;
156 
159  Eigen::MatrixXd ComputeDeviatoricTensor() const;
160 
163  Eigen::MatrixXd ComputeIdentityTensor() const;
164 
167  Eigen::MatrixXd ComputeIdentityOperator() const;
168 
172  double ComputeTensorTrace(const Eigen::VectorXd &T);
173 
177  double ComputeTensorNorm(const Eigen::VectorXd &T);
178 
179  private:
181  double K;
182 
184  double G;
185 
187  double Rho;
188 
190  double H;
191 
193  double beta;
194 
196  double SigmaY;
197 
199  double alpha;
200 
202  double alpha_n;
203 
205  Eigen::VectorXd Strain;
206 
208  Eigen::VectorXd Stress;
209 
211  Eigen::VectorXd BackStress;
212 
214  Eigen::VectorXd PlasticStrain;
215 
217  Eigen::MatrixXd TangentStiffness;
218 
220  Eigen::VectorXd BackStress_n;
221 
223  Eigen::VectorXd PlasticStrain_n;
224 
226  Eigen::VectorXd Strain_n;
227 
229  Eigen::VectorXd Stress_n;
230 
232  Eigen::MatrixXd TangentStiffness_n;
233 };
234 
235 #endif
std::unique_ptr< Material > CopyMaterial()
Eigen::MatrixXd TangentStiffness_n
Consistent tangent stiffness at time n.
Definition: Plastic3DJ2.hpp:232
Class for creating a triaxial J2 plastic material for three-dimensional elements. ...
Definition: Plastic3DJ2.hpp:49
Eigen::VectorXd GetStrain() const
Returns the material strain.
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.
Eigen::MatrixXd ComputeDeviatoricTensor() const
Constructs deviatoric tensor.
Plastic3DJ2(const double K, const double G, const double rho, const double Hbar, const double beta, const double SigmaY)
Creates a Plastic3DJ2 material to be specified at a Gauss-point in an Element.
double alpha_n
Internal hardening variable at time n.
Definition: Plastic3DJ2.hpp:202
void ReverseState()
Reverse the material states to previous converged state.
Eigen::VectorXd PlasticStrain_n
Plastic strain at time n.
Definition: Plastic3DJ2.hpp:223
double ComputeTensorNorm(const Eigen::VectorXd &T)
Computes the norm of a tensor of second-rank.
This file contains the abstract "Material object" declarations, which computes the strain...
void CommitState()
Perform converged material state update.
Eigen::VectorXd BackStress
Trial Back stress.
Definition: Plastic3DJ2.hpp:211
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
double beta
Linear hardening parameter.
Definition: Plastic3DJ2.hpp:193
double GetElasticityModulus() const
Access modulus of elasticity.
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
double SigmaY
Yield stress.
Definition: Plastic3DJ2.hpp:196
Eigen::VectorXd BackStress_n
Back stress at time n.
Definition: Plastic3DJ2.hpp:220
Eigen::VectorXd Strain_n
Strain at time n.
Definition: Plastic3DJ2.hpp:226
Eigen::MatrixXd TangentStiffness
Trial Consistent tangent stiffness.
Definition: Plastic3DJ2.hpp:217
double H
Kinematic hardening modulus.
Definition: Plastic3DJ2.hpp:190
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
Eigen::VectorXd PlasticStrain
Trial Plastic strain.
Definition: Plastic3DJ2.hpp:214
double alpha
Trial Internal hardening variable.
Definition: Plastic3DJ2.hpp:199
Eigen::VectorXd ComputeIdentityVector() const
Constructs a Second Order Identity Tensor.
Eigen::MatrixXd ComputeIdentityOperator() const
Constructs fourth-rank symmetric identity operator.
Eigen::VectorXd Stress_n
Stress vector at time n.
Definition: Plastic3DJ2.hpp:229
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the material state for this iteration.
Virtual class for creating a material object.
Definition: Material.hpp:45
double ComputeTensorTrace(const Eigen::VectorXd &T)
Computes the trace of a tensor of second-rank.
double G
Shear modulus.
Definition: Plastic3DJ2.hpp:184
double GetDensity() const
Access material density.
double GetShearModulus() const
Access shear modulus.
Eigen::VectorXd GetStress() const
Returns the material stress.
double Rho
Material density.
Definition: Plastic3DJ2.hpp:187
~Plastic3DJ2()
Destroys this Plastic3DJ2 material.
Eigen::VectorXd Strain
Trial Strain vector.
Definition: Plastic3DJ2.hpp:205
double GetEnergy() const
Access the material&#39;s energy at current strain.
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
void InitialState()
Brings the material states to its initial state in the element.
double K
Bulk modulus.
Definition: Plastic3DJ2.hpp:181
Eigen::MatrixXd ComputeIdentityTensor() const
Constructs fourth-rank identity tensor.
double GetBulkModulus() const
Access bulk modulus.
Eigen::VectorXd Stress
Trial Stress vector.
Definition: Plastic3DJ2.hpp:208