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

2D FemElement based on a triangle with a constant thickness More...

#include <Fem2DElementTriangle.h>

Inheritance diagram for SurgSim::Physics::Fem2DElementTriangle:
SurgSim::Physics::FemElement MockFem2DElement

Public Member Functions

 Fem2DElementTriangle ()
 Constructor.
 
 Fem2DElementTriangle (std::array< size_t, 3 > nodeIds)
 Constructor. More...
 
 Fem2DElementTriangle (std::shared_ptr< FemElementStructs::FemElementParameter > elementData)
 Constructor for FemElement object factory. More...
 
void setThickness (double thickness)
 Sets the triangle's thickness. More...
 
double getThickness () const
 Gets the triangle's thickness. 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.
 
SurgSim::Math::Matrix33d computeRotation (const SurgSim::Math::OdeState &state)
 Computes the triangle element's rotation given a state. More...
 
virtual void computeLocalStiffness (const SurgSim::Math::OdeState &state, Eigen::Matrix< double, 18, 18 > *localStiffnessMatrix)
 Computes the triangle's local stiffness matrix. More...
 
void computeStiffness (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *stiffnessMatrix)
 Computes the triangle's stiffness matrix. More...
 
virtual void computeLocalMass (const SurgSim::Math::OdeState &state, Eigen::Matrix< double, 18, 18 > *localMassMatrix)
 Computes the triangle's local mass matrix. More...
 
void computeMass (const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *massMatrix)
 Computes the triangle's mass matrix. More...
 
void doUpdateFMDK (const Math::OdeState &state, int options) override
 Update the FemElement based on the given state. More...
 
void computeShapeFunctionsParameters (const SurgSim::Math::OdeState &restState)
 Compute the various shape functions (membrane and plate deformations) parameters. More...
 
std::array< double, 9 > batozDhxDxi (double xi, double eta) const
 Batoz derivative dHx/dxi. More...
 
std::array< double, 9 > batozDhxDeta (double xi, double eta) const
 Batoz derivative dHx/deta. More...
 
std::array< double, 9 > batozDhyDxi (double xi, double eta) const
 Batoz derivative dHy/dxi. More...
 
std::array< double, 9 > batozDhyDeta (double xi, double eta) const
 Batoz derivative dHy/deta. More...
 
Matrix39Type batozStrainDisplacement (double xi, double eta) const
 Batoz strain displacement matrix evaluated at a given point. 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, 18, 1 > m_x0
 The element's rest state.
 
SurgSim::Math::Matrix33d m_initialRotation
 Initial rotation matrix for the element.
 
Eigen::Matrix< double, 18, 18 > m_MLocal
 Stiffness matrix (in local coordinate frame)
 
Eigen::Matrix< double, 18, 18 > m_KLocal
 Stiffness matrix (in local coordinate frame)
 
double m_restArea
 The triangle rest area.
 
double m_thickness
 Thickness of the element.
 
SurgSim::Math::Matrix33d m_membraneShapeFunctionsParameters
 Membrane (in-plane) deformation. More...
 
SurgSim::Math::Vector3d m_xij
 Batoz variable \(x_{ij} = x_i - x_j\).
 
SurgSim::Math::Vector3d m_yij
 Batoz variable \(y_{ij} = y_i - y_j\).
 
SurgSim::Math::Vector3d m_lij_sqr
 Batoz variable \(l_{ij}^2 = x_{ij}^2 + y_{ij}^2\).
 
SurgSim::Math::Vector3d m_ak
 Batoz variable \(a_{k} = -x_{ij}/l_i^2\).
 
SurgSim::Math::Vector3d m_bk
 Batoz variable \(b_{k} = 3/4 x_{ij} y_{ij}/l_{ij}2\).
 
SurgSim::Math::Vector3d m_ck
 Batoz variable \(c_{k} = (1/4 x_{ij}^2 - 1/2 y_{ij}^2)/l_{ij}^2\).
 
