Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
EQlin2DQuad4 Class Reference

Class for creating a 2D linearized equivalent linear four-node quadrilateral element in a mesh. More...

#include <EQlin2DQuad4.hpp>

Inheritance diagram for EQlin2DQuad4:
Collaboration diagram for EQlin2DQuad4:

Public Member Functions

 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. More...
 
 ~EQlin2DQuad4 ()
 Destroys this EQlin2DQuad4 object. More...
 
void CommitState ()
 Save the material states in the element. More...
 
void ReverseState ()
 Reverse the material/section states to previous converged state in this element. More...
 
void InitialState ()
 Brings the material/section state to its initial state in this element. More...
 
void UpdateState ()
 Update the material states in the element. More...
 
void SetDomain (std::map< unsigned int, std::shared_ptr< Node > > &nodes)
 Sets the finite element dependance among objects. More...
 
void SetDamping (const std::shared_ptr< Damping > &damping)
 Sets the damping model. More...
 
std::vector< unsigned int > GetTotalDegreeOfFreedom () const
 Gets the list of total-degree of freedom of this Element. More...
 
Eigen::MatrixXd GetStrain () const
 Gets the material/section (generalised) strain. More...
 
Eigen::MatrixXd GetStress () const
 Gets the material/section (generalised) stress. More...
 
Eigen::MatrixXd GetStrainRate () const
 Gets the material/section (generalised) strain-rate. More...
 
Eigen::MatrixXd GetStrainAt (double x3, double x2) const
 Gets the material strain in section at coordinate (x3,x2). More...
 
Eigen::MatrixXd GetStressAt (double x3, double x2) const
 Gets the material stress in section at coordinate (x3,x2). More...
 
Eigen::VectorXd GetVTKResponse (std::string response) const
 Gets the element internal response in VTK format for Paraview display. More...
 
double ComputeEnergy ()
 Computes the element energy for a given deformation. More...
 
Eigen::MatrixXd ComputeMassMatrix ()
 Compute the lumped/consistent mass matrix of the element. More...
 
Eigen::MatrixXd ComputeStiffnessMatrix ()
 Compute the stiffness matrix of the element using gauss-integration. More...
 
Eigen::MatrixXd ComputeDampingMatrix ()
 Compute the damping matrix of the element using gauss-integration. More...
 
Eigen::MatrixXd ComputePMLMatrix ()
 Compute the PML history matrix using gauss-integration. More...
 
Eigen::VectorXd ComputeInternalForces ()
 Compute the internal (elastic) forces acting on the element. More...
 
Eigen::VectorXd ComputeInternalDynamicForces ()
 Compute the elastic, inertial, and viscous forces acting on the element. More...
 
Eigen::VectorXd ComputeSurfaceForces (const std::shared_ptr< Load > &surface, unsigned int face)
 Compute the surface forces acting on the element. More...
 
Eigen::VectorXd ComputeBodyForces (const std::shared_ptr< Load > &body, unsigned int k=0)
 Compute the body forces acting on the element. More...
 
Eigen::VectorXd ComputeDomainReductionForces (const std::shared_ptr< Load > &drm, unsigned int k)
 Compute the domain reduction forces acting on the element. More...
 
- Public Member Functions inherited from Element
 Element (std::string name, const std::vector< unsigned int > nodes, unsigned int ndofs, unsigned int VTKcell, unsigned int SVLcell)
 Creates an Element in a finite element Mesh. More...
 
virtual ~Element ()=0
 Destroys this Element object. More...
 
std::string GetName () const
 Gets the Element Name. More...
 
unsigned int GetVTKCellType () const
 Gets the Element VTK cell type. More...
 
unsigned int GetVTKGroupType () const
 Gets the Element VTK group type. More...
 
unsigned int GetNumberOfNodes () const
 Returns the number of nodes in element. More...
 
unsigned int GetNumberOfDegreeOfFreedom () const
 Returns total number of degree of freedom in the element. More...
 
const std::vector< unsigned int > & GetNodes () const
 Returns the Node Connectivity Indexes. More...
 
bool HasFixedNode (const std::vector< std::shared_ptr< Node > > &nodes) const
 Returns if the element has fixed nodes. More...
 

Private Member Functions

Eigen::VectorXd ComputeStrain (const Eigen::MatrixXd &Bij) const
 Update strain in the element. More...
 
Eigen::VectorXd ComputeStrainRate (const Eigen::MatrixXd &Bij) const
 Update strain rate in the element. More...
 
