opensurgsim
Fem3DElementCube.h
1 // This file is a part of the OpenSurgSim project.
2 // Copyright 2013, SimQuest Solutions Inc.
3 //
4 // Licensed under the Apache License, Version 2.0 (the "License");
5 // you may not use this file except in compliance with the License.
6 // You may obtain a copy of the License at
7 //
8 // http://www.apache.org/licenses/LICENSE-2.0
9 //
10 // Unless required by applicable law or agreed to in writing, software
11 // distributed under the License is distributed on an "AS IS" BASIS,
12 // WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied.
13 // See the License for the specific language governing permissions and
14 // limitations under the License.
15 
16 #ifndef SURGSIM_PHYSICS_FEM3DELEMENTCUBE_H
17 #define SURGSIM_PHYSICS_FEM3DELEMENTCUBE_H
18 
19 #include <array>
20 
21 #include "SurgSim/Physics/FemElement.h"
23 
24 namespace SurgSim
25 {
26 
27 namespace Physics
28 {
29 SURGSIM_STATIC_REGISTRATION(Fem3DElementCube);
30 
43 {
44 public:
47 
56  explicit Fem3DElementCube(std::array<size_t, 8> nodeIds);
57 
67  explicit Fem3DElementCube(std::shared_ptr<FemElementStructs::FemElementParameter> elementData);
68 
69  SURGSIM_CLASSNAME(SurgSim::Physics::Fem3DElementCube)
70 
71 
72  void initialize(const SurgSim::Math::OdeState& state) override;
82 
83  double getVolume(const SurgSim::Math::OdeState& state) const override;
84 
86  const SurgSim::Math::Vector& naturalCoordinate) const override;
87 
89  const SurgSim::Math::Vector& cartesianCoordinate) const override;
90 
91 protected:
93  void initializeMembers();
94 
97  void buildConstitutiveMaterialMatrix(Eigen::Matrix<double, 6, 6>* constitutiveMatrix);
98 
102  void computeStiffness(const SurgSim::Math::OdeState& state,
103  Eigen::Matrix<double, 6, 24>* strain,
104  Eigen::Matrix<double, 6, 24>* stress,
106 
111 
112  void doUpdateFMDK(const Math::OdeState& state, int options) override;
113 
123  Eigen::Matrix<double, 6, 24>* strain,
124  Eigen::Matrix<double, 6, 24>* stress,
126 
136 
142  void evaluateJ(const SurgSim::Math::OdeState& state, double epsilon, double eta, double mu,
145  double* detJ) const;
146 
152  void evaluateStrainDisplacement(double epsilon, double eta, double mu, const SurgSim::Math::Matrix33d& Jinv,
153  Eigen::Matrix<double, 6, 24>* B) const;
154 
156  double m_restVolume;
157 
166  std::array<double, 8> m_shapeFunctionsEpsilonSign;
167  std::array<double, 8> m_shapeFunctionsEtaSign;
168  std::array<double, 8> m_shapeFunctionsMuSign;
170 
173 
188  double shapeFunction(size_t i, double epsilon, double eta, double mu) const;
195 
203  double dShapeFunctiondepsilon(size_t i, double epsilon, double eta, double mu) const;
204 
212  double dShapeFunctiondeta(size_t i, double epsilon, double eta, double mu) const;
213 
221  double dShapeFunctiondmu(size_t i, double epsilon, double eta, double mu) const;
222 
224  Eigen::Matrix<double, 24, 1> m_elementRestPosition;
225 
227  Eigen::Matrix<double, 6, 24> m_strain;
229  Eigen::Matrix<double, 6, 24> m_stress;
231  Eigen::Matrix<double, 6, 6> m_constitutiveMaterial;
232 };
233 
234 } // namespace Physics
235 
236 } // namespace SurgSim
237 
238 #endif // SURGSIM_PHYSICS_FEM3DELEMENTCUBE_H
Wraps glewInit() to separate the glew opengl definitions from the osg opengl definitions only imgui n...
Definition: AddRandomSphereBehavior.cpp:36
void computeStiffness(const SurgSim::Math::OdeState &state, Eigen::Matrix< double, 6, 24 > *strain, Eigen::Matrix< double, 6, 24 > *stress, SurgSim::Math::Matrix *k)
Computes the cube stiffness matrix along with the strain and stress matrices.
Definition: Fem3DElementCube.cpp:143
void addStrainStressStiffnessAtPoint(const SurgSim::Math::OdeState &state, const SurgSim::Math::gaussQuadraturePoint &epsilon, const SurgSim::Math::gaussQuadraturePoint &eta, const SurgSim::Math::gaussQuadraturePoint &mu, Eigen::Matrix< double, 6, 24 > *strain, Eigen::Matrix< double, 6, 24 > *stress, SurgSim::Math::Matrix *k)
Helper method to evaluate strain-stress and stiffness integral terms with a discrete sum using a Gaus...
Definition: Fem3DElementCube.cpp:267
void addMassMatrixAtPoint(const SurgSim::Math::OdeState &state, const SurgSim::Math::gaussQuadraturePoint &epsilon, const SurgSim::Math::gaussQuadraturePoint &eta, const SurgSim::Math::gaussQuadraturePoint &mu, SurgSim::Math::Matrix *m)
Helper method to evaluate mass integral terms with a discrete sum using a Gauss quadrature rule...
Definition: Fem3DElementCube.cpp:290
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: Fem3DElementCube.cpp:419
Eigen::Matrix< double, 24, 1 > m_elementRestPosition
The cube rest state (nodes ordered by m_nodeIds)
Definition: Fem3DElementCube.h:224
void evaluateStrainDisplacement(double epsilon, double eta, double mu, const SurgSim::Math::Matrix33d &Jinv, Eigen::Matrix< double, 6, 24 > *B) const
Helper method to evaluate the strain-displacement matrix at a given 3D parametric location c...
Definition: Fem3DElementCube.cpp:213
The state of an ode of 2nd order of the form with boundary conditions.
Definition: OdeState.h:38
Fem3DElementCube()
Constructor.
Definition: Fem3DElementCube.cpp:34
void initializeMembers()
Initializes variables needed before Initialize() is called.
Definition: Fem3DElementCube.cpp:59
Base class for all Fem Element (1D, 2D, 3D) It handles the node ids to which it is connected and requ...
Definition: FemElement.h:45
double dShapeFunctiondeta(size_t i, double epsilon, double eta, double mu) const
Shape functions derivative .
Definition: Fem3DElementCube.cpp:384
Definitions of a n-point Gaussian quadrature rule (a.k.a.
Class for Fem Element 3D based on a cube volume discretization.
Definition: Fem3DElementCube.h:42
void doUpdateFMDK(const Math::OdeState &state, int options) override
Update the FemElement based on the given state.
Definition: Fem3DElementCube.cpp:129
double shapeFunction(size_t i, double epsilon, double eta, double mu) const
Shape functions .
Definition: Fem3DElementCube.cpp:368
std::array< double, 8 > m_shapeFunctionsEpsilonSign
Definition: Fem3DElementCube.h:166
Eigen::Matrix< double, 6, 6 > m_constitutiveMaterial
Constitutive material matrix (Hooke&#39;s law in this case) defines the relationship between stress and s...
Definition: Fem3DElementCube.h:231
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
A dynamic size column vector.
Definition: Vector.h:68
void buildConstitutiveMaterialMatrix(Eigen::Matrix< double, 6, 6 > *constitutiveMatrix)
Build the constitutive material 6x6 matrix.
Definition: Fem3DElementCube.cpp:253
void initialize(const SurgSim::Math::OdeState &state) override
Initializes the element once everything has been set.
Definition: Fem3DElementCube.cpp:65
double dShapeFunctiondepsilon(size_t i, double epsilon, double eta, double mu) const
Shape functions derivative .
Definition: Fem3DElementCube.cpp:376
double getVolume(const SurgSim::Math::OdeState &state) const override
Gets the element volume based on the input state (in m-3)
Definition: Fem3DElementCube.cpp:321
void computeMass(const SurgSim::Math::OdeState &state, SurgSim::Math::Matrix *m)
Computes the cube mass matrix.
Definition: Fem3DElementCube.cpp:97
Eigen::Matrix< double, 3, 3, Eigen::RowMajor > Matrix33d
A 3x3 matrix of doubles.
Definition: Matrix.h:51
void evaluateJ(const SurgSim::Math::OdeState &state, double epsilon, double eta, double mu, SurgSim::Math::Matrix33d *J, SurgSim::Math::Matrix33d *Jinv, double *detJ) const
Helper method to evaluate matrix J = d(x,y,z)/d(epsilon,eta,mu) at a given 3D parametric location J e...
Definition: Fem3DElementCube.cpp:169
double dShapeFunctiondmu(size_t i, double epsilon, double eta, double mu) const
Shape functions derivative .
Definition: Fem3DElementCube.cpp:392
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: Fem3DElementCube.cpp:400
double m_restVolume
Cube rest volume.
Definition: Fem3DElementCube.h:156
Eigen::Matrix< double, 6, 24 > m_stress
Stress matrix (usually noted )
Definition: Fem3DElementCube.h:229
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > Matrix
A dynamic size matrix.
Definition: Matrix.h:65
1D Gauss-Legendre quadrature
Definition: GaussLegendreQuadrature.h:32
Eigen::Matrix< double, 6, 24 > m_strain
Strain matrix (usually noted )
Definition: Fem3DElementCube.h:227