SurgSim::Math::Vector3d m_dk
 Batoz variable \(d_{k} = -y_{ij}/l_{ij}^2\).
 
SurgSim::Math::Vector3d m_ek
 Batoz variable \(e_{k} = (1/4 y_{ij}^2 - 1/2 x_{ij}^2)/l_{ij}^2\).
 
SurgSim::Math::Vector3d m_Pk
 Batoz variable \(P_{k} = -6x_{ij}/l_{ij}^2 = 6.\textrm{m_a}_k\).
 
SurgSim::Math::Vector3d m_qk
 Batoz variable \(q_{k} = 3x_{ij}y_{ij}/l_{ij}^2 = 4.\textrm{m_b}_k\).
 
SurgSim::Math::Vector3d m_tk
 Batoz variable \(t_{k} = -6y_{ij}/l_{ij}^2 = 6.\textrm{m_d}_k\).
 
SurgSim::Math::Vector3d m_rk
 Batoz variable \(r_{k} = 3y_{ij}^2/l_{ij}^2\).
 
SurgSim::Math::Matrix m_integral_dT_d
 Plate mass matrix: integral terms related to the dof \((z)\).
 
SurgSim::Math::Matrix m_integralHyiHyj
 Plate mass matrix: integral terms related to the dof \((\theta_x)\).
 
SurgSim::Math::Matrix m_integralHxiHxj
 Plate mass matrix: integral terms related to the dof \((\theta_y)\).
 
- 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

2D FemElement based on a triangle with a constant thickness

The triangle is modeled as a shell (6DOF) which is decomposed into a membrane (in-plane 2DOF \((x, y)\)) and a plate (bending/twisting 3DOF \((z, \theta_x, \theta_y)\)). The thin-plate assumption does not consider the drilling dof \(\theta_z\). The system includes the drilling DOF for completeness but does not assign any mass or stiffness to it.

The membrane (in-plane) equations (mass and stiffness) are following "Theory of Matrix Structural Analysis" from J.S. Przemieniecki.

The thin-plate (bending) equations (mass and stiffness) are following "A Study Of Three-Node Triangular Plate Bending Elements", Jean-Louis Batoz Numerical Methods in Engineering, vol 15, 1771-1812 (1980).
The plate mass matrix is not detailed in the above paper, but the analytical equations have been derived from it. Moreover, to account for contribution of the displacement along z to the plate mass matrix, we used a cubic expression of this displacement given in: "Shell elements: modelizations DKT, DST, DKTG and Q4g", Code_Aster, 2013, Thomas De Soza.

Note
The element is considered to have a constant thickness.
The element uses linear elasticity (not visco-elasticity), so it does not have any damping.

Constructor & Destructor Documentation

§ Fem2DElementTriangle() [1/2]

SurgSim::Physics::Fem2DElementTriangle::Fem2DElementTriangle ( std::array< size_t, 3 >  nodeIds)
explicit

Constructor.

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

§ Fem2DElementTriangle() [2/2]

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

Constructor for FemElement object factory.

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

Member Function Documentation

§ batozDhxDeta()

std::array< double, 9 > SurgSim::Physics::Fem2DElementTriangle::batozDhxDeta ( double  xi,
double  eta 
) const
protected

Batoz derivative dHx/deta.

Parameters
xi,etaThe parametric coordinate (in [0 1] and xi+eta<1.0)
Returns
The vector dHx/deta evaluated at (xi, eta)

§ batozDhxDxi()

std::array< double, 9 > SurgSim::Physics::Fem2DElementTriangle::batozDhxDxi ( double  xi,
double  eta 
) const
protected

Batoz derivative dHx/dxi.

Parameters
xi,etaThe parametric coordinate (in [0 1] and xi+eta<1.0)
Returns
The vector dHx/dxi evaluated at (xi, eta)

§ batozDhyDeta()

