34 #ifndef _PML2DQUAD8_HPP_ 35 #define _PML2DQUAD8_HPP_ 40 #include <Eigen/Dense> 67 PML2DQuad8(
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=27);
92 void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
97 void SetDamping(
const std::shared_ptr<Damping> &damping);
123 Eigen::MatrixXd
GetStrainAt(
double x3,
double x2)
const;
130 Eigen::MatrixXd
GetStressAt(
double x3,
double x2)
const;
192 Eigen::VectorXd
ComputeBodyForces(
const std::shared_ptr<Load> &body,
unsigned int k=0);
242 Eigen::VectorXd
ComputeStrain(
const Eigen::MatrixXd &Bij)
const;
275 Eigen::VectorXd
ComputePMLStretchingFactors(
const double ri,
const double si,
const double rho,
const double mu,
const double lambda)
const;
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const
Computes the jacobian of the transformation.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
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.
Class for creating a 2D linearized eight-node perfectly matched layer element in a mesh...
Definition: PML2DQuad8.hpp:56
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
void InitialState()
Brings the material/section state to its initial state in this element.
void UpdateState()
Update the material states in the element.
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
PML2DQuad8(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=27)
Creates a PML2DQuad8 in a finite element Mesh.
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: PML2DQuad8.hpp:237
double t
Element thickness.
Definition: PML2DQuad8.hpp:204
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
double m_pml
Polynomial degree of thePerfectly Match Layer.
Definition: PML2DQuad8.hpp:207
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::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.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
double nx_pml
Horizontal normal components for outward normal.
Definition: PML2DQuad8.hpp:222
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
double ComputeEnergy()
Computes the element energy for a given deformation.
void CommitState()
Save the material states in the element.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of 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.
double x0_pml
Coordinates x0 of reference point.
Definition: PML2DQuad8.hpp:216
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: PML2DQuad8.hpp:231
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).
void ReverseState()
Reverse the material/section states to previous converged state in this element.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
~PML2DQuad8()
Destroys this PML2DQuad8 object.
double L_pml
Thickness of the PML region.
Definition: PML2DQuad8.hpp:210
std::vector< std::unique_ptr< Material > > theMaterial
The Element's material.
Definition: PML2DQuad8.hpp:234
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
double ny_pml
Vertical normal components for outward normal.
Definition: PML2DQuad8.hpp:225
double R_pml
Reflection coefficient for the PML.
Definition: PML2DQuad8.hpp:213
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
double y0_pml
Coordinates y0 of reference point.
Definition: PML2DQuad8.hpp:219
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: PML2DQuad8.hpp:228
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
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...
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
This file contains the "Element" object declarations, which defines an element in a finite element me...