Class SecondOrderDS#

Defined in Program listing for file kernel/src/modelingTools/SecondOrderDS.hpp

class SecondOrderDS : public DynamicalSystem#

Second Order non linear dynamical systems - \( M(q,z) \dot v = F(v, q, t, z) + p \).

This class defines and computes a generic ndof-dimensional second order Non Linear Dynamical System of the form :

\[\begin{split} M(q,z) \dot v = F(v, q, t, z) + p \\ \dot q = G(q,v) \end{split}\]

where

  • \( q \in R^{ndof} \) is the set of the coordinates,

  • \( \dot q =v \in R^{ndof} \) the velocity,

  • \( \ddot q = \dot v \in R^{ndof} \) the acceleration, i. e. the second time derivative of the generalized coordinates.

  • \( p \in R^{ndof} \) the reaction forces due to the Non Smooth Interaction.

  • \( M(q) \in R^{ndof \times ndof} \) is the inertia term (access : mass() method).

  • \( F( \dot q , q , t) \in R^{ndof} \) are the forces (access forces() method).

  • \( z \in R^{zSize} \) is a vector of arbitrary algebraic variables, some sort of discrete state.

q[i] is the derivative number i of q. Thus: q[0]= \( q \), global coordinates, q[1]= \( \dot q \), velocity, q[2]= \( \ddot q \), acceleration.

The following operators (and their jacobians) can be plugged, in the usual way (see User Guide, ‘User-defined plugins’)

If required (e.g. for Event-Driven like simulation), formulation as a first-order system is also available, and writes:

  • \( n= 2 ndof \)

  • \( x = \left[\begin{array}{c}q \\ \dot q\end{array}\right] \)

  • rhs given by:

    \[\begin{split} \dot x = \left[\begin{array}{c} \dot q\\ \ddot q = M^{-1}(q)\left[F(v, q , t, z) + p \right]\\ \end{array}\right] \end{split}\]

  • jacobian of the rhs, with respect to x

    \[\begin{split} \nabla_{x}rhs(x,t) = \left[\begin{array}{cc} 0 & I \\ \nabla_{q}(M^{-1}(q)F(v, q , t, z)) & \nabla_{\dot q}(M^{-1}(q)F(v, q , t, z)) \\ \end{array}\right] \end{split}\]
    with the input due to the non smooth law:

    \[\begin{split} \left[\begin{array}{c} 0 \\ p \end{array}\right] \end{split}\]

In that case, use the following methods:

Subclassed by LagrangianDS, NewtonEulerDS

Public Functions

inline virtual ~SecondOrderDS()#

destructor

inline SP::SiconosVector p(unsigned int level = 2) const#

get p

Parameters:

level – unsigned int, required level for p, default = 2

Returns:

pointer on a SiconosVector

inline SP::SiconosMatrix mass() const#

get mass matrix (pointer link)

Returns:

SP::SiconosMatrix

inline SP::SimpleMatrix inverseMass() const#

get (pointer) inverse or LU-factorization of the mass, used for LU-forward-backward computation

Returns:

pointer SP::SimpleMatrix

inline void setMassPtr(SP::SimpleMatrix newPtr)#

set mass to pointer newPtr

Parameters:

newPtr – a plugged matrix SP

inline virtual void setRhs(const SiconosVector &newValue) override#

set the value of the right-hand side, \( \dot x \)

Parameters:

newValueSiconosVector

inline virtual void setRhsPtr(SP::SiconosVector newPtr) override#

set right-hand side, \( \dot x \) (pointer link)

Parameters:

newPtr – SP::SiconosVector

virtual void computeForces(double time, SP::SiconosVector q, SP::SiconosVector velocity) = 0#

Compute \( F(v,q,t,z) \).

Parameters:
  • time – the current time

  • q – SP::SiconosVector: pointers on q

  • velocity – SP::SiconosVector: pointers on velocity

virtual void computeJacobianqForces(double time) = 0#

Compute \( \nabla_qF(v,q,t,z) \) for current \( q,v \) Default function to compute forces.

Parameters:

time – the current time

virtual void computeJacobianvForces(double time) = 0#

Compute \( \nabla_{\dot q}F(v,q,t,z) \) for current \( q,v \).

Parameters:

