opensurgsim
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 <Fem1DElementBeam.h>

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

Public Member Functions

 Fem1DElementBeam ()
 Constructor.
 
 Fem1DElementBeam (std::array< size_t, 2 > nodeIds)
 Constructor. More...
 
 Fem1DElementBeam (std::shared_ptr< FemElementStructs::FemElementParameter > elementData)
 Constructor for FemElement object factory. 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...
 
void initialize (const SurgSim::Math::OdeState &state) override
 Initialize the FemElement once everything has been set. More...
 
double getVolume (const SurgSim::Math::OdeState &state) const override
 Gets the element volume based on the input state (in m-3) 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...
 
SurgSim::Math::Vector computeCartesianCoordinate (const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &naturalCoordinate) const override
 Computes a given natural coordinate in cartesian coordinates. More...
 
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...
 
const SurgSim::Math::Matrix33dgetInitialRotation () const
 
- Public Member Functions inherited from SurgSim::Physics::FemElement
 FemElement ()
 Constructor.
 
virtual ~FemElement ()
 Virtual destructor.
 
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 void addForce (SurgSim::Math::Vector *F, double scale) const
 Adds the element force (computed for a given state) to a complete system force vector F (assembly) More...
 
virtual void addForce (SurgSim::Math::Vector *F) const
 
virtual void addMass (SurgSim::Math::SparseMatrix *M, double scale) const
 Adds the element mass matrix M (computed for a given state) to a complete system mass matrix M (assembly) More...
 
virtual void addMass (SurgSim::Math::SparseMatrix *M) const
 
virtual void addDamping (SurgSim::Math::SparseMatrix *D, double scale) const
 Adds the element damping matrix D (= -df/dv) (comuted for a given state) to a complete system damping matrix D (assembly) More...
 
virtual void addDamping (SurgSim::Math::SparseMatrix *D) const
 
virtual void addStiffness (SurgSim::Math::SparseMatrix *K, double scale) const
 Adds the element stiffness matrix K (= -df/dx) (computed for a given state) to a complete system stiffness matrix K (assembly) More...
 
virtual void addStiffness (SurgSim::Math::SparseMatrix *K) const
 
virtual void addFMDK (SurgSim::Math::Vector *F, SurgSim::Math::SparseMatrix *M, SurgSim::Math::SparseMatrix *D, SurgSim::Math::SparseMatrix *K) const
 Adds the element 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 (double alphaM, double alphaD, double alphaK, const SurgSim::Math::Vector &x, SurgSim::Math::Vector *F, SurgSim::Math::Vector *extractedX, SurgSim::Math::Vector *acumulator) const
 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) More...
 
bool isValidCoordinate (const SurgSim::Math::Vector &naturalCoordinate) const
 Determines whether a given natural coordinate is valid. More...
 
template<typename T , int Opt, typename StorageIndex >
void assembleMatrixBlocks (const Eigen::Ref< const Math::Matrix > &subMatrix, const std::vector< size_t > &blockIds, size_t blockSize, Eigen::SparseMatrix< T, Opt, StorageIndex > *matrix) const
 Add a sub-matrix made of squared-blocks into a matrix that could be un-initialized. More...
 
template<typename T , int Opt, typename StorageIndex >
void assembleMatrixBlocksNoInitialize (const Eigen::Ref< const Math::Matrix > &subMatrix, const std::vector< size_t > &blockIds, size_t blockSize, Eigen::SparseMatrix< T, Opt, StorageIndex > *matrix) const
 Add a sub-matrix made of squared-blocks into a matrix that is already initialized. More...
 
void updateFMDK (const Math::OdeState &state, int options)
 Update the FemElement based on the given state. More...
 

Protected Member Functions

void initializeMembers ()
 Initializes variables needed before Initialize() is called.
 
void computeInitialRotation (const SurgSim::Math::OdeState &state)
 Computes the beam element's initial rotation. More...
 
void computeStiffness (const SurgSim::Math::OdeState &state)
 Computes the beam's stiffness matrix. More...
 
void computeMass (const SurgSim::Math::OdeState &state)
 Computes the beam's mass matrix. More...
 
void doUpdateFMDK (const Math::OdeState &state, int options) override
 Update the FemElement based on the given state. More...
 
- Protected Member Functions inherited from SurgSim::Physics::FemElement
void setNumDofPerNode (size_t numDofPerNode)
 Sets the number of degrees of freedom per node. More...
 
void initializeFMDK ()
 Initialize f, M, D, K variables.
 
virtual void doInitializeFMDK ()
 Function to be overridden by the derived classes to initialize the f, M, D, K variables.
 

Protected Attributes

