33 #ifndef _HDRBYAMAMOTO3DLINK_HPP_ 34 #define _HDRBYAMAMOTO3DLINK_HPP_ 39 #include <Eigen/Dense> 65 HDRBYamamoto3DLink(
const std::vector<unsigned int> nodes,
double de,
double di,
double hr,
unsigned int dim);
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);
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
double Hr
Rubber bearing height.
Definition: HDRBYamamoto3DLink.hpp:211
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: HDRBYamamoto3DLink.hpp:247
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
Eigen::MatrixXd ComputeLocalAxes() const
Compute local axes.
Eigen::VectorXd Pn
Displacement Trajectory Vector.
Definition: HDRBYamamoto3DLink.hpp:229
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
Eigen::MatrixXd ComputeRotationMatrix() const
Compute rotation matrix.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::VectorXd Fn
Non-linear Restoring Force Vector.
Definition: HDRBYamamoto3DLink.hpp:235
Class for creating a bidirectional 3D two-node link element using Yamamoto HDRB formulation.
Definition: HDRBYamamoto3DLink.hpp:54
Eigen::VectorXd ComputeRelativeDeformation() const
Update strain in the element:
void CommitState()
Save the material states in the element.
~HDRBYamamoto3DLink()
Destroys this HDRBYamamoto3DLink object.
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
HDRBYamamoto3DLink(const std::vector< unsigned int > nodes, double de, double di, double hr, unsigned int dim)
Creates a HDRBYamamoto3DLink in a finite element Mesh.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
double alpha
Yamamoto's model parameter.
Definition: HDRBYamamoto3DLink.hpp:205
Eigen::VectorXd Fc
Non-linear Commited Restoring Force Vector.
Definition: HDRBYamamoto3DLink.hpp:238
Eigen::VectorXd Qaux
Non-linear Displacement Trajectory Vector.
Definition: HDRBYamamoto3DLink.hpp:244
Eigen::VectorXd Qn
Non-linear Displacement Trajectory Vector.
Definition: HDRBYamamoto3DLink.hpp:232
double cr
coefficient of shear stress
Definition: HDRBYamamoto3DLink.hpp:217
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction 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.
double cs
coefficient of shear stress
Definition: HDRBYamamoto3DLink.hpp:220
double n
Yamamoto's model parameter.
Definition: HDRBYamamoto3DLink.hpp:202
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
Eigen::VectorXd Paux
Trial Displacement Trajectory Vector.
Definition: HDRBYamamoto3DLink.hpp:241
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 GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
void InitialState()
Brings the material/section state to its initial state in this element.
double Ar
Rubber bearing area.
Definition: HDRBYamamoto3DLink.hpp:208
void UpdateState()
Update the material states in the element.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
double Krb
Rubber Bearing linear stiffness.
Definition: HDRBYamamoto3DLink.hpp:214
unsigned int Dimension
HDRB link dimension in 3D.
Definition: HDRBYamamoto3DLink.hpp:223
unsigned int nMax
Newton-Raphson maximum number of iterations.
Definition: HDRBYamamoto3DLink.hpp:226
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 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...