Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
SurgSim::Physics::Fem1DElementBeam Class Reference

1D FemElement based on a beam volume discretization with a fixed cross section More...

#include <SurgSim/Physics/Fem1DElementBeam.h>

Inheritance diagram for SurgSim::Physics::Fem1DElementBeam:
SurgSim::Physics::FemElement

Public Member Functions

 Fem1DElementBeam (std::array< size_t, 2 > nodeIds)
 Constructor. More...
 
void setRadius (double radius)
 Sets the beam's circular cross-section radius. More...
 
double getRadius () const
 Gets the beam's circular cross-section radius. More...
 
virtual void initialize (const SurgSim::Math::OdeState &state) override
 Initializes the FemElement once everything has been set. More...
 
virtual double getVolume (const SurgSim::Math::OdeState &state) const override
 Gets the element's volume based on the input state. More...
 
bool getShearingEnabled () const
 Gets whether shearing is enabled for the element. More...
 
void setShearingEnabled (bool enabled)
 Enables or disables shearing for the element. More...
 
virtual void addForce (const SurgSim::Math::OdeState &state, SurgSim::Math::Vector *F, double scale=1.0) override
 Adds the element's force (computed for a given state) to a complete system force vector F (assembly) More...
 
virtual void addMass (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *M, double scale=1.0) override
 Adds the element's mass matrix M (computed for a given state) to a complete system mass matrix M (assembly) More...
 
virtual void addDamping (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *D, double scale=1.0) override
 Adds the element's damping matrix D (= -df/dv) (computed for a given state) to a complete system damping matrix D (assembly) More...
 
virtual void addStiffness (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *K, double scale=1.0) override
 Adds the element's stiffness matrix K (= -df/dx) (computed for a given state) to a complete system stiffness matrix K (assembly) More...
 
virtual void addFMDK (const SurgSim::Math::OdeState &state, SurgSim::Math::Vector *F, SurgSim::Math::Matrix *M, SurgSim::Math::Matrix *D, SurgSim::Math::Matrix *K) override
 Adds the element's force vector, mass, stiffness and damping matrices (computed for a given state) into a complete system data structure F, M, D, K (assembly) More...
 
virtual void addMatVec (const SurgSim::Math::OdeState &state, double alphaM, double alphaD, double alphaK, const SurgSim::Math::Vector &x, SurgSim::Math::Vector *F)
 Adds the element's matrix-vector contribution F += (alphaM.M + alphaD.D + alphaK.K).x (computed for a given state) into a complete system data structure F (assembly) More...
 
virtual SurgSim::Math::Vector computeCartesianCoordinate (const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &naturalCoordinate) const
 Computes a given natural coordinate in cartesian coordinates. More...
 
virtual SurgSim::Math::Vector computeNaturalCoordinate (const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &cartesianCoordinate) const override
 Computes a natural coordinate given a global coordinate. More...
 
- Public Member Functions inherited from SurgSim::Physics::FemElement
 FemElement ()
 Constructor. More...
 
virtual ~FemElement ()
 Virtual destructor. More...
 
size_t getNumDofPerNode () const
 Gets the number of degree of freedom per node. More...
 
size_t getNumNodes () const
 Gets the number of nodes connected by this element. More...
 
size_t getNodeId (size_t elementNodeId) const
 Gets the elementNodeId-th node id. More...
 
const std::vector< size_t > & getNodeIds () const
 Gets the node ids for this element. More...
 
void setYoungModulus (double E)
 Sets the Young modulus (in N.m-2) More...
 
double getYoungModulus () const
 Gets the Young modulus (in N.m-2) More...
 
void setPoissonRatio (double nu)
 Sets the Poisson ratio (unitless) More...
 
double getPoissonRatio () const
 Gets the Poisson ratio (unitless) More...
 
void setMassDensity (double rho)
 Sets the mass density (in Kg.m-3) More...
 
double getMassDensity () const
 Gets the mass density (in Kg.m-3) More...
 
double getMass (const SurgSim::Math::OdeState &state) const
 Gets the element mass based on the input state (in Kg) More...
 
virtual bool update (const SurgSim::Math::OdeState &state)
 Update the element based on a given state. More...
 
bool isValidCoordinate (const SurgSim::Math::Vector &naturalCoordinate) const
 Determines whether a given natural coordinate is valid. More...
 

Protected Member Functions

