32 #ifndef _KIN2DQUAD4_HPP_ 33 #define _KIN2DQUAD4_HPP_ 38 #include <Eigen/Dense> 65 kin2DQuad4(
const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material,
const double th,
const std::string quadrature=
"GAUSS",
const unsigned int nGauss=4);
90 void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
95 void SetDamping(
const std::shared_ptr<Damping> &damping);
121 Eigen::MatrixXd
GetStrainAt(
double x3,
double x2)
const;
128 Eigen::MatrixXd
GetStressAt(
double x3,
double x2)
const;
190 Eigen::VectorXd
ComputeBodyForces(
const std::shared_ptr<Load> &body,
unsigned int k=0);
232 Eigen::VectorXd
ComputeStrain(
const double ri,
const double si,
const Eigen::MatrixXd &Jij)
const;
239 Eigen::VectorXd
ComputeStrainRate(
const double ri,
const double si,
const Eigen::MatrixXd &Bij)
const;
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
Eigen::MatrixXd ComputeSecondPiolaKirchhoffMatrix(const Eigen::VectorXd &Stress) const
Compute the Second Piola-Kirchhoff Stress Tensor.
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::MatrixXd ComputeDeformationGradientMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Computes the deformation gradient.
Eigen::MatrixXd TransformVectorToTensor(const Eigen::VectorXd &vector) const
Transform vector components into tensor.
kin2DQuad4(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 kin2DQuad4 in a finite element Mesh.
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the initial stiffness matrix of the element using gauss-integration.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
void UpdateState()
Update the material states in the element.
Eigen::MatrixXd ComputeNonLinearStrainDisplacementMatrix(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 abstract "Material object" declarations, which computes the strain...
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Class for creating a 2D finite kinematic four-node quadrilateral element in a mesh.
Definition: kin2DQuad4.hpp:54
Eigen::VectorXd ComputeStrainRate(const double ri, const double si, const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::MatrixXd ComputeLinearStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Evaluates the strain-displacement matrix at a given Gauss point.
~kin2DQuad4()
Destroys this kin2DQuad4 object.
std::vector< std::unique_ptr< Material > > theMaterial
The Element's material.
Definition: kin2DQuad4.hpp:222
Eigen::VectorXd TransformTensorToVector(const Eigen::MatrixXd &vector) const
Transform tensor components into vector.
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si) const
Evaluates the shape function matrix at a given Gauss point.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
double ComputeEnergy()
Computes the element energy for a given deformation.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si) const
Computes the jacobian of the transformation.
double t
Element thickness.
Definition: kin2DQuad4.hpp:213
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-...
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: kin2DQuad4.hpp:216
void ReverseState()
Reverse the material/section states to previous converged state in this element.
void CommitState()
Save the material states in the element.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
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.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: kin2DQuad4.hpp:219
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
Eigen::MatrixXd ComputeLinearizedStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Evaluates the linearized (small strain) strain-displacement matrix at a given Gauss point...
Eigen::VectorXd ComputeStrain(const double ri, const double si, const Eigen::MatrixXd &Jij) const
Update strain in the element.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: kin2DQuad4.hpp:225
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
This file contains the "Element" object declarations, which defines an element in a finite element me...