Documentation


animal::OdeSolver< t_Positions, t_Vector, t_Real > Class Template Reference

#include <odeSolver.h>

Inheritance diagram for animal::OdeSolver< t_Positions, t_Vector, t_Real >:

animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real > animal::OdeImplicitSolver< Points, Velocities, value_type< SpringStiffness >::type > animal::OdeImplicitSolver< t_Positions, t_Vector, value_type< t_VecReal >::type > animal::OdeImplicitSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, value_type< X3DTK::MFFloat >::type > animal::SimpleSolver< t_Positions, t_Vector, t_InvMasses, t_Real > animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex > animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings > animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 > animal::MassSpringEngine< X3DTK::MFFloat, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFInt32, X3DTK::MFVec3f, float > animal::MassSpringEngine< Inv_Masses, Velocities, SpringStiffness, IndexedSprings, Points, Real > X3DTK::X3D::AnimalEngineNode< AMassSpringEngine > X3DTK::X3D::MassSpringNode List of all members.

Detailed Description

template<class t_Positions, class t_Vector, class t_Real>
class animal::OdeSolver< t_Positions, t_Vector, t_Real >

Abstract class for ODE solvers.

Solves an ODE over a time step. Four explicit integration schemes are available. To derive your own ode solver from this class, you must overload the pure virtual method computeAcceleration( Vector&, const Positions&, const Vector&).

Choose the integration scheme using method setMethod( int ) with one of the following parameters:

Author:
François Faure, Oct. 2003

Definition at line 23 of file odeSolver.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


Public Member Functions

 OdeSolver ()
 Constructor.

virtual ~OdeSolver ()
 Destructor.

virtual void solveODE (Positions &pos, Vector &vel, Real dt)
 Update positions and velocities by solving an ODE over time step dt.

void setMethod (int m)
 Set the method.

int method () const
 Get the method.

void setMMIDsteps (int n)
 Set the number of internal time steps in the modified midpoint method.

int MMIDsteps () const
 The number of internal time steps in the modified midpoint method.


Public Attributes

int _method
 current integration method

int _MMIDsteps
 Number of substeps in MMID.


Static Public Attributes

const int EULER = 0
 Symbol of explicit Euler integration.

const int RK2 = 1
 Symbol of second-order Runge-Kutta integration.

const int RK4 = 2
 Symbol of fourth-order Runge-Kutta integration.

const int MODMID = 3
 Symbol of modified midpoint integration.

const int VVERLET = 4
 Symbol of velocity Verlet integration.


Protected Types

typedef Derivs< Vector, VectorDer
 Derivs.

typedef States< Positions,
Vector
Sta
 States.


Protected Member Functions

virtual void computeAccelerations (Vector &acceleration, const Positions &pos, const Vector &vel)=0
 Compute the acceleration for given positions and velocities.

void computeDerivative (Der &d, const Sta &s, Real)
 Compute the derivative(velocity, acceleration) of a given state(positions, velocities).

void euler (Positions &pos, Vector &vel, Real h, Real t=0)
 Apply euler time integration.

void rk2 (Positions &pos, Vector &vel, Real h, Real t=0)
 Apply second-order Runge-Kutta time integration.

void rk4 (Positions &pos, Vector &vel, Real h, Real t=0)
 Apply fourth-order Runge-Kutta time integration.

void modmid (Positions &pos, Vector &vel, Real h, int n, Real t=0)
 Apply modified midpoint time integration.

void velocityVerlet (Positions &pos, Vector &vel, Real h)
 Apply velocity Verlet time integration.

void integrate_euler (Sta &s, Real h, Real t, Der &d)
 Perform one Euler integration step: x += h*x'(x).

void integrate_rk2 (Sta &s, Real h, Real t, Der &d, Sta &s1)
 Perform one second-order Runge-Kutta integration step.

void integrate_rk4 (Sta &s, Real h, Real t, Der &d1, Der &d2, Der &d3, Der &d4, Sta &s1)
 Perform one fourth-order Runge-Kutta integration step.

void integrate_modmid (Sta &s, Real h, Real t, int n, Der &d1, Sta &s1, Sta &s2, Sta &s3)
 Perform one modified midpoint integration step.

void integrate_VVerlet (Sta &s, Real h, Der &d)

Protected Attributes

Vector dv1
 Auxiliary value to store a derivative (velocity).

