Go to the documentation of this file.
16 #ifndef SURGSIM_PHYSICS_FEM3DELEMENTTETRAHEDRON_H
17 #define SURGSIM_PHYSICS_FEM3DELEMENTTETRAHEDRON_H
31 SURGSIM_STATIC_REGISTRATION(Fem3DElementTetrahedron);
91 std::array<double, 4>* ai,
92 std::array<double, 4>* bi,
93 std::array<double, 4>* ci,
94 std::array<double, 4>* di)
const;
114 Eigen::Matrix<double, 12, 1>
m_x0;
117 Eigen::Matrix<double, 6, 6>
m_Em;
128 #endif // SURGSIM_PHYSICS_FEM3DELEMENTTETRAHEDRON_H
SurgSim::Math::Vector computeCartesianCoordinate(const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &naturalCoordinate) const override
Computes a given natural coordinate in cartesian coordinates.
Definition: Fem3DElementTetrahedron.cpp:353
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.
Definition: Fem3DElementTetrahedron.cpp:229
double m_restVolume
Shape functions: Tetrahedron rest volume.
Definition: Fem3DElementTetrahedron.h:109
Definitions of small fixed-size vector types.
The state of an ode of 2nd order of the form with boundary conditions.
Definition: OdeState.h:39
Eigen::Matrix< double, 6, 6 > m_Em
Elasticity material matrix (contains the elastic properties of the material)
Definition: Fem3DElementTetrahedron.h:117
Definitions of small fixed-size square matrix types.
Definition: CompoundShapeToGraphics.cpp:30
void computeStiffness(const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *k)
Computes the tetrahedron stiffness matrix.
Definition: Fem3DElementTetrahedron.cpp:171
std::array< double, 4 > m_di
Definition: Fem3DElementTetrahedron.h:111
Eigen::Matrix< double, 12, 1 > m_x0
The tetrahedon rest state.
Definition: Fem3DElementTetrahedron.h:114
void initializeMembers()
Initializes variables needed before Initialize() is called.
Definition: Fem3DElementTetrahedron.cpp:73
void computeMass(const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *m)
Computes the tetrahedron mass matrix.
Definition: Fem3DElementTetrahedron.cpp:126
Eigen::Matrix< double, 6, 12 > m_strain
Strain matrix.
Definition: Fem3DElementTetrahedron.h:119
std::array< double, 4 > m_bi
Definition: Fem3DElementTetrahedron.h:111
Eigen::Matrix< double, 6, 12 > m_stress
Stress matrix.
Definition: Fem3DElementTetrahedron.h:121
void doUpdateFMDK(const Math::OdeState &state, int options) override
Update the FemElement based on the given state.
Definition: Fem3DElementTetrahedron.cpp:78
void initialize(const SurgSim::Math::OdeState &state) override
Initialize the FemElement once everything has been set.
Definition: Fem3DElementTetrahedron.cpp:92
std::array< double, 4 > m_ci
Definition: Fem3DElementTetrahedron.h:111
double getVolume(const SurgSim::Math::OdeState &state) const override
Gets the element volume based on the input state (in m-3)
Definition: Fem3DElementTetrahedron.cpp:210
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
A dynamic size column vector.
Definition: Vector.h:68
SurgSim::Math::Vector computeNaturalCoordinate(const SurgSim::Math::OdeState &state, const SurgSim::Math::Vector &cartesianCoordinate) const override
Computes a natural coordinate given a global coordinate.
Definition: Fem3DElementTetrahedron.cpp:372
std::array< double, 4 > m_ai
Shape functions coefficients Ni(x,y,z) = 1/6V ( ai + x.bi + y.ci + z.di )
Definition: Fem3DElementTetrahedron.h:111
Base class for all Fem Element (1D, 2D, 3D) It handles the node ids to which it is connected and requ...
Definition: FemElement.h:46
Class for Fem Element 3D based on a tetrahedron volume discretization.
Definition: Fem3DElementTetrahedron.h:42
Fem3DElementTetrahedron()
Constructor.
Definition: Fem3DElementTetrahedron.cpp:50
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > Matrix
A dynamic size matrix.
Definition: Matrix.h:65
#define SURGSIM_CLASSNAME(ClassName)
Declare the class name of a class with the appropriate function header, do not use quotes.
Definition: Macros.h:21