16 #ifndef SURGSIM_PHYSICS_FEMELEMENT_H
17 #define SURGSIM_PHYSICS_FEMELEMENT_H
65 size_t getNodeId(
size_t elementNodeId)
const;
129 double scale = 1.0) = 0;
140 double scale = 1.0) = 0;
224 #endif // SURGSIM_PHYSICS_FEMELEMENT_H
Definition: DriveElementFromInputBehavior.cpp:27
double getPoissonRatio() const
Gets the Poisson ratio (unitless)
Definition: FemElement.cpp:83
size_t m_numDofPerNode
Number of degree of freedom per node for this element.
Definition: FemElement.h:205
double getMassDensity() const
Gets the mass density (in Kg.m-3)
Definition: FemElement.cpp:93
const std::vector< size_t > & getNodeIds() const
Gets the node ids for this element.
Definition: FemElement.cpp:63
virtual SurgSim::Math::Vector computeNaturalCoordinate(const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &cartesianCoordinate) const =0
Computes a natural coordinate given a global coordinate.
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > Matrix
A dynamic size matrix.
Definition: Matrix.h:65
double m_nu
Poisson ratio (unitless)
Definition: FemElement.h:217
virtual void addMatVec(const SurgSim::Math::OdeState &state, double alphaM, double alphaD, double alphaK, const SurgSim::Math::Vector &x, SurgSim::Math::Vector *F)=0
Adds the element matrix-vector contribution F += (alphaM.M + alphaD.D + alphaK.K).x (computed for a given state) into a complete system data structure F (assembly)
void setYoungModulus(double E)
Sets the Young modulus (in N.m-2)
Definition: FemElement.cpp:68
virtual void addMass(const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *M, double scale=1.0)=0
Adds the element mass matrix M (computed for a given state) to a complete system mass matrix M (assem...
OdeState defines the state y of an ode of 2nd order of the form M(x,v).a = F(x, v) with boundary cond...
Definition: OdeState.h:34
Base class for all Fem Element (1D, 2D, 3D) It handles the node ids to which it is connected and requ...
Definition: FemElement.h:42
double m_rho
Mass density (in Kg.m-3)
Definition: FemElement.h:211
virtual void addForce(const SurgSim::Math::OdeState &state, SurgSim::Math::Vector *F, double scale=1.0)=0
Adds the element force (computed for a given state) to a complete system force vector F (assembly) ...
size_t getNumDofPerNode() const
Gets the number of degree of freedom per node.
Definition: FemElement.cpp:43
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
A dynamic size column vector.
Definition: Vector.h:67
virtual void initialize(const SurgSim::Math::OdeState &state)
Initialize the FemElement once everything has been set.
Definition: FemElement.cpp:33
double m_E
Young modulus (in N.m-2)
Definition: FemElement.h:214
bool isValidCoordinate(const SurgSim::Math::Vector &naturalCoordinate) const
Determines whether a given natural coordinate is valid.
Definition: FemElement.cpp:108
double getYoungModulus() const
Gets the Young modulus (in N.m-2)
Definition: FemElement.cpp:73
virtual SurgSim::Math::Vector computeCartesianCoordinate(const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &naturalCoordinate) const =0
Computes a given natural coordinate in cartesian coordinates.
virtual void addStiffness(const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *K, double scale=1.0)=0
Adds the element stiffness matrix K (= -df/dx) (computed for a given state) to a complete system stif...
void setMassDensity(double rho)
Sets the mass density (in Kg.m-3)
Definition: FemElement.cpp:88
Definitions of small fixed-size square matrix types.
double getMass(const SurgSim::Math::OdeState &state) const
Gets the element mass based on the input state (in Kg)
Definition: FemElement.cpp:98
Definitions of small fixed-size vector types.
std::vector< size_t > m_nodeIds
Node ids connected by this element.
Definition: FemElement.h:208
virtual double getVolume(const SurgSim::Math::OdeState &state) const =0
Gets the element volume based on the input state (in m-3)
virtual bool update(const SurgSim::Math::OdeState &state)
Update the element based on a given state.
Definition: FemElement.cpp:103
size_t getNumNodes() const
Gets the number of nodes connected by this element.
Definition: FemElement.cpp:53
virtual void addFMDK(const SurgSim::Math::OdeState &state, SurgSim::Math::Vector *F, SurgSim::Math::Matrix *M, SurgSim::Math::Matrix *D, SurgSim::Math::Matrix *K)=0
Adds the element force vector, mass, stiffness and damping matrices (computed for a given state) into...
void setNumDofPerNode(size_t numDofPerNode)
Sets the number of degrees of freedom per node.
Definition: FemElement.cpp:48
virtual void addDamping(const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *D, double scale=1.0)=0
Adds the element damping matrix D (= -df/dv) (comuted for a given state) to a complete system damping...
size_t getNodeId(size_t elementNodeId) const
Gets the elementNodeId-th node id.
Definition: FemElement.cpp:58
void setPoissonRatio(double nu)
Sets the Poisson ratio (unitless)
Definition: FemElement.cpp:78
virtual ~FemElement()
Virtual destructor.
Definition: FemElement.cpp:30
FemElement()
Constructor.
Definition: FemElement.cpp:27