Vector da1
 Auxiliary value to store a derivative (acceleration).

Der d1
 Auxiliary value to store a derivative (velocity,acceleration).

Vector dv2
 Auxiliary value to store a derivative (velocity).

Vector da2
 Auxiliary value to store a derivative (acceleration).

Der d2
 Auxiliary value to store a derivative (velocity,acceleration).

Vector dv3
 Auxiliary value to store a derivative (velocity).

Vector da3
 Auxiliary value to store a derivative (acceleration).

Der d3
 Auxiliary value to store a derivative (velocity,acceleration).

Vector dv4
 Auxiliary value to store a derivative (velocity).

Vector da4
 Auxiliary value to store a derivative (acceleration).

Der d4
 Auxiliary value to store a derivative (velocity,acceleration).

Positions sp1
 Auxiliary value to store a state (positions).

Vector sv1
 Auxiliary value to store a state (velocities).

Sta s1
 Auxiliary value to store a state (positions, velocities).

Positions sp2
 Auxiliary value to store a state (positions).

Vector sv2
 Auxiliary value to store a state (velocities).

Sta s2
 Auxiliary value to store a state (positions, velocities).

Positions sp3
 Auxiliary value to store a state (positions).

Vector sv3
 Auxiliary value to store a state (velocities).

Sta s3
 Auxiliary value to store a state (positions, velocities).

Vector velHalf


Member Typedef Documentation

template<class t_Positions, class t_Vector, class t_Real>
typedef Derivs<Vector,Vector> animal::OdeSolver< t_Positions, t_Vector, t_Real >::Der [protected]
 

Derivs.

Definition at line 68 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
typedef t_Positions animal::OdeSolver< t_Positions, t_Vector, t_Real >::Positions
 

Vector of positions.

Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >, animal::SimpleSolver< t_Positions, t_Vector, t_InvMasses, t_Real >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >, animal::OdeImplicitSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, value_type< X3DTK::MFFloat >::type >, animal::OdeImplicitSolver< Points, Velocities, value_type< SpringStiffness >::type >, and animal::OdeImplicitSolver< t_Positions, t_Vector, value_type< t_VecReal >::type >.

Definition at line 26 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
typedef t_Real animal::OdeSolver< t_Positions, t_Vector, t_Real >::Real
 

time and other scalar values

Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >, animal::SimpleSolver< t_Positions, t_Vector, t_InvMasses, t_Real >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >, animal::OdeImplicitSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, value_type< X3DTK::MFFloat >::type >, animal::OdeImplicitSolver< Points, Velocities, value_type< SpringStiffness >::type >, and animal::OdeImplicitSolver< t_Positions, t_Vector, value_type< t_VecReal >::type >.

Definition at line 28 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
typedef States<Positions,Vector> animal::OdeSolver< t_Positions, t_Vector, t_Real >::Sta [protected]
 

States.

Definition at line 87 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
typedef t_Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::Vector
 

velocities, forces, accelerations

Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >, animal::SimpleSolver< t_Positions, t_Vector, t_InvMasses, t_Real >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >, animal::OdeImplicitSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, value_type< X3DTK::MFFloat >::type >, animal::OdeImplicitSolver< Points, Velocities, value_type< SpringStiffness >::type >, and animal::OdeImplicitSolver< t_Positions, t_Vector, value_type< t_VecReal >::type >.

Definition at line 27 of file odeSolver.h.


Constructor & Destructor Documentation

template<class P, class V, class R>
animal::OdeSolver< P, V, R >::OdeSolver  ) 
 

Constructor.

Definition at line 14 of file odeSolver.inl.

template<class t_Positions, class t_Vector, class t_Real>
virtual animal::OdeSolver< t_Positions, t_Vector, t_Real >::~OdeSolver  )  [inline, virtual]
 

Destructor.

Definition at line 39 of file odeSolver.h.


Member Function Documentation

template<class t_Positions, class t_Vector, class t_Real>
virtual void animal::OdeSolver< t_Positions, t_Vector, t_Real >::computeAccelerations Vector acceleration,
const Positions pos,
const Vector vel
[protected, pure virtual]
 

Compute the acceleration for given positions and velocities.

Parameters:
acceleration the result
pos the positions
vel the velocities

Implemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >, animal::SimpleSolver< t_Positions, t_Vector, t_InvMasses, t_Real >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >, animal::OdeImplicitSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, value_type< X3DTK::MFFloat >::type >, animal::OdeImplicitSolver< Points, Velocities, value_type< SpringStiffness >::type >, and animal::OdeImplicitSolver< t_Positions, t_Vector, value_type< t_VecReal >::type >.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::computeDerivative(), and animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_VVerlet().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::computeDerivative Der d,
const Sta s,
Real 
[protected]
 

Compute the derivative(velocity, acceleration) of a given state(positions, velocities).

Definition at line 84 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::computeAccelerations(), and animal::v_eq().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::euler Positions pos,
Vector vel,
Real  h,
Real  t = 0
[protected]
 

Apply euler time integration.

Class States is a container with entry type equal or compatible with State.

Parameters:
pos positions
vel velocities
h time step
t current time

Definition at line 102 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_euler(), and animal::size().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::integrate_euler Sta s,
Real  h,
Real  t,
Der d
[protected]
 

Perform one Euler integration step: x += h*x'(x).

Parameters:
s initial and final state
h time step
t time at the beginning of the time step
d an auxiliary value to store the derivative Auxiliary value need not be initialized, but must have correct size.

Definition at line 184 of file odeSolver.inl.

References animal::v_eq_euler_step().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::euler().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::integrate_modmid Sta s,
Real  h,
Real  t,
int  n,
Der d1,
Sta s1,
Sta s2,
Sta s3
[protected]
 

Perform one modified midpoint integration step.

Parameters:
s initial and final state
h time step
t time at the beginning of the time step
n number of substeps
d1 an auxiliary value to store the derivative
s1 an auxiliary value to store a state
s2 an auxiliary value to store a state
s3 an auxiliary value to store a state Auxiliary values need not be initialized, but they must have correct size.

Definition at line 245 of file odeSolver.inl.

References animal::v_eq(), animal::v_eq_ab(), animal::v_eq_euler_step(), and animal::v_peq().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::modmid().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::integrate_rk2 Sta s,
Real  h,
Real  t,
Der d,
Sta s1
[protected]
 

Perform one second-order Runge-Kutta integration step.

Parameters:
s initial and final state
h time step
t time at the beginning of the time step
d an auxiliary value to store the derivative
s1 an auxiliary value to store a state Auxiliary values need not be initialized, but they must have correct size.

Definition at line 197 of file odeSolver.inl.

References animal::v_eq_euler_step().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::rk2().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::integrate_rk4 Sta s,
Real  h,
Real  t,
Der d1,
Der d2,
Der d3,
Der d4,
Sta s1
[protected]
 

Perform one fourth-order Runge-Kutta integration step.

Parameters:
s initial and final state
h time step
t time at the beginning of the time step
d1 an auxiliary value to store the derivative
d2 an auxiliary value to store the derivative
d3 an auxiliary value to store the derivative
d4 an auxiliary value to store the derivative
s1 an auxiliary value to store a state Auxiliary values need not be initialized, but they must have correct size.

Definition at line 214 of file odeSolver.inl.

References animal::v_eq_euler_step(), animal::v_peq(), and animal::v_teq().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::rk4().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::integrate_VVerlet Sta s,
Real  h,
Der d
[protected]
 

Definition at line 289 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::computeAccelerations(), animal::v_eq(), and animal::OdeSolver< t_Positions, t_Vector, t_Real >::velHalf.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::velocityVerlet().

template<class P, class V, class R>
int animal::OdeSolver< P, V, R >::method  )  const
 

Get the method.

Definition at line 66 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::_method.

template<class t_Positions, class t_Vector, class t_Real>
int animal::OdeSolver< t_Positions, t_Vector, t_Real >::MMIDsteps  )  const [inline]
 

The number of internal time steps in the modified midpoint method.

Definition at line 54 of file odeSolver.h.

Referenced by X3DTK::Qt::MassSpringEngineQt::MassSpringEngineQt().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::modmid Positions pos,
Vector vel,
Real  h,
int  n,
Real  t = 0
[protected]
 

Apply modified midpoint time integration.

Class States is a container with entry type equal or compatible with State.

Parameters:
pos positions
vel velocities
h time step
n number of substeps
t current time

Definition at line 152 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_modmid(), and animal::size().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::rk2 Positions pos,
Vector vel,
Real  h,
Real  t = 0
[protected]
 

