32 #ifndef _NULL2DFRAME2_HPP_ 33 #define _NULL2DFRAME2_HPP_ 38 #include <Eigen/Dense> 84 void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
89 void SetDamping(
const std::shared_ptr<Damping> &damping);
115 Eigen::MatrixXd
GetStrainAt(
double x3,
double x2)
const;
122 Eigen::MatrixXd
GetStressAt(
double x3,
double x2)
const;
184 Eigen::VectorXd
ComputeBodyForces(
const std::shared_ptr<Load> &body,
unsigned int k=0);
void InitialState()
Brings the material/section state to its initial state in this 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.
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: null2DFrame2.hpp:199
Class for creating a 2D null frame element (no properties) in a mesh.
Definition: null2DFrame2.hpp:52
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
null2DFrame2(const std::vector< unsigned int > nodes)
Creates a null2DFrame2 in a finite element Mesh.
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: null2DFrame2.hpp:196
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
~null2DFrame2()
Destroys this null2DFrame2 object.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain 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.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
void UpdateState()
Update the material states in the element.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on 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 ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on 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-...
void CommitState()
Save the material states in the element.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
void ReverseState()
Reverse the material/section states to previous converged state in this element.
double ComputeEnergy()
Computes the element energy for a given deformation.
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
This file contains the "Element" object declarations, which defines an element in a finite element me...