32 #ifndef _LIN3DTETRA10_HPP_ 33 #define _LIN3DTETRA10_HPP_ 38 #include <Eigen/Dense> 64 lin3DTetra10(
const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material,
const std::string quadrature=
"GAUSS",
const unsigned int nGauss=8);
89 void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
94 void SetDamping(
const std::shared_ptr<Damping> &damping);
120 Eigen::MatrixXd
GetStrainAt(
double x3,
double x2)
const;
127 Eigen::MatrixXd
GetStressAt(
double x3,
double x2)
const;
189 Eigen::VectorXd
ComputeBodyForces(
const std::shared_ptr<Load> &body,
unsigned int k=0);
215 Eigen::VectorXd
ComputeStrain(
const Eigen::MatrixXd &Bij)
const;
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 "Load" class declarations, which defines a load that can act in a node or an e...
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: lin3DTetra10.hpp:210
Eigen::MatrixXd ComputeJacobianMatrix(const double ri, const double si, const double ti) const
Computes the jacobian of the transformation.
void UpdateState()
Update the material states in the element.
double ComputeEnergy()
Computes the element energy for a given deformation.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
void InitialState()
Brings the material/section state to its initial state in this element.
lin3DTetra10(const std::vector< unsigned int > nodes, std::unique_ptr< Material > &material, const std::string quadrature="GAUSS", const unsigned int nGauss=8)
Creates a lin3DTetra10 in a finite element Mesh.
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
std::vector< std::unique_ptr< Material > > theMaterial
The Element's material.
Definition: lin3DTetra10.hpp:207
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping 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).
Class for creating a 3D linearized eight-node hexahedron element in a mesh.
Definition: lin3DTetra10.hpp:54
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
~lin3DTetra10()
Destroys this lin3DTetra10 object.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
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.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
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.
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: lin3DTetra10.hpp:201
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri, const double si, const double ti) const
Evaluates the shape function matrix at a given Gauss point.
void CommitState()
Save the material states in the element.
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 GetStrain() const
Gets the material/section (generalised) strain.
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the stiffness matrix of the element using gauss-integration.
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: lin3DTetra10.hpp:204
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
This file contains the "Element" object declarations, which defines an element in a finite element me...
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.