Apply second-order Runge-Kutta time integration.

Class States is a container with entry type equal or compatible with State.

Parameters:
pos positions
vel velocities
h time step
t current time

Definition at line 118 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_rk2(), and animal::size().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::rk4 Positions pos,
Vector vel,
Real  h,
Real  t = 0
[protected]
 

Apply fourth-order Runge-Kutta time integration.

Class States is a container with entry type equal or compatible with State.

Parameters:
pos positions
vel velocities
h time step
t current time

Definition at line 133 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_rk4(), and animal::size().

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::setMethod int  m  ) 
 

Set the method.

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 60 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::_method.

template<class t_Positions, class t_Vector, class t_Real>
void animal::OdeSolver< t_Positions, t_Vector, t_Real >::setMMIDsteps int  n  )  [inline]
 

Set the number of internal time steps in the modified midpoint method.

Definition at line 51 of file odeSolver.h.

Referenced by X3DTK::Qt::MassSpringEngineQt::setMmidSubsteps().

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::solveODE Positions pos,
Vector vel,
Real  dt
[virtual]
 

Update positions and velocities by solving an ODE over time step dt.

Reimplemented in animal::MassSpringSolver< t_Positions, t_Vector, t_InvMasses, t_VecReal, t_VecIndex >, animal::OdeImplicitSolver< t_Positions, t_Vector, t_Real >, animal::SimpleSolver< t_Positions, t_Vector, t_InvMasses, t_Real >, animal::MassSpringSolver< Points, Velocities, Inv_Masses, SpringStiffness, IndexedSprings >, animal::MassSpringSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, X3DTK::MFFloat, X3DTK::MFFloat, X3DTK::MFInt32 >, animal::OdeImplicitSolver< X3DTK::MFVec3f, X3DTK::MFVec3f, value_type< X3DTK::MFFloat >::type >, animal::OdeImplicitSolver< Points, Velocities, value_type< SpringStiffness >::type >, and animal::OdeImplicitSolver< t_Positions, t_Vector, value_type< t_VecReal >::type >.

Definition at line 27 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::_method, animal::OdeSolver< t_Positions, t_Vector, t_Real >::_MMIDsteps, animal::OdeSolver< t_Positions, t_Vector, t_Real >::euler(), animal::OdeSolver< t_Positions, t_Vector, t_Real >::EULER, animal::OdeSolver< t_Positions, t_Vector, t_Real >::modmid(), animal::OdeSolver< t_Positions, t_Vector, t_Real >::MODMID, animal::OdeSolver< t_Positions, t_Vector, t_Real >::rk2(), animal::OdeSolver< t_Positions, t_Vector, t_Real >::RK2, animal::OdeSolver< t_Positions, t_Vector, t_Real >::rk4(), animal::OdeSolver< t_Positions, t_Vector, t_Real >::RK4, animal::OdeSolver< t_Positions, t_Vector, t_Real >::velocityVerlet(), and animal::OdeSolver< t_Positions, t_Vector, t_Real >::VVERLET.

template<class P, class V, class R>
void animal::OdeSolver< P, V, R >::velocityVerlet Positions pos,
Vector vel,
Real  h
[protected]
 

Apply velocity Verlet time integration.

Parameters:
pos positions
vel velocities
h time step
t current time

Definition at line 169 of file odeSolver.inl.

References animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_VVerlet(), animal::size(), and animal::OdeSolver< t_Positions, t_Vector, t_Real >::velHalf.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().


Member Data Documentation

template<class t_Positions, class t_Vector, class t_Real>
int animal::OdeSolver< t_Positions, t_Vector, t_Real >::_method
 

current integration method

Definition at line 56 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::method(), animal::OdeSolver< t_Positions, t_Vector, t_Real >::setMethod(), and animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class t_Positions, class t_Vector, class t_Real>
int animal::OdeSolver< t_Positions, t_Vector, t_Real >::_MMIDsteps
 

Number of substeps in MMID.

Definition at line 57 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class t_Positions, class t_Vector, class t_Real>
Der animal::OdeSolver< t_Positions, t_Vector, t_Real >::d1 [protected]
 

Auxiliary value to store a derivative (velocity,acceleration).

Definition at line 72 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Der animal::OdeSolver< t_Positions, t_Vector, t_Real >::d2 [protected]
 

