Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Plastic3DBA.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] Multiaxial Cyclic Plasticity Model for Clays, Ronaldo I. Borja and Alexander
25 // P. Amies, Journal of Geotechnical Engineering, Vol. 120, No. 6, June, 1994
26 //
27 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _PLASTIC3DBA_HPP_
33 #define _PLASTIC3DBA_HPP_
34 
35 #include <string>
36 #include <memory>
37 #include <Eigen/Dense>
38 
39 #include "Material.hpp"
40 
48 class Plastic3DBA : public Material{
49 
50  public:
61  Plastic3DBA(const double K, const double G, const double rho, const double H0, const double h, const double m, const double Su, const double beta);
62 
64  ~Plastic3DBA();
65 
66  //Clone the 'Plastic3DBA' material.
69  std::unique_ptr<Material> CopyMaterial();
70 
74  double GetDensity() const;
75 
78  double GetPoissonRatio() const;
79 
83  double GetBulkModulus() const;
84 
88  double GetShearModulus() const;
89 
92  double GetElasticityModulus() const;
93 
96  double GetEnergy() const;
97 
100  Eigen::MatrixXd GetDamping() const;
101 
106  Eigen::VectorXd GetStrain() const;
107 
112  Eigen::VectorXd GetStress() const;
113 
117  Eigen::VectorXd GetStrainRate() const;
118 
122  Eigen::VectorXd GetTotalStress() const;
123 
128  Eigen::MatrixXd GetTangentStiffness() const;
129 
133  Eigen::MatrixXd GetInitialTangentStiffness() const;
134 
138  void CommitState();
139 
142  void ReverseState();
143 
146  void InitialState();
147 
152  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
153 
154  protected:
157  Eigen::VectorXd ComputeIdentityVector() const;
158 
161  Eigen::MatrixXd ComputeDeviatoricTensor() const;
162 
165  Eigen::MatrixXd ComputeIdentityTensor() const;
166 
169  Eigen::MatrixXd ComputeIdentityOperator() const;
170 
174  double ComputeTensorTrace(const Eigen::VectorXd &T);
175 
179  double ComputeTensorNorm(const Eigen::VectorXd &T);
180 
185  double ComputeInnerProduct(const Eigen::VectorXd &V1, const Eigen::VectorXd &V2);
186 
191  Eigen::VectorXd ComputeHardening(Eigen::VectorXd const &DeviatoricIncrStrain, const Eigen::VectorXd &DeviatoricStress);
192 
197  Eigen::VectorXd ComputeHardeningBisection(Eigen::VectorXd const &DeviatoricIncrStrain, const Eigen::VectorXd &DeviatoricStress);
198 
199  private:
201  double K;
202 
204  double G;
205 
207  double Rho;
208 
210  double H0;
211 
213  double h;
214 
216  double m;
217 
219  double Su;
220 
222  double R;
223 
225  double beta;
226 
228  double psi;
229 
231  double kappa;
232 
234  double psi_n;
235 
237  double kappa_n;
238 
240  Eigen::VectorXd Strain;
241 
243  Eigen::VectorXd Strain_n;
244 
246  Eigen::VectorXd Stress;
247 
249  Eigen::VectorXd Stress_n;
250 
252  Eigen::VectorXd DeviatoricStress0;
253 
255  Eigen::VectorXd DeviatoricStress0_n;
256 
258  Eigen::MatrixXd TangentStiffness;
259 
261  Eigen::MatrixXd TangentStiffness_n;
262 
265 
268 
270  int rootFlag;
271 
274 
276  double Hn;
277 };
278 
279 #endif
double H0
Bounding surface hardening.
Definition: Plastic3DBA.hpp:210
double h
Exponential hardening parameter.
Definition: Plastic3DBA.hpp:213
~Plastic3DBA()
Destroys this Plastic3DBA material.
int FirstLoadFlag
First loading flag.
Definition: Plastic3DBA.hpp:264
Eigen::VectorXd ComputeHardening(Eigen::VectorXd const &DeviatoricIncrStrain, const Eigen::VectorXd &DeviatoricStress)
Compute hardening using the Newton-Rhapson method.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the material state for this iteration.
double m
Exponential hardening parameter.
Definition: Plastic3DBA.hpp:216
double ComputeInnerProduct(const Eigen::VectorXd &V1, const Eigen::VectorXd &V2)
Compute inner product of two tensor of first-rank.
double GetShearModulus() const
Access shear modulus.
Eigen::VectorXd ComputeIdentityVector() const
Constructs a Second Order Identity Tensor.
Class for creating an triaxial plastic bounding surface material for three-dimensional elements...
Definition: Plastic3DBA.hpp:48
Eigen::VectorXd Strain_n
Strain at time n.
Definition: Plastic3DBA.hpp:243
Eigen::VectorXd GetStrain() const
Returns the material strain.
void CommitState()
Perform converged material state update.
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
double psi_n
Internal hardening parameter at time n.
Definition: Plastic3DBA.hpp:234
std::unique_ptr< Material > CopyMaterial()
double GetDensity() const
Access material density.
Plastic3DBA(const double K, const double G, const double rho, const double H0, const double h, const double m, const double Su, const double beta)
Creates a Plastic3DBA material to be specified at a Gauss-point in an Element.
void InitialState()
Brings the material states to its initial state in the element.
Eigen::VectorXd GetStress() const
Returns the material stress.
int rootFlag
Flag to control the linear system solution.
Definition: Plastic3DBA.hpp:270
double K
Bulk modulus.
Definition: Plastic3DBA.hpp:201
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
double kappa_n
Internal hardening parameter at time n.
Definition: Plastic3DBA.hpp:237
double G
Shear modulus.
Definition: Plastic3DBA.hpp:204
Eigen::VectorXd Strain
Strain vector.
Definition: Plastic3DBA.hpp:240
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.
double R
Bounding surface radius.
Definition: Plastic3DBA.hpp:222
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
double kappa
Internal hardening parameters.
Definition: Plastic3DBA.hpp:231
Eigen::VectorXd DeviatoricStress0
Stress tensor at F0.
Definition: Plastic3DBA.hpp:252
Eigen::VectorXd Stress
Trial stress vector.
Definition: Plastic3DBA.hpp:246
double Hn
The hardening function.
Definition: Plastic3DBA.hpp:276
double GetEnergy() const
Access the material&#39;s energy at current strain.
Eigen::MatrixXd ComputeIdentityOperator() const
Constructs fourth-rank symmetric identity operator.
Eigen::MatrixXd TangentStiffness
Consistent tangent stiffness.
Definition: Plastic3DBA.hpp:258
Virtual class for creating a material object.
Definition: Material.hpp:45
Eigen::VectorXd ComputeHardeningBisection(Eigen::VectorXd const &DeviatoricIncrStrain, const Eigen::VectorXd &DeviatoricStress)
Compute hardening using the bi-section method.
Eigen::VectorXd DeviatoricStress0_n
Stress tensor at F0 at time n.
Definition: Plastic3DBA.hpp:255
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
double Su
Undrained soil strength.
Definition: Plastic3DBA.hpp:219
double ComputeTensorNorm(const Eigen::VectorXd &T)
Computes the norm of a tensor of second-rank.
double ComputeTensorTrace(const Eigen::VectorXd &T)
Computes the trace of a tensor of second-rank.
double Rho
Material density.
Definition: Plastic3DBA.hpp:207
int rootFlag_n
Flag to control the linear system solution at time n.
Definition: Plastic3DBA.hpp:273
Eigen::VectorXd Stress_n
Stress vector at time n.
Definition: Plastic3DBA.hpp:249
Eigen::MatrixXd TangentStiffness_n
Consistent tangent stiffness at time n.
Definition: Plastic3DBA.hpp:261
int FirstLoadFlag_n
First loading flag at time n.
Definition: Plastic3DBA.hpp:267
double GetElasticityModulus() const
Access modulus of elasticity.
Eigen::MatrixXd ComputeDeviatoricTensor() const
Constructs deviatoric tensor.
void ReverseState()
Reverse the material states to previous converged state.
Eigen::MatrixXd ComputeIdentityTensor() const
Constructs fourth-rank identity tensor.
double psi
Internal hardening parameters.
Definition: Plastic3DBA.hpp:228
double beta
The integration factor.
Definition: Plastic3DBA.hpp:225
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
double GetBulkModulus() const
Access bulk modulus.