34 #ifndef _UNXBOUCWEN2DLINK_HPP_ 35 #define _UNXBOUCWEN2DLINK_HPP_ 40 #include <Eigen/Dense> 68 UnxBoucWen2DLink(
const std::vector<unsigned int> nodes, std::vector<double> parameters, std::vector<double> variables,
const unsigned int dim,
const unsigned int dir,
double tol=1E-06,
unsigned int nmax=50);
93 void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
98 void SetDamping(
const std::shared_ptr<Damping> &damping);
124 Eigen::MatrixXd
GetStrainAt(
double x3,
double x2)
const;
131 Eigen::MatrixXd
GetStressAt(
double x3,
double x2)
const;
193 Eigen::VectorXd
ComputeBodyForces(
const std::shared_ptr<Load> &body,
unsigned int k=0);
262 double sign(
double x)
const;
double eta
Yielding exponent (sharpness of hysteresis loop corners)
Definition: UnxBoucWen2DLink.hpp:214
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
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.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
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.
double alpha
Yielding exponent (sharpness of hysteresis loop corners)
Definition: UnxBoucWen2DLink.hpp:211
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
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.
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
This file contains the abstract "Material object" declarations, which computes the strain...
Class for creating an uniaxial 2D two-node link element using Bouc-Wen formulation.
Definition: UnxBoucWen2DLink.hpp:55
double beta
First hysteretic shape parameter.
Definition: UnxBoucWen2DLink.hpp:217
double gamma
second hysteretic shape parameter
Definition: UnxBoucWen2DLink.hpp:220
Eigen::MatrixXd ComputeLocalAxes() const
Compute local axes.
double qbw
Non-linear Bouc-Wen internal force.
Definition: UnxBoucWen2DLink.hpp:238
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: UnxBoucWen2DLink.hpp:259
double kbw
Consistent Bouc-Wen stiffness.
Definition: UnxBoucWen2DLink.hpp:244
double Un
Displacements in local coordinates.
Definition: UnxBoucWen2DLink.hpp:235
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
Eigen::MatrixXd ComputeRotationMatrix() const
Compute rotation matrix.
void CommitState()
Save the material states in the element.
double sign(double x) const
Sign function.
Eigen::VectorXd ComputeRelativeDeformation() const
Update strain in the element:
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
double U
Trial displacements in local coordinates.
Definition: UnxBoucWen2DLink.hpp:232
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
double zn
Hysteretic evolution parameters.
Definition: UnxBoucWen2DLink.hpp:229
unsigned int Direction
The element direction in global coordinates.
Definition: UnxBoucWen2DLink.hpp:256
double qbc
Non-linear commited Bouc-Wen internal force.
Definition: UnxBoucWen2DLink.hpp:241
unsigned int nMax
Newton-Raphson maximum number of iterations.
Definition: UnxBoucWen2DLink.hpp:250
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
void UpdateState()
Update the material states in the element.
unsigned int Dimension
The element dimension in global coordinates.
Definition: UnxBoucWen2DLink.hpp:253
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
double Tol
Newton-Raphson Tolerance.
Definition: UnxBoucWen2DLink.hpp:223
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.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
double ComputeEnergy()
Computes the element energy for a given deformation.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
double Ko
Initial stiffness of hysteretic component.
Definition: UnxBoucWen2DLink.hpp:205
~UnxBoucWen2DLink()
Destroys this UnxBoucWen2DLink object.
double kbc
Consistent commited Bouc-Wen stiffness.
Definition: UnxBoucWen2DLink.hpp:247
double Fy
Yield force.
Definition: UnxBoucWen2DLink.hpp:208
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
double z
Trial Hysteretic evolution parameters.
Definition: UnxBoucWen2DLink.hpp:226
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
UnxBoucWen2DLink(const std::vector< unsigned int > nodes, std::vector< double > parameters, std::vector< double > variables, const unsigned int dim, const unsigned int dir, double tol=1E-06, unsigned int nmax=50)
Creates a UnxBoucWen2DLink in a finite element Mesh.
This file contains the "Element" object declarations, which defines an element in a finite element me...