Auxiliary value to store a derivative (velocity,acceleration).

Definition at line 76 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Der animal::OdeSolver< t_Positions, t_Vector, t_Real >::d3 [protected]
 

Auxiliary value to store a derivative (velocity,acceleration).

Definition at line 80 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Der animal::OdeSolver< t_Positions, t_Vector, t_Real >::d4 [protected]
 

Auxiliary value to store a derivative (velocity,acceleration).

Definition at line 84 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::da1 [protected]
 

Auxiliary value to store a derivative (acceleration).

Definition at line 71 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::da2 [protected]
 

Auxiliary value to store a derivative (acceleration).

Definition at line 75 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::da3 [protected]
 

Auxiliary value to store a derivative (acceleration).

Definition at line 79 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::da4 [protected]
 

Auxiliary value to store a derivative (acceleration).

Definition at line 83 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::dv1 [protected]
 

Auxiliary value to store a derivative (velocity).

Definition at line 70 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::dv2 [protected]
 

Auxiliary value to store a derivative (velocity).

Definition at line 74 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::dv3 [protected]
 

Auxiliary value to store a derivative (velocity).

Definition at line 78 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::dv4 [protected]
 

Auxiliary value to store a derivative (velocity).

Definition at line 82 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
const int animal::OdeSolver< t_Positions, t_Vector, t_Real >::EULER = 0 [static]
 

Symbol of explicit Euler integration.

Definition at line 29 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class t_Positions, class t_Vector, class t_Real>
const int animal::OdeSolver< t_Positions, t_Vector, t_Real >::MODMID = 3 [static]
 

Symbol of modified midpoint integration.

Definition at line 32 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class t_Positions, class t_Vector, class t_Real>
const int animal::OdeSolver< t_Positions, t_Vector, t_Real >::RK2 = 1 [static]
 

Symbol of second-order Runge-Kutta integration.

Definition at line 30 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class t_Positions, class t_Vector, class t_Real>
const int animal::OdeSolver< t_Positions, t_Vector, t_Real >::RK4 = 2 [static]
 

Symbol of fourth-order Runge-Kutta integration.

Definition at line 31 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().

template<class t_Positions, class t_Vector, class t_Real>
Sta animal::OdeSolver< t_Positions, t_Vector, t_Real >::s1 [protected]
 

Auxiliary value to store a state (positions, velocities).

Definition at line 91 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Sta animal::OdeSolver< t_Positions, t_Vector, t_Real >::s2 [protected]
 

Auxiliary value to store a state (positions, velocities).

Definition at line 95 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Sta animal::OdeSolver< t_Positions, t_Vector, t_Real >::s3 [protected]
 

Auxiliary value to store a state (positions, velocities).

Definition at line 99 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Positions animal::OdeSolver< t_Positions, t_Vector, t_Real >::sp1 [protected]
 

Auxiliary value to store a state (positions).

Definition at line 89 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Positions animal::OdeSolver< t_Positions, t_Vector, t_Real >::sp2 [protected]
 

Auxiliary value to store a state (positions).

Definition at line 93 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Positions animal::OdeSolver< t_Positions, t_Vector, t_Real >::sp3 [protected]
 

Auxiliary value to store a state (positions).

Definition at line 97 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::sv1 [protected]
 

Auxiliary value to store a state (velocities).

Definition at line 90 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::sv2 [protected]
 

Auxiliary value to store a state (velocities).

Definition at line 94 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::sv3 [protected]
 

Auxiliary value to store a state (velocities).

Definition at line 98 of file odeSolver.h.

template<class t_Positions, class t_Vector, class t_Real>
Vector animal::OdeSolver< t_Positions, t_Vector, t_Real >::velHalf [protected]
 

Definition at line 101 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::integrate_VVerlet(), and animal::OdeSolver< t_Positions, t_Vector, t_Real >::velocityVerlet().

template<class t_Positions, class t_Vector, class t_Real>
const int animal::OdeSolver< t_Positions, t_Vector, t_Real >::VVERLET = 4 [static]
 

Symbol of velocity Verlet integration.

Definition at line 33 of file odeSolver.h.

Referenced by animal::OdeSolver< t_Positions, t_Vector, t_Real >::solveODE().


The documentation for this class was generated from the following files:
Generated on Thu Dec 23 13:52:30 2004 by doxygen 1.3.6