16 #ifndef SURGSIM_PHYSICS_FEMREPRESENTATION_H 17 #define SURGSIM_PHYSICS_FEMREPRESENTATION_H 22 #include "SurgSim/DataStructures/IndexedLocalCoordinate.h" 27 #include "SurgSim/Physics/Fem.h" 58 virtual void loadFem(
const std::string& filename) = 0;
70 void addFemElement(
const std::shared_ptr<FemElement> element);
81 std::shared_ptr<FemElement>
getFemElement(
size_t femElementId);
112 void update(
double dt)
override;
163 bool useGlobalMassMatrix =
false,
bool useGlobalStiffnessMatrix =
false,
242 std::future<SurgSim::Math::Matrix> m_task;
251 double massCoefficient;
252 double stiffnessCoefficient;
255 bool m_useMassLumping;
261 #endif // SURGSIM_PHYSICS_FEMREPRESENTATION_H bool m_useComplianceWarping
Are we using Compliance Warping or not ?
Definition: FemRepresentation.h:227
Wraps glewInit() to separate the glew opengl definitions from the osg opengl definitions only imgui n...
Definition: AddRandomSphereBehavior.cpp:36
bool doInitialize() override
Interface to be implemented by derived classes.
Definition: FemRepresentation.cpp:85
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:234
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:38
const SurgSim::Math::Matrix & getComplianceMatrix() const override
Gets the compliance matrix associated with motion.
Definition: FemRepresentation.cpp:342
size_t getNumFemElements() const
Gets the number of FemElement.
Definition: FemRepresentation.cpp:175
Eigen::SparseMatrix< double > m_complianceWarpingTransformationForCalculation
Definition: FemRepresentation.h:241
void update(double dt) override
Update the representation state to the current time step.
Definition: FemRepresentation.cpp:232
bool m_isInitialComplianceMatrixComputed
For compliance warping: Is the initial compliance matrix computed ?
Definition: FemRepresentation.h:232
void setIsInitialComplianceMatrixComputed(bool flag)
Sets the flag keeping track of the initial compliance matrix calculation (compliance warping case) ...
Definition: FemRepresentation.cpp:669
virtual void calculateComplianceWarpingTransformation(const SurgSim::Math::OdeState &state)
Calculates and stores the compliance warping transformation matrix.
Definition: FemRepresentation.cpp:362
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:631
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:212
bool isInitialComplianceMatrixComputed() const
Gets the flag keeping track of the initial compliance matrix calculation (compliance warping case) ...
Definition: FemRepresentation.cpp:664
bool isComplianceWarpingSynchronous() const
Get whether or not any compliance warping will be done synchronously.
Definition: FemRepresentation.cpp:314
Definitions of useful sparse matrix functions.
virtual bool isValidCoordinate(const SurgSim::DataStructures::IndexedLocalCoordinate &coordinate) const
Determines whether the associated coordinate is valid.
Definition: FemRepresentation.cpp:186
void updateFMDK(const SurgSim::Math::OdeState &state, int options) override
Update the OdeEquation (and support data) based on the given state.
Definition: FemRepresentation.cpp:547
void setMassLumping(bool useMassLumping)
Enable mass lumping for the mass matrix, currently we use the row sum i.e.
Definition: FemRepresentation.cpp:320
Math::Matrix applyCompliance(const Math::OdeState &state, const Math::Matrix &b) override
Calculate the product where is the compliance matrix with boundary conditions applied.
Definition: FemRepresentation.cpp:329
double getRayleighDampingMass() const
Gets the Rayleigh mass parameter.
Definition: FemRepresentation.cpp:207
void setComplianceWarpingSynchronous(bool complianceWarpingSynchronous)
If using compliance warping, set the calculation to synchronous or asynchronous.
Definition: FemRepresentation.cpp:308
virtual void setFemElementType(const std::string &type)
Sets the FemElement type pulled from the object factory.
Definition: FemRepresentation.cpp:74
Eigen::Matrix< double, Eigen::Dynamic, 1 > Vector
A dynamic size column vector.
Definition: Vector.h:68
std::vector< std::shared_ptr< FemElement > > m_femElements
FemElements.
Definition: FemRepresentation.h:219
void setRayleighDampingMass(double massCoef)
Sets the Rayleigh mass parameter.
Definition: FemRepresentation.cpp:217
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:412
Definitions of small fixed-size vector types.
const std::string & getFemElementType() const
Gets the FemElement type pulled from the object factory.
Definition: FemRepresentation.cpp:80
void computeFMDK(const SurgSim::Math::OdeState &state) override
Evaluation of , , and .
Definition: FemRepresentation.cpp:479
void computeD(const SurgSim::Math::OdeState &state) override
Evaluation of for a given state.
Definition: FemRepresentation.cpp:423
std::vector< double > m_massPerNode
Useful information per node.
Definition: FemRepresentation.h:216
double getTotalMass() const
Gets the total mass of the fem.
Definition: FemRepresentation.cpp:192
void beforeUpdate(double dt) override
Preprocessing done before the update call.
Definition: FemRepresentation.cpp:222
virtual ~FemRepresentation()
Destructor.
Definition: FemRepresentation.cpp:66
void setComplianceWarping(bool useComplianceWarping)
Set the compliance warping flag.
Definition: FemRepresentation.cpp:296
bool getMassLumping() const
Definition: FemRepresentation.cpp:325
bool getComplianceWarping() const
Get the compliance warping flag (default = false)
Definition: FemRepresentation.cpp:303
void addFemElement(const std::shared_ptr< FemElement > element)
Adds a FemElement.
Definition: FemRepresentation.cpp:170
double getRayleighDampingStiffness() const
Gets the Rayleigh stiffness parameter.
Definition: FemRepresentation.cpp:202
Eigen::SparseMatrix< double > m_complianceWarpingTransformation
The system-size transformation matrix. It contains nodes transformation on the diagonal blocks...
Definition: FemRepresentation.h:237
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:565
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > Matrix
A dynamic size matrix.
Definition: Matrix.h:65
void computeF(const SurgSim::Math::OdeState &state) override
Evaluation of the RHS function for a given state.
Definition: FemRepresentation.cpp:396
bool m_isComplianceWarpingSynchronous
Is the compliance warping computation be synchronous (slower but up-to-date results)?
Definition: FemRepresentation.h:230
virtual SurgSim::Math::Matrix getNodeTransformation(const SurgSim::Math::OdeState &state, size_t nodeId) const
Retrieves a specific node transformation (useful for compliance warping)
Definition: FemRepresentation.cpp:353
std::shared_ptr< FemElement > getFemElement(size_t femElementId)
Retrieves a given FemElement from its id.
Definition: FemRepresentation.cpp:180
std::string m_femElementType
The FemElement factory parameter invoked in doInitialize() when using a MeshAsset either with LoadFem...
Definition: FemRepresentation.h:225
void computeK(const SurgSim::Math::OdeState &state) override
Evaluation of for a given state.
Definition: FemRepresentation.cpp:462
void updateComplianceMatrix(const SurgSim::Math::OdeState &state)
Updates the compliance matrix using nodes transformation (useful for compliance warping) ...
Definition: FemRepresentation.cpp:368
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:641