void computeInitialRotation (const SurgSim::Math::OdeState &state)
 Computes the beam element's initial rotation. More...
 
void computeStiffness (const SurgSim::Math::OdeState &state, Eigen::Matrix< double, 12, 12 > *k)
 Computes the beam's stiffness matrix. More...
 
void computeMass (const SurgSim::Math::OdeState &state, Eigen::Matrix< double, 12, 12 > *m)
 Computes the beam's mass matrix. More...
 
- Protected Member Functions inherited from SurgSim::Physics::FemElement
void setNumDofPerNode (size_t numDofPerNode)
 Sets the number of degrees of freedom per node. More...
 

Protected Attributes

Eigen::Matrix< double, 12, 1 > m_x0
 The element's rest state. More...
 
Eigen::Matrix< double, 12, 12 > m_R0
 Initial rotation matrix for the element. More...
 
Eigen::Matrix< double, 12, 12 > m_M
 Mass matrix (in global coordinate frame) More...
 
Eigen::Matrix< double, 12, 12 > m_MLocal
 Stiffness matrix (in local coordinate frame) More...
 
Eigen::Matrix< double, 12, 12 > m_K
 Stiffness matrix (in global coordinate frame) More...
 
Eigen::Matrix< double, 12, 12 > m_KLocal
 Stiffness matrix (in local coordinate frame) More...
 
double m_G
 Physical shear modulus G = E/( 2(1+mu) ) More...
 
double m_restLength
 Rest length. More...
 
double m_radius
 radius for a circular Beam More...
 
double m_A
 Cross sectional area = PI.radius.radius if circular. More...
 
bool m_haveShear
 Does this beam element have shear. More...
 
double m_shearFactor
 Shear factor (usually 5/8) More...
 
double m_Asy
 The shear area in the y and z directions (=0 => no shear) http://en.wikipedia.org/wiki/Timoshenko_beam_theory. More...
 
double m_Asz
 
double m_Phi_y
 Shear deformation parameters Phi_y=12.E.Iz/(G.Asy.L^2) or 0 if As?=0 Phi_z=12.E.Iy/(G.Asz.L^2) or 0 if As?=0. More...
 
double m_Phi_z
 
double m_Iy
 Cross sectional moment of inertia. More...
 
double m_Iz
 
double m_J
 Polar moment of inertia. More...
 
- Protected Attributes inherited from SurgSim::Physics::FemElement
size_t m_numDofPerNode
 Number of degree of freedom per node for this element. More...
 
std::vector< size_t > m_nodeIds
 Node ids connected by this element. More...
 
double m_rho
 Mass density (in Kg.m-3) More...
 
double m_E
 Young modulus (in N.m-2) More...
 
double m_nu
 Poisson ratio (unitless) More...
 

Detailed Description

1D FemElement based on a beam volume discretization with a fixed cross section

The inertia property (mass) and the stiffness matrices are derived from "Theory of Matrix Structural Analysis" from J.S. Przemieniecki. The deformation is based on linear elasticity theory and not on visco-elasticity theory; therefore, the element does not have any damping components.

Note
The element is considered to have a circular cross section.

Constructor & Destructor Documentation

SurgSim::Physics::Fem1DElementBeam::Fem1DElementBeam ( std::array< size_t, 2 >  nodeIds)
explicit

Constructor.

Parameters
nodeIdsAn array of 2 node ids (A, B) defining this beam element with respect to a DeformableRepresentaitonState which is passed to the initialize method.

Member Function Documentation

void SurgSim::Physics::Fem1DElementBeam::addDamping ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix D,
double  scale = 1.0 
)
overridevirtual

Adds the element's damping matrix D (= -df/dv) (computed for a given state) to a complete system damping matrix D (assembly)

Parameters
stateThe state to compute the damping matrix with
[in,out]DThe complete system damping matrix to add the element damping matrix into
scaleA factor to scale the added damping matrix with
Note
The element's damping matrix is a square matrix of size getNumDofPerNode() x getNumNodes().
This method supposes that the incoming state contains information with the same number of dof per node as getNumDofPerNode().
The beam uses linear elasticity (not visco-elasticity), so it does not have any damping.

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::addFMDK ( const SurgSim::Math::OdeState state,
SurgSim::Math::Vector F,
SurgSim::Math::Matrix M,
SurgSim::Math::Matrix D,
SurgSim::Math::Matrix K 
)
overridevirtual

