AnimaL |
Tutorial |
Documentation |
#include <odeImplicitSolver.h>
Inheritance diagram for animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >:
Solves an ODE over a time step. Several explicit integration schemes inherited from parent class OdeSolver are available, along with implicit Euler (method number 4).
This class models a simplified system where damping forces are proportional with elastic forces, applying a uniform factor defined by the damping ratio.
To derive your own ode solver from this class, you must overload the following pure virtual methods:
In this class, it is assumed that matrix df/dx is updated by method computeForce( Vector&, const Positions&, const Vector&).
The solution of the equation system uses a simplified version of damping: the damping matrix is considered strictly proportional with the stiffness matrix. The factor is set using method setDampingRatio( Real ).
Definition at line 33 of file odeImplicitSolver.h.
Public Types | |
typedef t_Positions | Positions |
Vector of positions. | |
typedef t_Vector | Vector |
velocities, forces, accelerations | |
typedef t_Real | Real |
time and other scalar values | |
typedef container_traits< Vector >::value_type | Vec |
Public Member Functions | |
OdeImplicitSolver () | |
Constructor. | |
virtual void | solveODE (Positions &pos, Vector &vel, Real dt) |
Update positions and velocities by solving an ODE over time step dt. | |
void | setDampingRatio (Real d) |
Set the damping ratio used in the equation system. | |
Real | dampingRatio () const |
damping ratio used in the equation system | |
void | setMaxCGiterations (int n) |
Set the maximum number of CG iterations. | |
int | maxCGiterations () const |
Maximum number of CG iterations. | |
void | setSmallDenominatorThreshold (Real s) |
Set the other CG stop criterion. | |
Real | smallDenominatorThreshold () const |
CG stop criterion. | |
virtual void | applyForces (Vector &a, const Vector &f)=0 |
Compute the acceleration of the particles when given forces are applied. | |
virtual void | computeAccelerations (Vector &a, const Positions &p, const Vector &v)=0 |
Compute the acceleration of the particles for given positions and velocities. | |
virtual void | computeAccelerationsAndStiffness (Vector &a, const Positions &p, const Vector &v)=0 |
Compute the acceleration of the particles for given positions and velocities, and updates the gradient matrix. | |
virtual void | v_eq_h_dfdx_x (Vector &v, Real h, const Vector &x)=0 |
Compute ![]() | |
virtual void | v_peq_h_dfdx_x (Vector &v, Real h, const Vector &x)=0 |
Compute ![]() | |
void | compute_implicit_euler_step (Positions &dpos, Vector &dvel, const Positions &pos, const Vector &vel, Real dt) |
Compute position and velocity updates using implicit Euler time integration. | |
void | implicit_euler (Positions &pos, Vector &vel, Real dt) |
Update positions and velocities using implicit Euler time integration. | |
void | multImplicitMatrix (Vector &y, const Vector &x, Real h) |
Compute ![]() | |
Public Attributes | |
Real | _dampingRatio |
damping ratio used in the equation system | |
int | _maxIter |
Maximum number of CG iterations. | |
Real | _smallDenominatorThreshold |
the other CG stop criterion | |
Static Public Attributes | |
const int | IMPLICIT_EULER = 5 |
Symbol of implicit Euler integration. |
|
Vector of positions.
Reimplemented from animal::OdeSolver< t_Positions, t_Vector, t_Real >. Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. Definition at line 37 of file odeImplicitSolver.h. |
|
time and other scalar values
Reimplemented from animal::OdeSolver< t_Positions, t_Vector, t_Real >. Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. Definition at line 39 of file odeImplicitSolver.h. |
|
|
velocities, forces, accelerations
Reimplemented from animal::OdeSolver< t_Positions, t_Vector, t_Real >. Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. Definition at line 38 of file odeImplicitSolver.h. |
|
Constructor.
Definition at line 16 of file odeImplicitSolver.inl. References animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::IMPLICIT_EULER. |
|
Compute the acceleration of the particles when given forces are applied.
Implemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::compute_implicit_euler_step(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::implicit_euler(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::multImplicitMatrix(). |
|
Compute position and velocity updates using implicit Euler time integration.
Definition at line 126 of file odeImplicitSolver.inl. References animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::_maxIter, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::_smallDenominatorThreshold, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::applyForces(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::computeAccelerationsAndStiffness(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::multImplicitMatrix(), animal::v_assign(), animal::v_dot(), animal::v_eq(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::v_eq_h_dfdx_x(), animal::v_peq(), and animal::v_teq(). |
|
Compute the acceleration of the particles for given positions and velocities.
Implements animal::OdeSolver< t_Positions, t_Vector, t_Real >. Implemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. |
|
Compute the acceleration of the particles for given positions and velocities, and updates the gradient matrix.
Implemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::compute_implicit_euler_step(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::implicit_euler(). |
|
damping ratio used in the equation system
Definition at line 58 of file odeImplicitSolver.h. Referenced by X3DTK::Qt::MassSpringEngineQt::MassSpringEngineQt(). |
|
Update positions and velocities using implicit Euler time integration.
Definition at line 46 of file odeImplicitSolver.inl. References animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::_maxIter, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::_smallDenominatorThreshold, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::applyForces(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::computeAccelerationsAndStiffness(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::multImplicitMatrix(), animal::v_assign(), animal::v_dot(), animal::v_eq(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::v_eq_h_dfdx_x(), animal::v_peq(), and animal::v_teq(). Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::solveODE(). |
|
Maximum number of CG iterations.
Definition at line 66 of file odeImplicitSolver.h. Referenced by X3DTK::Qt::MassSpringEngineQt::MassSpringEngineQt(). |
|
Compute Used internally by the CG solution.
Definition at line 202 of file odeImplicitSolver.inl. References animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::_dampingRatio, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::applyForces(), animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::v_eq_h_dfdx_x(), and animal::v_peq(). Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::compute_implicit_euler_step(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::implicit_euler(). |
|
Set the damping ratio used in the equation system. The damping factor is the product of the stiffness with the damping ratio.
Definition at line 55 of file odeImplicitSolver.h. Referenced by X3DTK::Qt::MassSpringEngineQt::setImplicitEulerDampingRatio(). |
|
Set the maximum number of CG iterations.
Definition at line 63 of file odeImplicitSolver.h. Referenced by X3DTK::Qt::MassSpringEngineQt::setMaxCGiterations(). |
|
Set the other CG stop criterion.
Definition at line 71 of file odeImplicitSolver.h. Referenced by X3DTK::Qt::MassSpringEngineQt::setImplicitEulerDividerThreshold(). |
|
CG stop criterion.
Definition at line 74 of file odeImplicitSolver.h. Referenced by X3DTK::Qt::MassSpringEngineQt::MassSpringEngineQt(). |
|
Update positions and velocities by solving an ODE over time step dt.
Reimplemented from animal::OdeSolver< t_Positions, t_Vector, t_Real >. Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. Definition at line 31 of file odeImplicitSolver.inl. References animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::implicit_euler(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::IMPLICIT_EULER. |
|
|
Compute
Implemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, and animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >. |
|
damping ratio used in the equation system
Definition at line 80 of file odeImplicitSolver.h. Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::multImplicitMatrix(). |
|
Maximum number of CG iterations.
Definition at line 81 of file odeImplicitSolver.h. Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::compute_implicit_euler_step(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::implicit_euler(). |
|
the other CG stop criterion
Definition at line 82 of file odeImplicitSolver.h. Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::compute_implicit_euler_step(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::implicit_euler(). |
|
Symbol of implicit Euler integration.
Definition at line 40 of file odeImplicitSolver.h. Referenced by animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::OdeImplicitSolver(), and animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >::solveODE(). |