Eigen::MatrixXd ComputeJacobianMatrix (const double ri, const double si) const
 Computes the jacobian of the transformation. More...
 
Eigen::MatrixXd ComputeShapeFunctionMatrix (const double ri, const double si) const
 Evaluates the shape function matrix at a given Gauss point. More...
 
Eigen::MatrixXd ComputeStrainDisplacementMatrix (const double ri, const double si, const Eigen::MatrixXd &Jij) const
 Evaluates the strain-displacement matrix at a given Gauss point. More...
 
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. More...
 

Private Attributes

double t
 Element thickness. More...
 
double cf1
 Lower corner frequencies for damping. More...
 
double cf2
 Upper corner frequencies for damping. More...
 
double zref
 Reference elevation to compute GGmax and damping. More...
 
std::string Type
 Determine the type of modulus reduction and damping curve. More...
 
std::shared_ptr< DampingtheDamping
 The Damping model. More...
 
std::vector< std::shared_ptr< Node > > theNodes
 The Element's Nodes. More...
 
std::vector< std::unique_ptr< Material > > theMaterial
 The Element's material. More...
 
std::unique_ptr< QuadratureRuleQuadraturePoints
 Coordinate of Gauss points. More...
 

Detailed Description

Class for creating a 2D linearized equivalent linear four-node quadrilateral element in a mesh.

See also
Element.hpp Mesh.hpp

Constructor & Destructor Documentation

◆ EQlin2DQuad4()

EQlin2DQuad4::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.

Parameters
nodesThe Node connectivity array of this Element.
materialPointer to the Material that this Element is made out of.
thThe thickness of the lin2DQuad4 Element.
quadratureThe integration rule to be employed.
nGaussNumber of Gauss points for Element integration.
typeThe type of modulus reduction and damping curve.
zrefReference elevation to compute GGmax and damping.
cf1Lower corner frequencies for damping.
cf2Upper corner frequencies for damping.
Note
More details can be found at EQlin2DQuad4.
See also
EQlin2DQuad4::theNodes, EQlin2DQuad4::theMaterial, EQlin2DQuad4::QuadraturePoints.

◆ ~EQlin2DQuad4()

EQlin2DQuad4::~EQlin2DQuad4 ( )

Destroys this EQlin2DQuad4 object.

Member Function Documentation

◆ CommitState()

void EQlin2DQuad4::CommitState ( )
virtual

Save the material states in the element.

Note
This function sets the trial states as converged ones in Material.

Implements Element.

◆ ComputeBodyForces()

Eigen::VectorXd EQlin2DQuad4::ComputeBodyForces ( const std::shared_ptr< Load > &  body,
unsigned int  k = 0 
)
virtual

Compute the body forces acting on the element.

Parameters
bodyPointer to the Load object that contains this information.
kThe time step at which the body load is evaluated.
Returns
Vector with the Element surface force.
Note
The body force vector can be revisited in ZeroLength1D.
See also
Assembler::ComputeExternalForceVector(), Integrator::ComputeEffectiveForce().

Implements Element.

◆ ComputeDampingMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputeDampingMatrix ( )
virtual

Compute the damping matrix of the element using gauss-integration.

Returns
Matrix with the Element damping matrix.
Note
The damping matrix can be revisited in EQlin2DQuad4.
See also
Assembler::ComputeDampingMatrix(), Integrator::ComputeEffectiveStiffness().

Implements Element.

◆ ComputeDomainReductionForces()

Eigen::VectorXd EQlin2DQuad4::ComputeDomainReductionForces ( const std::shared_ptr< Load > &  drm,
unsigned int  k 
)
virtual

Compute the domain reduction forces acting on the element.

Parameters
drmPointer to the DRM Load object that contains this information.
kThe time step at which the body load is evaluated.
Returns
Vector with the Element domain reduction forces.
Note
The DRM force vector can be revisited in ZeroLength1D.
See also
Assembler::ComputeExternalForceVector(), Integrator::ComputeEffectiveForce().

Implements Element.

◆ ComputeEnergy()

double EQlin2DQuad4::ComputeEnergy ( )
virtual

Computes the element energy for a given deformation.

Returns
Scalar with the element deformation energy.

Implements Element.

◆ ComputeGGmaxDamping()

Eigen::VectorXd EQlin2DQuad4::ComputeGGmaxDamping ( const double  vs,
const double  z,
const double  rho,
const double  gamma 
) const
private

Compute GGmax and Damping for nType = 1 ==> equivalent linear model.

Parameters
vsThe shear-wave velocity.
zThe center of the Element depth.
rhoThe material density.
gammaThe reference strain.
Note
The Gmax and damping are computed as described in EQlin2DQuad4.

