Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Hertzian1DLinear.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] Tuning of Surface-Acoustic-Wave Dispersion via Magnetically Modulated
25 // Contact Resonances Antonio Palermo, Yifan Wang, Paolo Celli, and Chiara
26 // Daraio Phys. Rev. Applied 11, 044057 – Published 18 April 2019
27 //
28 // Description:
32 //------------------------------------------------------------------------------
33 
34 #ifndef _HERTZIAN1DLINEAR_HPP_
35 #define _HERTZIAN1DLINEAR_HPP_
36 
37 #include <string>
38 #include <memory>
39 #include <Eigen/Dense>
40 
41 #include "Material.hpp"
42 
50 class Hertzian1DLinear : public Material{
51 
52  public:
59  Hertzian1DLinear(const double k1, const double k2, const double k3, const double rho = 0.0);
60 
63 
64  //Clone the 'Hertzian1DLinear' material.
65  std::unique_ptr<Material> CopyMaterial();
66 
70  double GetDensity() const;
71 
74  double GetPoissonRatio() const;
75 
78  double GetBulkModulus() const;
79 
82  double GetShearModulus() const;
83 
86  double GetElasticityModulus() const;
87 
90  double GetEnergy() const;
91 
94  Eigen::MatrixXd GetDamping() const;
95 
99  Eigen::VectorXd GetStrain() const;
100 
104  Eigen::VectorXd GetStress() const;
105 
109  Eigen::VectorXd GetStrainRate() const;
110 
113  Eigen::VectorXd GetTotalStress() const;
114 
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);
142 
143  private:
145  double k1;
146 
148  double k2;
149 
151  double k3;
152 
154  double Rho;
155 
157  Eigen::VectorXd Strain;
158 
160  Eigen::VectorXd newStrain;
161 };
162 
163 #endif
Eigen::VectorXd GetStrainRate() const
Returns the material strain rate.
double GetBulkModulus() const
Access bulk modulus.
void UpdateState(const Eigen::VectorXd strain, const unsigned int cond)
Update the material state for this iteration.
double k1
The linear Spring Constant.
Definition: Hertzian1DLinear.hpp:145
Eigen::VectorXd GetStress() const
Returns the material stress.
std::unique_ptr< Material > CopyMaterial()
Clone the selected material.
void CommitState()
Perform converged material state update.
double Rho
Material density.
Definition: Hertzian1DLinear.hpp:154
Eigen::MatrixXd GetTangentStiffness() const
Returns the material stiffness.
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::VectorXd GetStrain() const
Returns the material strain.
double GetShearModulus() const
Access shear modulus.
Hertzian1DLinear(const double k1, const double k2, const double k3, const double rho=0.0)
Creates a Hertzian1DLinear material to be specified at a Gauss-point in an Element.
Eigen::MatrixXd GetDamping() const
Returns the material viscous damping.
double GetElasticityModulus() const
Access modulus of elasticity.
double GetEnergy() const
Access the material&#39;s energy at current strain.
double k3
The cubic Spring Constant.
Definition: Hertzian1DLinear.hpp:151
double GetPoissonRatio() const
Returns the Poisson&#39;s ratio.
Virtual class for creating a material object.
Definition: Material.hpp:45
Eigen::VectorXd Strain
Strain vector.
Definition: Hertzian1DLinear.hpp:157
Eigen::MatrixXd GetInitialTangentStiffness() const
Returns the initial material stiffness.
~Hertzian1DLinear()
Destroys this Hertzian1DLinear material.
double GetDensity() const
Access material density.
void ReverseState()
Reverse the material states to previous converged state.
Class for creating an uniaxial hertzian contact material for one-dimensional elements.
Definition: Hertzian1DLinear.hpp:50
Eigen::VectorXd newStrain
Commited Strain vector.
Definition: Hertzian1DLinear.hpp:160
void InitialState()
Brings the material states to its initial state in the element.
double k2
The quadratic Spring Constant.
Definition: Hertzian1DLinear.hpp:148
Eigen::VectorXd GetTotalStress() const
Computes the material total stress.