Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
lin2DQuad8.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] Finite Element Procedures, Bathe, K.J., Chapter 5: pages 342-348.
25 // Prentice-Hall, 1996.
26 //
27 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _LIN2DQUAD8_HPP_
33 #define _LIN2DQUAD8_HPP_
34 
35 #include <map>
36 #include <memory>
37 #include <string>
38 #include <Eigen/Dense>
39 
40 #include "Node.hpp"
41 #include "Load.hpp"
42 #include "Material.hpp"
43 #include "Element.hpp"
44 #include "Damping.hpp"
45 #include "QuadratureRule.hpp"
46 
54 class lin2DQuad8 : public Element{
55 
56  public:
65  lin2DQuad8(const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material, const double th, const std::string quadrature="GAUSS", const unsigned int nGauss=4);
66 
68  ~lin2DQuad8();
69 
72  void CommitState();
73 
76  void ReverseState();
77 
80  void InitialState();
81 
84  void UpdateState();
85 
90  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
91 
95  void SetDamping(const std::shared_ptr<Damping> &damping);
96 
99  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
100 
104  Eigen::MatrixXd GetStrain() const;
105 
109  Eigen::MatrixXd GetStress() const;
110 
114  Eigen::MatrixXd GetStrainRate() const;
115 
121  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
122 
128  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
129 
134  Eigen::VectorXd GetVTKResponse(std::string response) const;
135 
138  double ComputeEnergy();
139 
144  Eigen::MatrixXd ComputeMassMatrix();
145 
150  Eigen::MatrixXd ComputeStiffnessMatrix();
151 
156  Eigen::MatrixXd ComputeDampingMatrix();
157 
162  Eigen::MatrixXd ComputePMLMatrix();
163 
168  Eigen::VectorXd ComputeInternalForces();
169 
174  Eigen::VectorXd ComputeInternalDynamicForces();
175 
182  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
183 
190  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
191 
198  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
199 
200  private:
202  double t;
203 
205  std::shared_ptr<Damping> theDamping;
206 
208  std::vector<std::shared_ptr<Node> > theNodes;
209 
211  std::vector<std::unique_ptr<Material> > theMaterial;
212 
214  std::unique_ptr<QuadratureRule> QuadraturePoints;
215 
219  Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const;
220 
224  Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const;
225 
230  Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const;
231 
236  Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const;
237 
243  Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const;
244 
248  Eigen::MatrixXd ComputeInitialStiffnessMatrix() const;
249 };
250 
251 #endif
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
lin2DQuad8(const std::vector< unsigned int > nodes, std::unique_ptr< Material > &material, const double th, const std::string quadrature="GAUSS", const unsigned int nGauss=4)
Creates a lin2DQuad8 in a finite element Mesh.
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Class for creating a 2D linearized eight-node quadrilateral element in a mesh.
Definition: lin2DQuad8.hpp:54
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
std::vector< std::unique_ptr< Material > > theMaterial
The Element&#39;s material.
Definition: lin2DQuad8.hpp:211
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
void UpdateState()
Update the material states in the element.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const
Computes the jacobian of the transformation.
This file contains the abstract "Material object" declarations, which computes the strain...
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: lin2DQuad8.hpp:205
void ReverseState()
Reverse the material/section states to previous converged state in this element.
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element&#39;s Nodes.
Definition: lin2DQuad8.hpp:208
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the initial stiffness matrix of the element using gauss-integration.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
~lin2DQuad8()
Destroys this lin2DQuad8 object.
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const
Evaluates the shape function matrix at a given Gauss point.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
double ComputeEnergy()
Computes the element energy for a given deformation.
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: lin2DQuad8.hpp:214
Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Evaluates the strain-displacement matrix at a given Gauss point.
void InitialState()
Brings the material/section state to its initial state in this 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.
void CommitState()
Save the material states in the element.
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const
Update strain in the element.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
This file contains the "Element" object declarations, which defines an element in a finite element me...
double t
Element thickness.
Definition: lin2DQuad8.hpp:202