32 #ifndef _KIN3DTRUSS2_HPP_    33 #define _KIN3DTRUSS2_HPP_    38 #include <Eigen/Dense>    62         kin3DTruss2(
const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material, 
const double area);
    87         void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
    92         void SetDamping(
const std::shared_ptr<Damping> &damping);
   118         Eigen::MatrixXd 
GetStrainAt(
double x3, 
double x2) 
const;
   125         Eigen::MatrixXd 
GetStressAt(
double x3, 
double x2) 
const;
   187         Eigen::VectorXd 
ComputeBodyForces(
const std::shared_ptr<Load> &body, 
unsigned int k=0);
 std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes. 
Definition: kin3DTruss2.hpp:208
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element. 
double ComputeLength() const
Compute the current length of the element. 
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 GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2). 
~kin3DTruss2()
Destroys this kin3DTruss2 object. 
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element. 
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element. 
double A
Area of cross-section. 
Definition: kin3DTruss2.hpp:202
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element. 
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain. 
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::MatrixXd ComputeGeometricAxes() const
Compute/update the geometric global axis of the element. 
Eigen::VectorXd ComputeStrain() const
Update strain in the element: 
std::unique_ptr< Material > theMaterial
The Element's material. 
Definition: kin3DTruss2.hpp:211
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model. 
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element. 
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element. 
std::shared_ptr< Damping > theDamping
The Damping model. 
Definition: kin3DTruss2.hpp:205
double ComputeEnergy()
Computes the element energy for a given deformation. 
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the initial stiffness matrix of the element. 
Class for creating a 3D finite kinematic two-node truss element in a mesh. 
Definition: kin3DTruss2.hpp:53
Eigen::VectorXd ComputeStrainRate() const
Update strain rate in the element: 
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element. 
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate. 
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2). 
Eigen::MatrixXd ComputeTransformationAxes() const
Compute/update the tranformation matrix from local to global axis. 
Virtual class for creating an element in a mesh. 
Definition: Element.hpp:51
kin3DTruss2(const std::vector< unsigned int > nodes, std::unique_ptr< Material > &material, const double area)
Creates a kin3DTruss2 in a finite element Mesh. 
void CommitState()
Save the material states in the element. 
void InitialState()
Brings the material/section state to its initial state in this 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. 
Eigen::VectorXd ComputeLocalAxes() const
Compute/update the local axis-1 of the element. 
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress. 
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration. 
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping 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. 
double ComputeAxialForce() const
Compute the axial force in element. 
double Lo
Length of the element undeformed configuration. 
Definition: kin3DTruss2.hpp:199
void UpdateState()
Update the material states in the element. 
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration. 
void ReverseState()
Reverse the material/section states to previous converged state in this element. 
This file contains the "Element" object declarations, which defines an element in a finite element me...