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

Class for Fem Element 3D based on a tetrahedron volume discretization. More...

#include <Fem3DElementTetrahedron.h>

Inheritance diagram for SurgSim::Physics::Fem3DElementTetrahedron:
SurgSim::Physics::FemElement MockFem3DElementTet SurgSim::Physics::Fem3DElementCorotationalTetrahedron MockFem3DElementCorotationalTet

Public Member Functions

 Fem3DElementTetrahedron ()
 Constructor.
 
 Fem3DElementTetrahedron (std::array< size_t, 4 > nodeIds)
 Constructor. More...
 
 Fem3DElementTetrahedron (std::shared_ptr< FemElementStructs::FemElementParameter > elementData)
 Constructor for FemElement object factory. 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...
 
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...
 
- 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 doUpdateFMDK (const Math::OdeState &state, int options) override
 Update the FemElement based on the given state. More...
 
void computeShapeFunctions (const SurgSim::Math::OdeState &state, double *volume, std::array< double, 4 > *ai, std::array< double, 4 > *bi, std::array< double, 4 > *ci, std::array< double, 4 > *di) const
 Computes the tetrahedron shape functions. More...
 
void computeStiffness (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *k)
 Computes the tetrahedron stiffness matrix. More...
 
void computeMass (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *m)
 Computes the tetrahedron 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...
 
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 tetrahedon rest state.
 
- 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

Class for Fem Element 3D based on a tetrahedron volume discretization.

Note
The inertia property (mass) of the tetrahedron is derived from
"Theory of Matrix Structural Analysis" from J.S. Przemieniecki
The force and stiffness matrix of the tetrahedron is derived from
http://www.colorado.edu/engineering/CAS/courses.d/AFEM.d/AFEM.Ch09.d/AFEM.Ch09.pdf
The deformation is based on the linear elasticity theory and not on the visco-elasticity theory.
Therefore the element does not have any damping component.
Fem3DElementTetrahedron uses linear elasticity (not visco-elasticity), so it does not give any damping.

Constructor & Destructor Documentation

§ Fem3DElementTetrahedron() [1/2]

SurgSim::Physics::Fem3DElementTetrahedron::Fem3DElementTetrahedron ( std::array< size_t, 4 >  nodeIds)
explicit

Constructor.

Parameters
nodeIdsAn array of 4 node ids defining this tetrahedron element in a overall mesh
Note
It is required that the triangle ABC is CCW looking from D (i.e. dot(cross(AB, AC), AD) > 0)
This is required from the signed volume calculation method getVolume()
A warning will be logged when the initialize function is called if this condition is not met, but the
simulation will keep running. Behavior will be undefined because of possible negative volume terms.

§ Fem3DElementTetrahedron() [2/2]

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

Constructor for FemElement object factory.

Parameters
elementDataA FemElement3D struct defining this tetrahedron element in a overall mesh
Note
It is required that the triangle ABC is CCW looking from D (i.e. dot(cross(AB, AC), AD) > 0)
This is required from the signed volume calculation method getVolume()
A warning will be logged when the initialize function is called if this condition is not met, but the
simulation will keep running. Behavior will be undefined because of possible negative volume terms.
Exceptions
SurgSim::Framework::AssertionFailureif nodeIds has a size different than 4

Member Function Documentation

§ computeCartesianCoordinate()

SurgSim::Math::Vector SurgSim::Physics::Fem3DElementTetrahedron::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.

§ computeMass()

void SurgSim::Physics::Fem3DElementTetrahedron::computeMass ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix m 
)
protected

Computes the tetrahedron mass matrix.

Parameters
stateThe state to compute the mass matrix from
[out]mThe mass matrix to store the result into

§ computeNaturalCoordinate()

SurgSim::Math::Vector SurgSim::Physics::Fem3DElementTetrahedron::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.

§ computeShapeFunctions()

void SurgSim::Physics::Fem3DElementTetrahedron::computeShapeFunctions ( const SurgSim::Math::OdeState state,
double *  volume,
std::array< double, 4 > *  ai,
std::array< double, 4 > *  bi,
std::array< double, 4 > *  ci,
std::array< double, 4 > *  di 
) const
protected

Computes the tetrahedron shape functions.

Parameters
stateThe deformable rest state to compute the shape function from
[out]volumethe volume calculated with the given state
[out]aifrom the shape function, Ni(x, y, z) = 1/6*volume (ai + bi.x + ci.y + di.z)
[out]bifrom the shape function, Ni(x, y, z) = 1/6*volume (ai + bi.x + ci.y + di.z)
[out]cifrom the shape function, Ni(x, y, z) = 1/6*volume (ai + bi.x + ci.y + di.z)
[out]difrom the shape function, Ni(x, y, z) = 1/6*volume (ai + bi.x + ci.y + di.z)

§ computeStiffness()

void SurgSim::Physics::Fem3DElementTetrahedron::computeStiffness ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix k 
)
protected

Computes the tetrahedron stiffness matrix.

Parameters
stateThe state to compute the stiffness matrix from
[out]kThe stiffness matrix to store the result into

Elasticity material matrix (contains the elastic properties of the material)

Shape functions: Tetrahedron rest volume

Shape functions coefficients Ni(x,y,z) = 1/6V ( ai + x.bi + y.ci + z.di )

§ doUpdateFMDK()

void SurgSim::Physics::Fem3DElementTetrahedron::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.

§ getVolume()

double SurgSim::Physics::Fem3DElementTetrahedron::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)

Computes the tetrahedron volume 1/6 * | 1 p0x p0y p0z | | 1 p1x p1y p1z | | 1 p2x p2y p2z | | 1 p3x p3y p3z |

Implements SurgSim::Physics::FemElement.

§ initialize()

void SurgSim::Physics::Fem3DElementTetrahedron::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.


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