Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
EQlin2DQuad4.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] M. Darendeli (2001), Development of a new family of normalized modulus
25 // reduction and material damping curves
26 // [2] J. Shi & D. Asimaki (2017), From stiffness to strength: Formulation and
27 // validation of a hybrid hyperbolic nonlinear soil model for siteā€response
28 // analyses
29 //
30 // Description:
33 //------------------------------------------------------------------------------
34 
35 #ifndef _EQLIN2DQUAD4_HPP_
36 #define _EQLIN2DQUAD4_HPP_
37 
38 #include <map>
39 #include <memory>
40 #include <string>
41 #include <Eigen/Dense>
42 
43 #include "Node.hpp"
44 #include "Load.hpp"
45 #include "Material.hpp"
46 #include "Element.hpp"
47 #include "Damping.hpp"
48 #include "QuadratureRule.hpp"
49 
57 class EQlin2DQuad4 : public Element{
58 
59  public:
72  EQlin2DQuad4(const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material, const double th, const std::string quadrature="GAUSS", const unsigned int nGauss=4, const std::string Type="DARENDELI", const double zref=0.0, const double cf1=0.0, const double cf2=0.0);
73 
75  ~EQlin2DQuad4();
76 
79  void CommitState();
80 
83  void ReverseState();
84 
87  void InitialState();
88 
91  void UpdateState();
92 
97  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
98 
102  void SetDamping(const std::shared_ptr<Damping> &damping);
103 
106  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
107 
111  Eigen::MatrixXd GetStrain() const;
112 
116  Eigen::MatrixXd GetStress() const;
117 
121  Eigen::MatrixXd GetStrainRate() const;
122 
128  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
129 
135  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
136 
141  Eigen::VectorXd GetVTKResponse(std::string response) const;
142 
145  double ComputeEnergy();
146 
151  Eigen::MatrixXd ComputeMassMatrix();
152 
157  Eigen::MatrixXd ComputeStiffnessMatrix();
158 
163  Eigen::MatrixXd ComputeDampingMatrix();
164 
169  Eigen::MatrixXd ComputePMLMatrix();
170 
175  Eigen::VectorXd ComputeInternalForces();
176 
181  Eigen::VectorXd ComputeInternalDynamicForces();
182 
189  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
190 
197  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
198 
205  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
206 
207  private:
209  double t;
210 
212  double cf1;
213 
215  double cf2;
216 
218  double zref;
219 
221  std::string Type;
222 
224  std::shared_ptr<Damping> theDamping;
225 
227  std::vector<std::shared_ptr<Node> > theNodes;
228 
230  std::vector<std::unique_ptr<Material> > theMaterial;
231 
233  std::unique_ptr<QuadratureRule> QuadraturePoints;
234 
238  Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const;
239 
243  Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const;
244 
249  Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const;
250 
255  Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const;
256 
262  Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const;
263 
270  Eigen::VectorXd ComputeGGmaxDamping(const double vs, const double z, const double rho, const double gamma) const;
271 };
272 
273 #endif
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const
Evaluates the shape function matrix at a given Gauss point.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
std::vector< std::unique_ptr< Material > > theMaterial
The Element&#39;s material.
Definition: EQlin2DQuad4.hpp:230
double cf2
Upper corner frequencies for damping.
Definition: EQlin2DQuad4.hpp:215
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const
Computes the jacobian of the transformation.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
This file contains the abstract "Material object" declarations, which computes the strain...
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
double zref
Reference elevation to compute GGmax and damping.
Definition: EQlin2DQuad4.hpp:218
Eigen::VectorXd ComputeGGmaxDamping(const double vs, const double z, const double rho, const double gamma) const
Compute GGmax and Damping for nType = 1 ==> equivalent linear model.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Evaluates the strain-displacement matrix at a given Gauss point.
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element&#39;s Nodes.
Definition: EQlin2DQuad4.hpp:227
double t
Element thickness.
Definition: EQlin2DQuad4.hpp:209
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
Class for creating a 2D linearized equivalent linear four-node quadrilateral element in a mesh...
Definition: EQlin2DQuad4.hpp:57
Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const
Update strain in the element.
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
double ComputeEnergy()
Computes the element energy for a given deformation.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
void UpdateState()
Update the material states in the element.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
~EQlin2DQuad4()
Destroys this EQlin2DQuad4 object.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: EQlin2DQuad4.hpp:224
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
void CommitState()
Save the material states in the element.
This file contains the "Node object" declarations, which is stores the coordinates, state variables, degrees-of-freedom, and total-degree of freedom lists of a node in a finite element mesh.
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
EQlin2DQuad4(const std::vector< unsigned int > nodes, std::unique_ptr< Material > &material, const double th, const std::string quadrature="GAUSS", const unsigned int nGauss=4, const std::string Type="DARENDELI", const double zref=0.0, const double cf1=0.0, const double cf2=0.0)
Creates a EQlin2DQuad4 in a finite element Mesh.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: EQlin2DQuad4.hpp:233
double cf1
Lower corner frequencies for damping.
Definition: EQlin2DQuad4.hpp:212
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
void InitialState()
Brings the material/section state to its initial state in this element.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
std::string Type
Determine the type of modulus reduction and damping curve.
Definition: EQlin2DQuad4.hpp:221
This file contains the "Element" object declarations, which defines an element in a finite element me...