32 #ifndef _LIN2DFRAME2_HPP_ 33 #define _LIN2DFRAME2_HPP_ 38 #include <Eigen/Dense> 65 lin2DFrame2(
const std::vector<unsigned int> nodes, std::unique_ptr<Section> §ion,
bool formulation=
false,
const std::string quadrature=
"GAUSS",
unsigned int nGauss=3);
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);
233 Eigen::VectorXd
ComputeStrain(
const Eigen::MatrixXd &Bij)
const;
void CommitState()
Save the material states in the element.
~lin2DFrame2()
Destroys this lin2DFrame2 object.
double ComputeLength() const
Compute the current length of the element.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri) const
Evaluates the strain-displacement matrix at a given Gauss point.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: lin2DFrame2.hpp:214
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
Eigen::MatrixXd ComputeLocalAxes() const
Compute/update the local axis-1 of the element.
Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const
Update strain in the element.
double L
Length of the element.
Definition: lin2DFrame2.hpp:202
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
This file contains the abstract "Section object" declarations, which computes the strain...
void ReverseState()
Reverse the material/section states to previous converged state in this element.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri) const
Evaluates the shape function matrix at a given Gauss point.
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: lin2DFrame2.hpp:220
std::vector< std::unique_ptr< Section > > theSection
The Element's Section.
Definition: lin2DFrame2.hpp:217
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
Class for creating a 2D linearized two-node frame element in a mesh.
Definition: lin2DFrame2.hpp:54
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: lin2DFrame2.hpp:211
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
lin2DFrame2(const std::vector< unsigned int > nodes, std::unique_ptr< Section > §ion, bool formulation=false, const std::string quadrature="GAUSS", unsigned int nGauss=3)
Creates a lin2DFrame2 in a finite element Mesh.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) 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-...
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the initial stiffness matrix of the element.
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
void InitialState()
Brings the material/section state to its initial state in this element.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
bool Formulation
The element formulation.
Definition: lin2DFrame2.hpp:205
void UpdateState()
Update the material states in the element.
double ComputeEnergy()
Computes the element energy for a given deformation.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress 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.
double Phi
Section Shear Factor.
Definition: lin2DFrame2.hpp:208
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
This file contains the "Element" object declarations, which defines an element in a finite element me...