Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Elastic1DGap.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]
25 //
26 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _ELASTIC1DGAP_HPP_
33 #define _ELASTIC1DGAP_HPP_
34 
35 #include "Material.hpp"
36 
44 class Elastic1DGap : public Material {
45 
46  public:
52  Elastic1DGap(double E, double gap, bool behavior=false);
53 
55  ~Elastic1DGap();
56 
57  //Clone the 'Elastic1DGap' material.
60  std::unique_ptr<Material> CopyMaterial();
61 
64  double GetDensity() const;
65 
68  double GetPoissonRatio() const;
69 
72  double GetBulkModulus() const;
73 
76  double GetShearModulus() const;
77 
81  double GetElasticityModulus() const;
82 
85  double GetEnergy() const;
86 
89  Eigen::MatrixXd GetDamping() const;
90 
95  Eigen::VectorXd GetStrain() const;
96 
101  Eigen::VectorXd GetStress() const;
102 
106  Eigen::VectorXd GetStrainRate() const;
107 
112  Eigen::VectorXd GetTotalStress() const;
113 
118  Eigen::MatrixXd GetTangentStiffness() const;
119 
123  Eigen::MatrixXd GetInitialTangentStiffness() const;
124 
127  void CommitState();
128 
131  void ReverseState();
132 
135  void InitialState();
136 
141  void UpdateState(const Eigen::VectorXd strain, const unsigned int cond=0);
142 
143  private:
145  double E;
146 
148  double Gap;
149 
151  bool Behavior;
152 
154  double oldStrain;
155 
157  double newStrain;
158 };
159 
160 #endif
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
double GetEnergy() const
Access the material&#39;s energy at current strain.
Eigen::VectorXd GetStress() const
Returns the material stress.
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond=0)
Update the material state for this iteration.
std::unique_ptr< Material > CopyMaterial()
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
double GetBulkModulus() const
Access bulk modulus.
bool Behavior
The material behavior (true: Tension, false: Compression)
Definition: Elastic1DGap.hpp:151
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
void CommitState()
Perform converged material state update.
Eigen::VectorXd GetStrain() const
Returns the material strain.
Elastic1DGap(double E, double gap, bool behavior=false)
Creates a Elastic1DGap material to be specified at a Gauss-point in an Element.
double newStrain
Commited Elastic Gap history variables:
Definition: Elastic1DGap.hpp:157
~Elastic1DGap()
Destroys this Elastic1DGap material.
Class for creating an uniaxial elastic gap material compatible with zero-length element.
Definition: Elastic1DGap.hpp:44
Virtual class for creating a material object.
Definition: Material.hpp:45
double GetElasticityModulus() const
Access modulus of elasticity.
double Gap
The separation from which the material start reacting.
Definition: Elastic1DGap.hpp:148
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
void InitialState()
Brings the material states to its initial state in the element.
void ReverseState()
Reverse the material states to previous converged state.
double GetShearModulus() const
Returns linear shear modulus.
double GetDensity() const
Access material density.
double E
The elastic moduli of this material.
Definition: Elastic1DGap.hpp:145
double oldStrain
Elastic Gap history variables:
Definition: Elastic1DGap.hpp:154
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.