time – the current time

inline virtual unsigned int dimension() const override#

return the number of degrees of freedom of the system

Returns:

an unsigned int.

virtual SP::SiconosVector q() const = 0#

generalized coordinates of the system (vector of size dimension())

Returns:

pointer on a SiconosVector

virtual void setQ(const SiconosVector &newValue) = 0#

set value of generalized coordinates vector (copy)

Parameters:

newValue

virtual void setQPtr(SP::SiconosVector newPtr) = 0#

set value of generalized coordinates vector (pointer link)

Parameters:

newPtr

inline SP::SiconosVector q0() const#

get initial state (pointer link)

Returns:

pointer on a SiconosVector

virtual void setQ0(const SiconosVector &newValue) = 0#

set initial state (copy)

Parameters:

newValue

virtual void setQ0Ptr(SP::SiconosVector newPtr) = 0#

set initial state (pointer link)

Parameters:

newPtr

virtual SP::SiconosVector velocity() const = 0#

get velocity vector (pointer link)

Returns:

pointer on a SiconosVector

virtual void setVelocity(const SiconosVector &newValue) = 0#

set velocity vector (copy)

Parameters:

newValue

virtual void setVelocityPtr(SP::SiconosVector newPtr) = 0#

set velocity vector (pointer link)

Parameters:

newPtr

virtual SP::SiconosVector velocity0() const = 0#

get initial velocity (pointer)

Returns:

pointer on a SiconosVector

virtual void setVelocity0(const SiconosVector &newValue) = 0#

set initial velocity (copy)

Parameters:

newValue

virtual void setVelocity0Ptr(SP::SiconosVector newPtr) = 0#

set initial velocity (pointer link)

Parameters:

newPtr

virtual SP::SiconosVector acceleration() const = 0#

get acceleration (pointer link)

Returns:

pointer on a SiconosVector

virtual SP::SiconosVector forces() const = 0#

get \( F(v,q,t,z) \) (pointer link)

Returns:

pointer on a SiconosVector

virtual SP::SiconosMatrix jacobianqForces() const = 0#
Returns:

\( \nabla_qF(v,q,t,z) \) (pointer link)

virtual SP::SiconosMatrix jacobianvForces() const = 0#

get \( \nabla_{\dot q}F(v,q,t,z) \) (pointer link)

Returns:

pointer on a SiconosMatrix

virtual const SiconosMemory &qMemory() = 0#

get all the values of the state vector q stored in memory.

note: not const due to SchatzmanPaoliOSI::initializeWorkVectorsForDS

Returns:

a memory

virtual const SiconosMemory &velocityMemory() = 0#

get all the values of the state vector velocity stored in memory.

note: not const due to SchatzmanPaoliOSI::initializeWorkVectorsForDS

Returns:

a memory

virtual const SiconosMemory &forcesMemory() = 0#

get forces in memory buff

Returns:

pointer on a SiconosMemory

virtual void initMemory(unsigned int size) override = 0#

initialize the SiconosMemory objects with a positive size.

Parameters:

size – the size of the SiconosMemory. must be >= 0

virtual void computeMass() = 0#

default function to compute the mass

virtual void computeMass(SP::SiconosVector position) = 0#

function to compute the mass

Parameters:

position – value used to evaluate the mass matrix

virtual void setBoundaryConditions(SP::BoundaryCondition newbd)#

set Boundary Conditions

Parameters:

newbd – BoundaryConditions

inline SP::BoundaryCondition boundaryConditions()#

get Boundary Conditions

Returns:

SP::BoundaryCondition pointer on a BoundaryConditions

inline void setReactionToBoundaryConditions(SP::SiconosVector newrbd)#

set Reaction to Boundary Conditions

Parameters:

newrbd – BoundaryConditions pointer

inline SP::SiconosVector reactionToBoundaryConditions()#

get Reaction to Boundary Conditions

Returns:

pointer on a BoundaryConditions

virtual void init_inverse_mass() = 0#

Allocate memory for the lu factorization of the mass of the system.

Useful for some integrators with system inversion involving the mass

virtual void update_inverse_mass() = 0#

Update the content of the lu factorization of the mass of the system, if required.

virtual void init_forces() = 0#

Allocate memory for forces and its jacobian.