Adds the element's force vector, mass, stiffness and damping matrices (computed for a given state) into a complete system data structure F, M, D, K (assembly)

Parameters
stateThe state to compute everything with
[in,out]FThe complete system force vector to add the element force into
[in,out]MThe complete system mass matrix to add the element mass matrix into
[in,out]DThe complete system damping matrix to add the element damping matrix into
[in,out]KThe complete system stiffness matrix to add the element stiffness matrix into
Note
This method supposes that the incoming state contains information with the same number of dof per node as getNumDofPerNode().

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::addForce ( const SurgSim::Math::OdeState state,
SurgSim::Math::Vector F,
double  scale = 1.0 
)
overridevirtual

Adds the element's force (computed for a given state) to a complete system force vector F (assembly)

Parameters
stateThe state to compute the force with
[in,out]FThe complete system force vector to add the element's force into
scaleA factor to scale the added force with
Note
The element's force is of size (getNumDofPerNode() x getNumNodes()).
This method supposes that the incoming state contains information with the same number of dof per node as getNumDofPerNode().

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::addMass ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix M,
double  scale = 1.0 
)
overridevirtual

Adds the element's mass matrix M (computed for a given state) to a complete system mass matrix M (assembly)

Parameters
stateThe state to compute the mass matrix with
[in,out]MThe complete system mass matrix to add the element's mass-matrix into
scaleA factor to scale the added mass matrix with
Note
The element's mass matrix is a square matrix of size getNumDofPerNode() x getNumNodes().
This method supposes that the incoming state contains information with the same number of dof per node as getNumDofPerNode()

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::addMatVec ( const SurgSim::Math::OdeState state,
double  alphaM,
double  alphaD,
double  alphaK,
const SurgSim::Math::Vector x,
SurgSim::Math::Vector F 
)
virtual

Adds the element's matrix-vector contribution F += (alphaM.M + alphaD.D + alphaK.K).x (computed for a given state) into a complete system data structure F (assembly)

Parameters
stateThe state to compute everything with
alphaMThe scaling factor for the mass contribution
alphaDThe scaling factor for the damping contribution
alphaKThe scaling factor for the stiffness contribution
xA complete system vector to use as the vector in the matrix-vector multiplication
[in,out]FThe complete system force vector to add the element matrix-vector contribution into
Note
This method supposes that the incoming state contains information with the same number of dof per node as getNumDofPerNode().

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::addStiffness ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix K,
double  scale = 1.0 
)
overridevirtual

Adds the element's stiffness matrix K (= -df/dx) (computed for a given state) to a complete system stiffness matrix K (assembly)

Parameters
stateThe state to compute the stiffness matrix with
[in,out]KThe complete system stiffness matrix to add the element stiffness matrix into
scaleA factor to scale the added stiffness matrix with
Note
The element stiffness matrix is square of size getNumDofPerNode() x getNumNodes().
This method supposes that the incoming state contains information with the same number of dof per node as getNumDofPerNode()

Implements SurgSim::Physics::FemElement.

SurgSim::Math::Vector SurgSim::Physics::Fem1DElementBeam::computeCartesianCoordinate ( const SurgSim::Math::OdeState state,
const SurgSim::Math::Vector naturalCoordinate 
) const
virtual

Computes a given natural coordinate in cartesian coordinates.

Parameters
stateThe state at which to transform coordinates
naturalCoordinateThe coordinates to transform
Returns
The resultant cartesian coordinates

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::computeInitialRotation ( const SurgSim::Math::OdeState state)
protected

Computes the beam element's initial rotation.

Parameters
stateThe state to compute the rotation from
Note
This method stores the result in m_R0
void SurgSim::Physics::Fem1DElementBeam::computeMass ( const SurgSim::Math::OdeState state,
Eigen::Matrix< double, 12, 12 > *  m 
)
protected

Computes the beam's mass matrix.

Parameters
stateThe state to compute the stiffness matrix from
[out]mThe mass matrix to store the result into
SurgSim::Math::Vector SurgSim::Physics::Fem1DElementBeam::computeNaturalCoordinate ( const SurgSim::Math::OdeState state,
const SurgSim::Math::Vector cartesianCoordinate 
) const
overridevirtual

Computes a natural coordinate given a global coordinate.