◆ ComputeInternalDynamicForces()

Eigen::VectorXd EQlin2DQuad4::ComputeInternalDynamicForces ( )
virtual

Compute the elastic, inertial, and viscous forces acting on the element.

Returns
Vector with the Element dynamic internal force.
Note
The internal force vector can be revisited in Element.
See also
Assembler::ComputeDynamicInternalForceVector().

Implements Element.

◆ ComputeInternalForces()

Eigen::VectorXd EQlin2DQuad4::ComputeInternalForces ( )
virtual

Compute the internal (elastic) forces acting on the element.

Returns
Vector with the Element internal force.
Note
The internal force vector can be revisited in EQlin2DQuad4.
See also
Assembler::ComputeInternalForceVector(), Integrator::ComputeEffectiveForce().

Implements Element.

◆ ComputeJacobianMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputeJacobianMatrix ( const double  ri,
const double  si 
) const
private

Computes the jacobian of the transformation.

Parameters
riThe i-th Gauss coordinate in the r-axis.
siThe i-th Gauss coordinate in the s-axis.
Returns
The Jacobian matrix evaluated at i-th Gauss point (ri,si).

◆ ComputeMassMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputeMassMatrix ( )
virtual

Compute the lumped/consistent mass matrix of the element.

Returns
Matrix with the Element mass matrix.
Note
The mass matrix can be revisited in EQlin2DQuad4.
See also
Assembler::ComputeMassMatrix(), Integrator::ComputeEffectiveStiffness().

Implements Element.

◆ ComputePMLMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputePMLMatrix ( )
virtual

Compute the PML history matrix using gauss-integration.

Returns
Matrix with the PML Element matrix.
Note
The PML matrix is none existent for this element.
See also
Assembler::ComputePMLHistoryMatrix(), Integrator::ComputeEffectiveStiffness().

Implements Element.

◆ ComputeShapeFunctionMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputeShapeFunctionMatrix ( const double  ri,
const double  si 
) const
private

Evaluates the shape function matrix at a given Gauss point.

Parameters
riThe i-th Gauss coordinate in the r-axis.
siThe i-th Gauss coordinate in the s-axis.
Returns
Matrix with the shape function at the Gauss coordinate.

◆ ComputeStiffnessMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputeStiffnessMatrix ( )
virtual

Compute the stiffness matrix of the element using gauss-integration.

Returns
Matrix with the Element stiffness matrix.
Note
The stiffness matrix can be revisited in EQlin2DQuad4.
See also
Assembler::ComputeStiffnessMatrix(), Integrator::ComputeEffectiveStiffness().

Implements Element.

◆ ComputeStrain()

Eigen::VectorXd EQlin2DQuad4::ComputeStrain ( const Eigen::MatrixXd &  Bij) const
private

Update strain in the element.

Parameters
BijThe strain-displacement matrix.
Returns
Vector with the strain at integration point.

◆ ComputeStrainDisplacementMatrix()

Eigen::MatrixXd EQlin2DQuad4::ComputeStrainDisplacementMatrix ( const double  ri,
const double  si,
const Eigen::MatrixXd &  Jij 
) const
private

Evaluates the strain-displacement matrix at a given Gauss point.

Parameters
riThe i-th Gauss coordinate in the r-axis.
siThe i-th Gauss coordinate in the s-axis.
JijThe Jacobian matrix evaluated at (ri,si).
Returns
Matrix with the strain-displacement operator.

◆ ComputeStrainRate()

Eigen::VectorXd EQlin2DQuad4::ComputeStrainRate ( const Eigen::MatrixXd &  Bij) const
private

Update strain rate in the element.

Parameters
BijThe strain-displacement matrix.
Returns
Vector with the strain-rate at integration point.

◆ ComputeSurfaceForces()

Eigen::VectorXd EQlin2DQuad4::ComputeSurfaceForces ( const std::shared_ptr< Load > &  surface,
unsigned int  face 
)
virtual

Compute the surface forces acting on the element.

Parameters
surfacePointer to the Load object that contains this information.
kThe time step at which the surface load is evaluated.
Returns
Vector with the Element surface force.
Note
The surface force vector can be revisited in ZeroLength1D.
See also
Assembler::ComputeExternalForceVector(), Integrator::ComputeEffectiveForce().

Implements Element.

◆ GetStrain()

Eigen::MatrixXd EQlin2DQuad4::GetStrain ( ) const
virtual

Gets the material/section (generalised) strain.

