Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
PML2DQuad4.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] Kucukcoban, S., and Loukas F. Kallivokas. "A symmetric hybrid formulation
25 // for transient wave simulations in PML-truncated heterogeneous media." Wave
26 // Motion Vol 50. Issue 1 (2013): 57-79.
27 //
28 // Description:
32 //------------------------------------------------------------------------------
33 
34 #ifndef _PML2DQUAD4_HPP_
35 #define _PML2DQUAD4_HPP_
36 
37 #include <map>
38 #include <memory>
39 #include <string>
40 #include <Eigen/Dense>
41 
42 #include "Node.hpp"
43 #include "Load.hpp"
44 #include "Material.hpp"
45 #include "Element.hpp"
46 #include "QuadratureRule.hpp"
47 #include "Damping.hpp"
48 
56 class PML2DQuad4 : public Element{
57 
58  public:
67  PML2DQuad4(const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material, const std::vector<double> parameters, const std::string quadrature="GAUSS", const unsigned int nGauss=4);
68 
70  ~PML2DQuad4();
71 
74  void CommitState();
75 
78  void ReverseState();
79 
82  void InitialState();
83 
86  void UpdateState();
87 
92  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
93 
97  void SetDamping(const std::shared_ptr<Damping> &damping);
98 
101  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
102 
106  Eigen::MatrixXd GetStrain() const;
107 
111  Eigen::MatrixXd GetStress() const;
112 
116  Eigen::MatrixXd GetStrainRate() const;
117 
123  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
124 
130  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
131 
136  Eigen::VectorXd GetVTKResponse(std::string response) const;
137 
140  double ComputeEnergy();
141 
146  Eigen::MatrixXd ComputeMassMatrix();
147 
152  Eigen::MatrixXd ComputeStiffnessMatrix();
153 
158  Eigen::MatrixXd ComputeDampingMatrix();
159 
164  Eigen::MatrixXd ComputePMLMatrix();
165 
170  Eigen::VectorXd ComputeInternalForces();
171 
176  Eigen::VectorXd ComputeInternalDynamicForces();
177 
184  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
185 
192  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
193 
200  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
201 
202  private:
204  double t;
205 
207  double m_pml;
208 
210  double L_pml;
211 
213  double R_pml;
214 
216  double x0_pml;
217 
219  double y0_pml;
220 
222  double nx_pml;
223 
225  double ny_pml;
226 
228  std::shared_ptr<Damping> theDamping;
229 
231  std::vector<std::shared_ptr<Node> > theNodes;
232 
234  std::vector<std::unique_ptr<Material> > theMaterial;
235 
237  std::unique_ptr<QuadratureRule> QuadraturePoints;
238 
242  Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const;
243 
247  Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const;
248 
253  Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const;
254 
259  Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const;
260 
266  Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const;
267 
275  Eigen::VectorXd ComputePMLStretchingFactors(const double ri, const double si, const double rho, const double mu, const double lambda) const;
276 };
277 
278 #endif
double R_pml
Reflection coefficient for the PML.
Definition: PML2DQuad4.hpp:213
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const
Evaluates the shape function matrix at a given Gauss point.
void InitialState()
Brings the material/section state to its initial state in this element.
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
double y0_pml
Coordinates y0 of reference point.
Definition: PML2DQuad4.hpp:219
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
~PML2DQuad4()
Destroys this PML2DQuad4 object.
This file contains the abstract "Material object" declarations, which computes the strain...
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: PML2DQuad4.hpp:237
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
Class for creating a 2D linearized four-node perfectly matched layer element in a mesh...
Definition: PML2DQuad4.hpp:56
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
double ComputeEnergy()
Computes the element energy for a given deformation.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: PML2DQuad4.hpp:228
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const
Computes the jacobian of the transformation.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Evaluates the strain-displacement matrix at a given Gauss point.
double ny_pml
Vertical normal components for outward normal.
Definition: PML2DQuad4.hpp:225
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const
Update strain in the element.
double m_pml
Polynomial degree of thePerfectly Match Layer.
Definition: PML2DQuad4.hpp:207
PML2DQuad4(const std::vector< unsigned int > nodes, std::unique_ptr< Material > &material, const std::vector< double > parameters, const std::string quadrature="GAUSS", const unsigned int nGauss=4)
Creates a PML2DQuad4 in a finite element Mesh.
double nx_pml
Horizontal normal components for outward normal.
Definition: PML2DQuad4.hpp:222
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
void UpdateState()
Update 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-...
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
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::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
double L_pml
Thickness of the PML region.
Definition: PML2DQuad4.hpp:210
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
void CommitState()
Save the material states in the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element&#39;s Nodes.
Definition: PML2DQuad4.hpp:231
double x0_pml
Coordinates x0 of reference point.
Definition: PML2DQuad4.hpp:216
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
Eigen::VectorXd ComputePMLStretchingFactors(const double ri, const double si, const double rho, const double mu, const double lambda) const
Evaluates the stretching parameters of PML.
double t
Element thickness.
Definition: PML2DQuad4.hpp:204
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
std::vector< std::unique_ptr< Material > > theMaterial
The Element&#39;s material.
Definition: PML2DQuad4.hpp:234
This file contains the "Element" object declarations, which defines an element in a finite element me...