Eigen::Matrix< double, 12, 1 > m_x0
 The element's rest state.
 
SurgSim::Math::Matrix33d m_R0
 Initial rotation matrix for the element.
 
Eigen::Matrix< double, 12, 12 > m_MLocal
 Stiffness matrix (in local coordinate frame)
 
Eigen::Matrix< double, 12, 12 > m_KLocal
 Stiffness matrix (in local coordinate frame)
 
double m_restLength
 Rest length.
 
double m_radius
 radius for a circular Beam
 
double m_A
 Cross sectional area = PI.radius.radius if circular.
 
bool m_haveShear
 Does this beam element have shear.
 
double m_shearFactor
 Shear factor (usually 5/8)
 
double m_J
 Polar moment of inertia.
 
- Protected Attributes inherited from SurgSim::Physics::FemElement
size_t m_numDofPerNode
 Number of degree of freedom per node for this element.
 
std::vector< size_t > m_nodeIds
 Node ids connected by this element.
 
SurgSim::Math::Vector m_f
 The force vector.
 
SurgSim::Math::Matrix m_K
 The stiffness matrix.
 
SurgSim::Math::Matrix m_M
 The mass matrix.
 
SurgSim::Math::Matrix m_D
 The damping matrix.
 
bool m_useDamping
 Flag to specify of the damping is used.
 
double m_rho
 Mass density (in Kg.m-3)
 
double m_E
 Young modulus (in N.m-2)
 
double m_nu
 Poisson ratio (unitless)
 

Additional Inherited Members

- Public Types inherited from SurgSim::Physics::FemElement
typedef SurgSim::Framework::ObjectFactory1< FemElement, std::shared_ptr< FemElementStructs::FemElementParameter > > FactoryType
 
- Static Public Member Functions inherited from SurgSim::Physics::FemElement
static FactoryTypegetFactory ()
 

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

§ Fem1DElementBeam() [1/2]

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

Constructor.

Parameters
nodeIdsAn array of 2 node ids defining this beam element with respect to a DeformableRepresentationState which is passed to the initialize method.

§ Fem1DElementBeam() [2/2]

SurgSim::Physics::Fem1DElementBeam::Fem1DElementBeam ( std::shared_ptr< FemElementStructs::FemElementParameter elementData)
explicit

Constructor for FemElement object factory.

Parameters
elementDataA FemElement struct defining this beam element with respect to a DeformableRepresentationState which is passed to the initialize method.
Exceptions
SurgSim::Framework::AssertionFailureif nodeIds has a size different than 2

Member Function Documentation

§ computeCartesianCoordinate()

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

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.

§ computeInitialRotation()

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

§ computeMass()

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

Computes the beam's mass matrix.

Parameters
stateThe state to compute the stiffness matrix from

§ computeNaturalCoordinate()

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.

§ computeStiffness()

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

Computes the beam's stiffness matrix.

Parameters
stateThe state to compute the stiffness matrix from

Cross sectional moment of inertia

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

Shear deformation parameters

§ doUpdateFMDK()

void SurgSim::Physics::Fem1DElementBeam::doUpdateFMDK ( const Math::OdeState state,
int  options 
)
overrideprotectedvirtual

Update the FemElement based on the given state.

Parameters
state\((x, v)\) the current position and velocity to evaluate the various terms with
optionsFlag to specify which of the F, M, D, K needs to be updated

Implements SurgSim::Physics::FemElement.

§ getInitialRotation()

const SurgSim::Math::Matrix33d & SurgSim::Physics::Fem1DElementBeam::getInitialRotation ( ) const
Returns
the initial beam rotation

§ getRadius()

double SurgSim::Physics::Fem1DElementBeam::getRadius ( ) const

Gets the beam's circular cross-section radius.

Returns
The radius of the beam

§ getShearingEnabled()

bool SurgSim::Physics::Fem1DElementBeam::getShearingEnabled ( ) const

Gets whether shearing is enabled for the element.

Returns
True if shearing is enabled

§ getVolume()

double SurgSim::Physics::Fem1DElementBeam::getVolume ( const SurgSim::Math::OdeState state) const
overridevirtual

Gets the element volume based on the input state (in m-3)

Parameters
stateThe state to compute the volume with
Returns
The volume of this element (in m-3)

Implements SurgSim::Physics::FemElement.

§ initialize()

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

Initialize the FemElement once everything has been set.

Parameters
stateThe state to initialize the FemElement with

Reimplemented from SurgSim::Physics::FemElement.

§ setRadius()

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

Sets the beam's circular cross-section radius.

Parameters
radiusThe radius of the beam

§ setShearingEnabled()

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

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