Returns
Matrix with the strain at each integration point.
Note
The index (i,j) are the strain and Gauss-point respectively.

Implements Element.

◆ GetStrainAt()

Eigen::MatrixXd EQlin2DQuad4::GetStrainAt ( double  x3,
double  x2 
) const
virtual

Gets the material strain in section at coordinate (x3,x2).

Parameters
x3Local coordinate along the x3-axis.
x2Local coordinate along the x2-axis.
Returns
Matrix with the strain at coordinate (x3,x2).
Note
The strains are interpolated at this coordinate.

Implements Element.

◆ GetStrainRate()

Eigen::MatrixXd EQlin2DQuad4::GetStrainRate ( ) const
virtual

Gets the material/section (generalised) strain-rate.

Returns
Matrix with the strain-rate at each integration point.
Note
The index (i,j) are the strain-rate and Gauss-point respectively.

Implements Element.

◆ GetStress()

Eigen::MatrixXd EQlin2DQuad4::GetStress ( ) const
virtual

Gets the material/section (generalised) stress.

Returns
Matrix with the stress at each integration point.
Note
The index (i,j) are the stress and Gauss-point respectively.

Implements Element.

◆ GetStressAt()

Eigen::MatrixXd EQlin2DQuad4::GetStressAt ( double  x3,
double  x2 
) const
virtual

Gets the material stress in section at coordinate (x3,x2).

Parameters
x3Local coordinate along the x3-axis.
x2Local coordinate along the x2-axis.
Returns
Matrix with the stresses at coordinate (x3,x2).
Note
The stresses are interpolated at this coordinate.

Implements Element.

◆ GetTotalDegreeOfFreedom()

std::vector<unsigned int> EQlin2DQuad4::GetTotalDegreeOfFreedom ( ) const
virtual

Gets the list of total-degree of freedom of this Element.

Returns
Vector with the list of degree-of-freedom of this Element.

Implements Element.

◆ GetVTKResponse()

Eigen::VectorXd EQlin2DQuad4::GetVTKResponse ( std::string  response) const
virtual

Gets the element internal response in VTK format for Paraview display.

Parameters
responseThe response to be display in Paraview.
Returns
Vector with the response at the Element center.
Note
The current responses are: "Strain", "Stress".

Implements Element.

◆ InitialState()

void EQlin2DQuad4::InitialState ( )
virtual

Brings the material/section state to its initial state in this element.

Note
This function returns the meterial states to the beginning.

Implements Element.

◆ ReverseState()

void EQlin2DQuad4::ReverseState ( )
virtual

Reverse the material/section states to previous converged state in this element.

Note
This function returns the trial states to previous converged states at the Material level.

Implements Element.

◆ SetDamping()

void EQlin2DQuad4::SetDamping ( const std::shared_ptr< Damping > &  damping)
virtual

Sets the damping model.

Parameters
dampingPointer to the damping model.
Note
Several Element objects can share the same damping model.

Implements Element.

◆ SetDomain()

void EQlin2DQuad4::SetDomain ( std::map< unsigned int, std::shared_ptr< Node > > &  nodes)
virtual

Sets the finite element dependance among objects.

Parameters
nodesThe Node list of the Mesh object.
Note
This function sets the relation between Node and Element objects.
See also
lin2DQuad4::theNodes.

Implements Element.

◆ UpdateState()

void EQlin2DQuad4::UpdateState ( )
virtual

Update the material states in the element.

Note
This function update the trial states at the Material level.

Implements Element.

Member Data Documentation

◆ cf1

double EQlin2DQuad4::cf1
private

Lower corner frequencies for damping.

◆ cf2

double EQlin2DQuad4::cf2
private

Upper corner frequencies for damping.

◆ QuadraturePoints

std::unique_ptr<QuadratureRule> EQlin2DQuad4::QuadraturePoints
private

Coordinate of Gauss points.

◆ t

double EQlin2DQuad4::t
private

Element thickness.

◆ theDamping

std::shared_ptr<Damping> EQlin2DQuad4::theDamping
private

The Damping model.

◆ theMaterial

std::vector<std::unique_ptr<Material> > EQlin2DQuad4::theMaterial
private

The Element's material.

◆ theNodes

std::vector<std::shared_ptr<Node> > EQlin2DQuad4::theNodes
private

The Element's Nodes.

◆ Type

std::string EQlin2DQuad4::Type
private

Determine the type of modulus reduction and damping curve.

◆ zref

double EQlin2DQuad4::zref
private

Reference elevation to compute GGmax and damping.