FemRepresentation.h
Go to the documentation of this file.
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_FEMREPRESENTATION_H
17 #define SURGSIM_PHYSICS_FEMREPRESENTATION_H
18 
19 #include <memory>
20 
22 #include "SurgSim/Math/Matrix.h"
24 #include "SurgSim/Math/Vector.h"
26 #include "SurgSim/Physics/Fem.h"
27 
28 namespace SurgSim
29 {
30 
31 namespace Physics
32 {
33 
34 class FemElement;
35 class FemPlyReaderDelegate;
36 
47 {
48 public:
51  explicit FemRepresentation(const std::string& name);
52 
54  virtual ~FemRepresentation();
55 
58  virtual void loadFem(const std::string& filename) = 0;
59 
62  virtual void setFemElementType(const std::string& type);
63 
66  const std::string& getFemElementType() const;
67 
70  void addFemElement(const std::shared_ptr<FemElement> element);
71 
74  size_t getNumFemElements() const;
75 
81  std::shared_ptr<FemElement> getFemElement(size_t femElementId);
82 
85  double getTotalMass() const;
86 
89  double getRayleighDampingStiffness() const;
90 
93  double getRayleighDampingMass() const;
94 
97  void setRayleighDampingStiffness(double stiffnessCoef);
98 
101  void setRayleighDampingMass(double massCoef);
102 
107 
110  void beforeUpdate(double dt) override;
111 
112  void update(double dt) override;
113 
118  void setComplianceWarping(bool useComplianceWarping);
119 
122  bool getComplianceWarping() const;
123 
128  Math::Matrix applyCompliance(const Math::OdeState& state, const Math::Matrix& b) override;
129 
130  const SurgSim::Math::Matrix& getComplianceMatrix() const override;
131 
132  void updateFMDK(const SurgSim::Math::OdeState& state, int options) override;
133 
134 protected:
147  bool useGlobalMassMatrix = false, bool useGlobalStiffnessMatrix = false,
148  double scale = 1.0);
149 
154  void addFemElementsForce(SurgSim::Math::Vector* f, const SurgSim::Math::OdeState& state, double scale = 1.0);
155 
161  void addGravityForce(SurgSim::Math::Vector* f, const SurgSim::Math::OdeState& state, double scale = 1.0);
162 
163  bool doInitialize() override;
164 
170 
175  virtual SurgSim::Math::Matrix getNodeTransformation(const SurgSim::Math::OdeState& state, size_t nodeId);
176 
180 
183  void setIsInitialComplianceMatrixComputed(bool flag);
184 
185  void computeF(const SurgSim::Math::OdeState& state) override;
186 
187  void computeM(const SurgSim::Math::OdeState& state) override;
188 
189  void computeD(const SurgSim::Math::OdeState& state) override;
190 
191  void computeK(const SurgSim::Math::OdeState& state) override;
192 
193  void computeFMDK(const SurgSim::Math::OdeState& state) override;
194 
196  std::vector<double> m_massPerNode;
197 
199  std::vector<std::shared_ptr<FemElement>> m_femElements;
200 
205  std::string m_femElementType;
206 
207 private:
211  struct
212  {
216 
218 
220 
222 
224  Eigen::SparseMatrix<double> m_complianceWarpingTransformation;
225 };
226 
227 } // namespace Physics
228 } // namespace SurgSim
229 
230 #endif // SURGSIM_PHYSICS_FEMREPRESENTATION_H
bool m_useComplianceWarping
Are we using Compliance Warping or not ?
Definition: FemRepresentation.h:217
Definition: CompoundShapeToGraphics.cpp:29
bool doInitialize() override
Interface to be implemented by derived classes.
Definition: FemRepresentation.cpp:75
A generic (size_t index, Vector coordinate) pair.
Definition: IndexedLocalCoordinate.h:29
virtual void loadFem(const std::string &filename)=0
Loads the FEM file into an Fem class data structure.
SurgSim::Math::Matrix m_complianceWarpingMatrix
The compliance warping matrix if compliance warping in use.
Definition: FemRepresentation.h:221
Finite Element Model (a.k.a FEM) is a deformable model (a set of nodes connected by FemElement)...
Definition: FemRepresentation.h:46
FemRepresentation(const std::string &name)
Constructor.
Definition: FemRepresentation.cpp:39
const SurgSim::Math::Matrix & getComplianceMatrix() const override
Gets the compliance matrix associated with motion.
Definition: FemRepresentation.cpp:293
size_t getNumFemElements() const
Gets the number of FemElement.
Definition: FemRepresentation.cpp:165
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > Matrix
A dynamic size matrix.
Definition: Matrix.h:65
void update(double dt) override
Update the representation state to the current time step.
Definition: FemRepresentation.cpp:222
bool m_isInitialComplianceMatrixComputed
For compliance warping: Is the initial compliance matrix computed ?
Definition: FemRepresentation.h:219
Base class for all deformable representations MassSprings, Finite Element Models,...
Definition: DeformableRepresentation.h:50
void setIsInitialComplianceMatrixComputed(bool flag)
Sets the flag keeping track of the initial compliance matrix calculation (compliance warping case) ...
Definition: FemRepresentation.cpp:565
Base class for all deformable representations (abstract class)
void addFemElementsForce(SurgSim::Math::Vector *f, const SurgSim::Math::OdeState &state, double scale=1.0)
Adds the FemElements forces to f (given a state)
Definition: FemRepresentation.cpp:527
double massCoefficient
Definition: FemRepresentation.h:213
The state of an ode of 2nd order of the form with boundary conditions.
Definition: OdeState.h:38
void setRayleighDampingStiffness(double stiffnessCoef)
Sets the Rayleigh stiffness parameter.
Definition: FemRepresentation.cpp:202
bool isInitialComplianceMatrixComputed() const
Gets the flag keeping track of the initial compliance matrix calculation (compliance warping case) ...
Definition: FemRepresentation.cpp:560
Definitions of useful sparse matrix functions.
bool isValidCoordinate(const SurgSim::DataStructures::IndexedLocalCoordinate &coordinate) const
Determines whether the associated coordinate is valid.
Definition: FemRepresentation.cpp:176
double stiffnessCoefficient
Definition: FemRepresentation.h:214
void updateFMDK(const SurgSim::Math::OdeState &state, int options) override
Update the OdeEquation (and support data) based on the given state.
Definition: FemRepresentation.cpp:467
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
A dynamic size column vector.
Definition: Vector.h:68
Math::Matrix applyCompliance(const Math::OdeState &state, const Math::Matrix &b) override
Calculate the product C.b where C is the compliance matrix with boundary conditions.
Definition: FemRepresentation.cpp:280
double getRayleighDampingMass() const
Gets the Rayleigh mass parameter.
Definition: FemRepresentation.cpp:197
virtual void setFemElementType(const std::string &type)
Sets the FemElement type pulled from the object factory.
Definition: FemRepresentation.cpp:64
std::vector< std::shared_ptr< FemElement > > m_femElements
FemElements.
Definition: FemRepresentation.h:199
void setRayleighDampingMass(double massCoef)
Sets the Rayleigh mass parameter.
Definition: FemRepresentation.cpp:207
Definitions of small fixed-size square matrix types.
void computeM(const SurgSim::Math::OdeState &state) override
Evaluation of the LHS matrix for a given state.
Definition: FemRepresentation.cpp:344
virtual SurgSim::Math::Matrix getNodeTransformation(const SurgSim::Math::OdeState &state, size_t nodeId)
Retrieves a specific node transformation (useful for compliance warping)
Definition: FemRepresentation.cpp:304
Definitions of small fixed-size vector types.
const std::string & getFemElementType() const
Gets the FemElement type pulled from the object factory.
Definition: FemRepresentation.cpp:70
void computeFMDK(const SurgSim::Math::OdeState &state) override
Evaluation of , , and .
Definition: FemRepresentation.cpp:411
void computeD(const SurgSim::Math::OdeState &state) override
Evaluation of for a given state.
Definition: FemRepresentation.cpp:355
std::vector< double > m_massPerNode
Useful information per node.
Definition: FemRepresentation.h:196
double getTotalMass() const
Gets the total mass of the fem.
Definition: FemRepresentation.cpp:182
void beforeUpdate(double dt) override
Preprocessing done before the update call.
Definition: FemRepresentation.cpp:212
virtual ~FemRepresentation()
Destructor.
Definition: FemRepresentation.cpp:60
void setComplianceWarping(bool useComplianceWarping)
Set the compliance warping flag.
Definition: FemRepresentation.cpp:268
bool getComplianceWarping() const
Get the compliance warping flag (default = false)
Definition: FemRepresentation.cpp:275
void addFemElement(const std::shared_ptr< FemElement > element)
Adds a FemElement.
Definition: FemRepresentation.cpp:160
double getRayleighDampingStiffness() const
Gets the Rayleigh stiffness parameter.
Definition: FemRepresentation.cpp:192
Eigen::SparseMatrix< double > m_complianceWarpingTransformation
The system-size transformation matrix. It contains nodes transformation on the diagonal blocks...
Definition: FemRepresentation.h:224
void addRayleighDampingForce(SurgSim::Math::Vector *f, const SurgSim::Math::OdeState &state, bool useGlobalMassMatrix=false, bool useGlobalStiffnessMatrix=false, double scale=1.0)
Adds the Rayleigh damping forces.
Definition: FemRepresentation.cpp:479
void computeF(const SurgSim::Math::OdeState &state) override
Evaluation of the RHS function for a given state.
Definition: FemRepresentation.cpp:328
struct SurgSim::Physics::FemRepresentation::@3 m_rayleighDamping
Rayleigh damping parameters (massCoefficient and stiffnessCoefficient) D = massCoefficient.M + stiffnessCoefficient.K Matrices: D = damping, M = mass, K = stiffness.
std::shared_ptr< FemElement > getFemElement(size_t femElementId)
Retrieves a given FemElement from its id.
Definition: FemRepresentation.cpp:170
std::string m_femElementType
The FemElement factory parameter invoked in doInitialize() when using a MeshAsset either with LoadFem...
Definition: FemRepresentation.h:205
void computeK(const SurgSim::Math::OdeState &state) override
Evaluation of for a given state.
Definition: FemRepresentation.cpp:394
void updateComplianceMatrix(const SurgSim::Math::OdeState &state)
Updates the compliance matrix using nodes transformation (useful for compliance warping) ...
Definition: FemRepresentation.cpp:312
void addGravityForce(SurgSim::Math::Vector *f, const SurgSim::Math::OdeState &state, double scale=1.0)
Adds the gravity force to f (given a state)
Definition: FemRepresentation.cpp:537