Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
PlasticPlaneStrainBA.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] R.J. Borja, A.P. Amies (1994), Multiaxial cyclic plasticity model for
25 // clays, ASCE.
26 //
27 // Description:
31 //------------------------------------------------------------------------------
32 
33 #ifndef _PLASTICPLANESTRAINBA_HPP_
34 #define _PLASTICPLANESTRAINBA_HPP_
35 
36 #include <string>
37 #include <Eigen/Dense>
38 #include <iostream>
39 
40 #include "Material.hpp"
41 
50 
51  public:
62  PlasticPlaneStrainBA(const double K, const double G, const double rho, const double H0, const double h, const double m, const double Su, const double beta);
63 
66 
67  //Clone the 'PlasticPlaneStrainBA' material.
70  std::unique_ptr<Material> CopyMaterial();
71 
75  double GetDensity() const;
76 
79  double GetPoissonRatio() const;
80 
84  double GetBulkModulus() const;
85 
89  double GetShearModulus() const;
90 
93  double GetElasticityModulus() const;
94 
97  double GetEnergy() const;
98 
101  Eigen::MatrixXd GetDamping() const;
102 
107  Eigen::VectorXd GetStrain() const;
108 
113  Eigen::VectorXd GetStress() const;
114 
118  Eigen::VectorXd GetStrainRate() const;
119 
123  Eigen::VectorXd GetTotalStress() const;
124 
129  Eigen::MatrixXd GetTangentStiffness() const;
130 
134  Eigen::MatrixXd GetInitialTangentStiffness() const;
135 
139  void CommitState();
140 
143  void ReverseState();
144 
147  void InitialState();
148 
153  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond);
154 
155  protected:
158  Eigen::VectorXd ComputeIdentityVector() const;
159 
162  Eigen::MatrixXd ComputeDeviatoricTensor() const;
163 
166  Eigen::MatrixXd ComputeIdentityTensor() const;
167 
170  Eigen::MatrixXd ComputeIdentityOperator() const;
171 
175  double ComputeTensorTrace(const Eigen::VectorXd &T);
176 
180  double ComputeTensorNorm(const Eigen::VectorXd &T);
181 
186  double ComputeInnerProduct(const Eigen::VectorXd &V1, const Eigen::VectorXd &V2);
187 
192  Eigen::VectorXd ComputeHardening(const Eigen::VectorXd &_DeviatoricIncrStrain, const Eigen::VectorXd &_DeviatoricStress);
193 
198  Eigen::VectorXd ComputeHardeningBisection(const Eigen::VectorXd &_DeviatoricIncrStrain, const Eigen::VectorXd &_DeviatoricStress);
199 
200  private:
202  double K;
203 
205  double G;
206 
208  double Rho;
209 
211  double H0;
212 
214  double h;
215 
217  double m;
218 
220  double Su;
221 
223  double R;
224 
226  double beta;
227 
229  double psi;
230 
232  double kappa;
233 
235  double psi_n;
236 
238  double kappa_n;
239 
241  Eigen::VectorXd Strain;
242 
244  Eigen::VectorXd Strain_n;
245 
247  Eigen::VectorXd Stress;
248 
250  Eigen::VectorXd Stress_n;
251 
253  Eigen::VectorXd DeviatoricStress0;
254 
256  Eigen::VectorXd DeviatoricStress0_n;
257 
259  Eigen::MatrixXd TangentStiffness;
260 
262  Eigen::MatrixXd TangentStiffness_n;
263 
266 
269 
271  int rootFlag;
272 
275 
277  double Hn;
278 };
279 
280 #endif
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
double ComputeInnerProduct(const Eigen::VectorXd &V1, const Eigen::VectorXd &V2)
Compute inner product of two tensor of first-rank.
double G
Shear modulus.
Definition: PlasticPlaneStrainBA.hpp:205
void ReverseState()
Reverse the material states to previous converged state.
Eigen::VectorXd ComputeHardeningBisection(const Eigen::VectorXd &_DeviatoricIncrStrain, const Eigen::VectorXd &_DeviatoricStress)
Compute hardening using the bi-section method.
double kappa_n
Internal hardening parameter at time n.
Definition: PlasticPlaneStrainBA.hpp:238
~PlasticPlaneStrainBA()
Destroys this PlasticPlaneStrainBA material.
Eigen::MatrixXd ComputeIdentityOperator() const
Constructs fourth-rank symmetric identity operator.
double m
Exponential hardening parameter.
Definition: PlasticPlaneStrainBA.hpp:217
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::VectorXd GetStress() const
Returns the material stress.
double GetBulkModulus() const
Access bulk modulus.
double Hn
The hardening function.
Definition: PlasticPlaneStrainBA.hpp:277
Eigen::VectorXd Stress_n
Stress vector at time n.
Definition: PlasticPlaneStrainBA.hpp:250
Eigen::VectorXd DeviatoricStress0_n
Stress tensor at F0 at time n.
Definition: PlasticPlaneStrainBA.hpp:256
double kappa
Internal hardening parameter.
Definition: PlasticPlaneStrainBA.hpp:232
int rootFlag_n
Flag to control the linear system solution at time n.
Definition: PlasticPlaneStrainBA.hpp:274
Eigen::VectorXd ComputeIdentityVector() const
Constructs a Second Order Identity Tensor.
Eigen::VectorXd ComputeHardening(const Eigen::VectorXd &_DeviatoricIncrStrain, const Eigen::VectorXd &_DeviatoricStress)
Compute hardening using the Newton-Rhapson method.
Eigen::VectorXd GetStrain() const
Returns the material strain.
double R
Bounding surface radius.
Definition: PlasticPlaneStrainBA.hpp:223
double psi_n
Internal hardening parameter at time n.
Definition: PlasticPlaneStrainBA.hpp:235
Eigen::VectorXd Strain
Strain vector.
Definition: PlasticPlaneStrainBA.hpp:241
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
double ComputeTensorTrace(const Eigen::VectorXd &T)
Computes the trace of a tensor of second-rank.
double Su
Undrained soil strength.
Definition: PlasticPlaneStrainBA.hpp:220
int rootFlag
Flag to control the linear system solution.
Definition: PlasticPlaneStrainBA.hpp:271
PlasticPlaneStrainBA(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 PlasticPlaneStrainBA material to be specified at a Gauss-point in an Element.
double GetElasticityModulus() const
Access modulus of elasticity.
Eigen::MatrixXd TangentStiffness
Consistent tangent stiffness.
Definition: PlasticPlaneStrainBA.hpp:259
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the material state for this iteration.
double beta
integration factor.
Definition: PlasticPlaneStrainBA.hpp:226
Virtual class for creating a material object.
Definition: Material.hpp:45
double Rho
Material density.
Definition: PlasticPlaneStrainBA.hpp:208
Eigen::MatrixXd ComputeDeviatoricTensor() const
Constructs deviatoric tensor.
Eigen::VectorXd Stress
Trial stress vector.
Definition: PlasticPlaneStrainBA.hpp:247
double GetShearModulus() const
Access shear modulus.
double GetDensity() const
Access material density.
double H0
Bounding surface hardening.
Definition: PlasticPlaneStrainBA.hpp:211
void InitialState()
Brings the material states to its initial state in the element.
Class for creating a biaxial plastic bounding surface material for two-dimensional elements...
Definition: PlasticPlaneStrainBA.hpp:49
double K
Bulk modulus.
Definition: PlasticPlaneStrainBA.hpp:202
double h
exponential hardening parameter.
Definition: PlasticPlaneStrainBA.hpp:214
Eigen::MatrixXd ComputeIdentityTensor() const
Constructs fourth-rank identity tensor.
Eigen::VectorXd DeviatoricStress0
Stress tensor at F0.
Definition: PlasticPlaneStrainBA.hpp:253
int FirstLoadFlag_n
First loading flag at time n.
Definition: PlasticPlaneStrainBA.hpp:268
Eigen::MatrixXd TangentStiffness_n
Consistent tangent stiffness at time n.
Definition: PlasticPlaneStrainBA.hpp:262
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
void CommitState()
Perform converged material state update.
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
double ComputeTensorNorm(const Eigen::VectorXd &T)
Computes the norm of a tensor of second-rank.
double psi
Internal hardening parameter.
Definition: PlasticPlaneStrainBA.hpp:229
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
Eigen::VectorXd Strain_n
Strain at time n.
Definition: PlasticPlaneStrainBA.hpp:244
std::unique_ptr< Material > CopyMaterial()
double GetEnergy() const
Access the material&#39;s energy at current strain.
int FirstLoadFlag
First loading flag.
Definition: PlasticPlaneStrainBA.hpp:265