Parameters
stateThe state at which to transform coordinates
cartesianCoordinateThe coordinates to transform
Returns
The resultant natural coordinates

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::computeStiffness ( const SurgSim::Math::OdeState state,
Eigen::Matrix< double, 12, 12 > *  k 
)
protected

Computes the beam's stiffness matrix.

Parameters
stateThe state to compute the stiffness matrix from
[out]kThe stiffness matrix to store the result into
double SurgSim::Physics::Fem1DElementBeam::getRadius ( ) const

Gets the beam's circular cross-section radius.

Returns
The radius of the beam
bool SurgSim::Physics::Fem1DElementBeam::getShearingEnabled ( ) const

Gets whether shearing is enabled for the element.

Returns
True if shearing is enabled
double SurgSim::Physics::Fem1DElementBeam::getVolume ( const SurgSim::Math::OdeState state) const
overridevirtual

Gets the element's volume based on the input state.

Parameters
stateThe state to compute the volume with
Returns
The element's volume

Implements SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::initialize ( const SurgSim::Math::OdeState state)
overridevirtual

Initializes the FemElement once everything has been set.

Parameters
stateThe state to initialize the FemElement with
Note
We use the theory of linear elasticity, so this method pre-computes the stiffness and mass matrices

Reimplemented from SurgSim::Physics::FemElement.

void SurgSim::Physics::Fem1DElementBeam::setRadius ( double  radius)

Sets the beam's circular cross-section radius.

Parameters
radiusThe radius of the beam
void SurgSim::Physics::Fem1DElementBeam::setShearingEnabled ( bool  enabled)

Enables or disables shearing for the element.

Shearing can only be meaningfully enabled or disabled before the element has had initialize called.

Parameters
enabledBoolean determining whether shearing is enabled

Member Data Documentation

double SurgSim::Physics::Fem1DElementBeam::m_A
protected

Cross sectional area = PI.radius.radius if circular.

double SurgSim::Physics::Fem1DElementBeam::m_Asy
protected

The shear area in the y and z directions (=0 => no shear) http://en.wikipedia.org/wiki/Timoshenko_beam_theory.

double SurgSim::Physics::Fem1DElementBeam::m_Asz
protected
double SurgSim::Physics::Fem1DElementBeam::m_G
protected

Physical shear modulus G = E/( 2(1+mu) )

bool SurgSim::Physics::Fem1DElementBeam::m_haveShear
protected

Does this beam element have shear.

double SurgSim::Physics::Fem1DElementBeam::m_Iy
protected

Cross sectional moment of inertia.

double SurgSim::Physics::Fem1DElementBeam::m_Iz
protected
double SurgSim::Physics::Fem1DElementBeam::m_J
protected

Polar moment of inertia.

Eigen::Matrix<double, 12, 12> SurgSim::Physics::Fem1DElementBeam::m_K
protected

Stiffness matrix (in global coordinate frame)

Eigen::Matrix<double, 12, 12> SurgSim::Physics::Fem1DElementBeam::m_KLocal
protected

Stiffness matrix (in local coordinate frame)

Eigen::Matrix<double, 12, 12> SurgSim::Physics::Fem1DElementBeam::m_M
protected

Mass matrix (in global coordinate frame)

Eigen::Matrix<double, 12, 12> SurgSim::Physics::Fem1DElementBeam::m_MLocal
protected

Stiffness matrix (in local coordinate frame)

double SurgSim::Physics::Fem1DElementBeam::m_Phi_y
protected

Shear deformation parameters Phi_y=12.E.Iz/(G.Asy.L^2) or 0 if As?=0 Phi_z=12.E.Iy/(G.Asz.L^2) or 0 if As?=0.

double SurgSim::Physics::Fem1DElementBeam::m_Phi_z
protected
Eigen::Matrix<double, 12, 12> SurgSim::Physics::Fem1DElementBeam::m_R0
protected

Initial rotation matrix for the element.

double SurgSim::Physics::Fem1DElementBeam::m_radius
protected

radius for a circular Beam

double SurgSim::Physics::Fem1DElementBeam::m_restLength
protected

Rest length.

double SurgSim::Physics::Fem1DElementBeam::m_shearFactor
protected

Shear factor (usually 5/8)

Eigen::Matrix<double, 12, 1> SurgSim::Physics::Fem1DElementBeam::m_x0
protected

The element's rest state.


The documentation for this class was generated from the following files: