Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
PML3DHexa8.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] Fathi, A., Poursartip, B., & Kallivokas, L. F. (2015). Time‐domain hybrid
25 // formulations for wave simulations in three‐dimensional PML‐truncated
26 // heterogeneous media. International Journal for Numerical Methods in
27 // Engineering, 101(3), 165-198.
28 //
29 // Description:
33 //------------------------------------------------------------------------------
34 
35 #ifndef _PML3DHEXA8_HPP_
36 #define _PML3DHEXA8_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 "QuadratureRule.hpp"
48 #include "Damping.hpp"
49 
57 class PML3DHexa8 : public Element{
58 
59  public:
68  PML3DHexa8(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=8);
69 
71  ~PML3DHexa8();
72 
75  void CommitState();
76 
79  void ReverseState();
80 
83  void InitialState();
84 
87  void UpdateState();
88 
93  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
94 
98  void SetDamping(const std::shared_ptr<Damping> &damping);
99 
102  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
103 
107  Eigen::MatrixXd GetStrain() const;
108 
112  Eigen::MatrixXd GetStress() const;
113 
117  Eigen::MatrixXd GetStrainRate() const;
118 
124  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
125 
131  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
132 
137  Eigen::VectorXd GetVTKResponse(std::string response) const;
138 
141  double ComputeEnergy();
142 
147  Eigen::MatrixXd ComputeMassMatrix();
148 
153  Eigen::MatrixXd ComputeStiffnessMatrix();
154 
159  Eigen::MatrixXd ComputeDampingMatrix();
160 
165  Eigen::MatrixXd ComputePMLMatrix();
166 
171  Eigen::VectorXd ComputeInternalForces();
172 
177  Eigen::VectorXd ComputeInternalDynamicForces();
178 
185  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
186 
193  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
194 
201  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
202 
203  private:
205  double m_pml;
206 
208  double L_pml;
209 
211  double R_pml;
212 
214  double x0_pml;
215 
217  double y0_pml;
218 
220  double z0_pml;
221 
223  double nx_pml;
224 
226  double ny_pml;
227 
229  double nz_pml;
230 
232  std::shared_ptr<Damping> theDamping;
233 
235  std::vector<std::shared_ptr<Node> > theNodes;
236 
238  std::vector<std::unique_ptr<Material> > theMaterial;
239 
241  std::unique_ptr<QuadratureRule> QuadraturePoints;
242 
246  Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const;
247 
251  Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const;
252 
258  Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si, const double ti) const;
259 
265  Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si, const double ti) const;
266 
273  Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const double ti, const Eigen::MatrixXd &Jij) const;
274 
283  Eigen::VectorXd ComputePMLStretchingFactors(const double ri, const double si, const double ti, const double rho, const double mu, const double lambda) const;
284 };
285 
286 #endif
double nx_pml
Horizontal normal components for x1.
Definition: PML3DHexa8.hpp:223
double y0_pml
Coordinates y0 of reference point.
Definition: PML3DHexa8.hpp:217
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
double z0_pml
Coordinates z0 of reference point.
Definition: PML3DHexa8.hpp:220
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
This file contains the abstract "Material object" declarations, which computes the strain...
void ReverseState()
Reverse the material/section states to previous converged state in this element.
double R_pml
Reflection coefficient for the PML.
Definition: PML3DHexa8.hpp:211
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
double x0_pml
Coordinates x0 of reference point.
Definition: PML3DHexa8.hpp:214
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: PML3DHexa8.hpp:232
Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri, const double si, const double ti, const Eigen::MatrixXd &Jij) const
Evaluates the strain-displacement matrix at a given Gauss point.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
std::vector< std::unique_ptr< Material > > theMaterial
The Element&#39;s material.
Definition: PML3DHexa8.hpp:238
Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const
Update strain in the element.
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si, const double ti) const
Computes the jacobian of the transformation.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
void UpdateState()
Update the material states in the element.
void InitialState()
Brings the material/section state to its initial state in this element.
double m_pml
Polynomial degree of thePerfectly Match Layer.
Definition: PML3DHexa8.hpp:205
Eigen::VectorXd ComputePMLStretchingFactors(const double ri, const double si, const double ti, const double rho, const double mu, const double lambda) const
Evaluates the stretching parameters of PML.
double ny_pml
Horizontal normal components for x2.
Definition: PML3DHexa8.hpp:226
void CommitState()
Save the material states in the element.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
std::vector< std::shared_ptr< Node > > theNodes
The Element&#39;s Nodes.
Definition: PML3DHexa8.hpp:235
double ComputeEnergy()
Computes the element energy for a given deformation.
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-...
double nz_pml
Vertical normal components for x3.
Definition: PML3DHexa8.hpp:229
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si, const double ti) const
Evaluates the shape function matrix at a given Gauss point.
PML3DHexa8(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=8)
Creates a PML3DHexa8 in a finite element Mesh.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
double L_pml
Thickness of the PML region.
Definition: PML3DHexa8.hpp:208
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: PML3DHexa8.hpp:241
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
~PML3DHexa8()
Destroys this PML3DHexa8 object.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
Class for creating a 3D linearized eight-node perfectly matched layer element in a mesh...
Definition: PML3DHexa8.hpp:57
This file contains the "Element" object declarations, which defines an element in a finite element me...