std::array< double, 9 > SurgSim::Physics::Fem2DElementTriangle::batozDhyDeta ( double  xi,
double  eta 
) const
protected

Batoz derivative dHy/deta.

Parameters
xi,etaThe parametric coordinate (in [0 1] and xi+eta<1.0)
Returns
The vector dHy/deta evaluated at (xi, eta)

§ batozDhyDxi()

std::array< double, 9 > SurgSim::Physics::Fem2DElementTriangle::batozDhyDxi ( double  xi,
double  eta 
) const
protected

Batoz derivative dHy/dxi.

Parameters
xi,etaThe parametric coordinate (in [0 1] and xi+eta<1.0)
Returns
The vector dHy/dxi evaluated at (xi, eta)

§ batozStrainDisplacement()

Fem2DElementTriangle::Matrix39Type SurgSim::Physics::Fem2DElementTriangle::batozStrainDisplacement ( double  xi,
double  eta 
) const
protected

Batoz strain displacement matrix evaluated at a given point.

Parameters
xi,etaThe parametric coordinate (in [0 1] and xi+eta<1.0)
Returns
The 3x9 strain displacement matrix evaluated at (xi, eta)

§ computeCartesianCoordinate()

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

§ computeLocalMass()

void SurgSim::Physics::Fem2DElementTriangle::computeLocalMass ( const SurgSim::Math::OdeState state,
Eigen::Matrix< double, 18, 18 > *  localMassMatrix 
)
protectedvirtual

Computes the triangle's local mass matrix.

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

§ computeLocalStiffness()

void SurgSim::Physics::Fem2DElementTriangle::computeLocalStiffness ( const SurgSim::Math::OdeState state,
Eigen::Matrix< double, 18, 18 > *  localStiffnessMatrix 
)
protectedvirtual

Computes the triangle's local stiffness matrix.

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

§ computeMass()

void SurgSim::Physics::Fem2DElementTriangle::computeMass ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix massMatrix 
)
protected

Computes the triangle's mass matrix.

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

§ computeNaturalCoordinate()

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

§ computeRotation()

SurgSim::Math::Matrix33d SurgSim::Physics::Fem2DElementTriangle::computeRotation ( const SurgSim::Math::OdeState state)
protected

Computes the triangle element's rotation given a state.

Parameters
stateThe state to compute the rotation from
Returns
The rotation matrix of the element in the given state

§ computeShapeFunctionsParameters()

void SurgSim::Physics::Fem2DElementTriangle::computeShapeFunctionsParameters ( const SurgSim::Math::OdeState restState)
protected

Compute the various shape functions (membrane and plate deformations) parameters.

Parameters
restStatethe rest state to compute the shape functions parameters from

§ computeStiffness()

void SurgSim::Physics::Fem2DElementTriangle::computeStiffness ( const SurgSim::Math::OdeState state,
SurgSim::Math::Matrix stiffnessMatrix 
)
protected

Computes the triangle's stiffness matrix.

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

§ doUpdateFMDK()

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

§ getThickness()

double SurgSim::Physics::Fem2DElementTriangle::getThickness ( ) const

Gets the triangle's thickness.

Returns
The thickness of the triangle

§ getVolume()

double SurgSim::Physics::Fem2DElementTriangle::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::Fem2DElementTriangle::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.

§ setThickness()

void SurgSim::Physics::Fem2DElementTriangle::setThickness ( double  thickness)

Sets the triangle's thickness.

Parameters
thicknessThe thickness of the triangle

Member Data Documentation

§ m_membraneShapeFunctionsParameters

SurgSim::Math::Matrix33d SurgSim::Physics::Fem2DElementTriangle::m_membraneShapeFunctionsParameters
protected

Membrane (in-plane) deformation.

DOF simulated: (x, y). "Theory of Matrix Structural Analysis" from J.S. Przemieniecki. Shape functions are \(f_i(x, y) = a_i + b_i.x + c_i.y\), here we store \((a_i, b_i, c_i)\) on each row.


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