Class LagrangianDS¶
Defined in Program listing for file kernel/src/modelingTools/LagrangianDS.hpp
-
class LagrangianDS : public SecondOrderDS¶
Lagrangian non linear dynamical systems - \( M(q,z) \dot v = F(v, q, t, z) + p \).
This class defines and computes a generic ndof-dimensional Lagrangian Non Linear Dynamical System of the form :
\[\begin{split} M(q,z) \dot v + F_{gyr}(v, q, z) + F_{int}(v , q , t, z) = F_{ext}(t, z) + p \\ \dot q = v \end{split}\]where
\( q \in R^{ndof} \) is the set of the generalized coordinates,
\( \dot q =v \in R^{ndof} \) the velocity, i. e. the time derivative of the generalized coordinates (Lagrangian systems).
\( \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_{gyr}(\dot q, q) \in R^{ndof} \) is the non linear inertia term (access fGyr() method).
\( F_{int}(\dot q , q , t) \in R^{ndof} \) are the internal forces (access fInt() method).
\( F_{ext}(t) \in R^{ndof} \) are the external forces (access fExt() method).
\( z \in R^{zSize} \) is a vector of arbitrary algebraic variables, some sort of discrete state.
The equation of motion is also shortly denoted as \( M(q,z) \dot v = F(v, q, t, z) + p \)
where \( F(v, q, t, z) \in R^{ndof} \) collects the total forces acting on the system, that is \( F(v, q, t, z) = F_{ext}(t, z) - F_{gyr}(v, q, z) + F_{int}(v, q , t, z) \)
This vector is saved and may be accessed using forces() method.
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’)
\( M(q) \) (computeMass())
\( F_{gyr}(v, q, z) \) (computeFGyr())
\( F_{int}(v , q , t, z) \) (computeFInt())
\( F_{ext}(t, z) \) (computeFExt())
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:
initRhs() to allocate/initialize memory for these new operators,
rhs() to get the rhs vector
computeRhs(), computeJacobianRhsx() …, to update the content of rhs, its jacobians …
Subclassed by LagrangianLinearDiagonalDS, LagrangianLinearTIDS
Public Functions
-
LagrangianDS(SP::SiconosVector position, SP::SiconosVector velocity)¶
constructor from initial state only, \( dv = p \)
- Parameters
position – SiconosVector : initial coordinates of this DynamicalSystem
velocity – SiconosVector : initial velocity of this DynamicalSystem
-
LagrangianDS(SP::SiconosVector position, SP::SiconosVector velocity, SP::SiconosMatrix mass)¶
constructor from initial state and mass, \( Mdv = p \)
- Parameters
position – SiconosVector : initial coordinates of this DynamicalSystem
velocity – SiconosVector : initial velocity of this DynamicalSystem
mass – SiconosMatrix : mass matrix
-
LagrangianDS(SP::SiconosVector position, SP::SiconosVector velocity, const std::string &plugin)¶
constructor from initial state and mass (plugin) \( Mdv = p \)
- Parameters
position – SiconosVector : initial coordinates of this DynamicalSystem
velocity – SiconosVector : initial velocity of this DynamicalSystem
plugin – std::string: plugin path to compute mass matrix
-
inline virtual ~LagrangianDS()¶
destructor
-
virtual void resetToInitialState() override¶
reset the state to the initial state
-
virtual void initRhs(double time) override¶
allocate (if needed) and compute rhs and its jacobian.
- Parameters
time – of initialization
-
virtual void initializeNonSmoothInput(unsigned int level) override¶
set nonsmooth input to zero
- Parameters
level – input-level to be initialized.
-
virtual void computeRhs(double time) override¶
update right-hand side for the current state
- Parameters
time – of interest
-
virtual void computeJacobianRhsx(double time) override¶
update \( \nabla_x rhs \) for the current state
- Parameters
time – of interest
-
virtual void resetAllNonSmoothParts() override¶
reset non-smooth part of the rhs (i.e.
p), for all ‘levels’
-
virtual void resetNonSmoothPart(unsigned int level) override¶
set nonsmooth part of the rhs (i.e.
p) to zero for a given level
- Parameters
level –
-
inline virtual void setRhs(const SiconosVector &newValue) override¶
set the value of the right-hand side, \( \dot x \)
- Parameters
newValue – SiconosVector
-
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) override¶
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) override¶
Compute \( \nabla_qF(v,q,t,z) \) for current \( q,v \) Default function to compute forces.
- Parameters
time – the current time
-
inline virtual void computeJacobianqDotForces(double time)¶
Compute \( \nabla_{\dot q}F(v,q,t,z) \) for current \( q,v \).
- Parameters
time – the current time
-
virtual void computeJacobianvForces(double time) override¶
Compute \( \nabla_{\dot q}F(v,q,t,z) \) for current \( q,v \).
- Parameters
time – the current time
-
inline virtual SP::SiconosVector q() const override¶
generalized coordinates of the system (vector of size dimension())
- Returns
pointer on a SiconosVector
-
virtual void setQ(const SiconosVector &newValue) override¶
set value of generalized coordinates vector (copy)
- Parameters
newValue –
-
virtual void setQPtr(SP::SiconosVector newPtr) override¶
set value of generalized coordinates vector (pointer link)
- Parameters
newPtr –
-
virtual void setQ0(const SiconosVector &newValue) override¶
set initial state (copy)
- Parameters
newValue –
-
virtual void setQ0Ptr(SP::SiconosVector newPtr) override¶
set initial state (pointer link)
- Parameters
newPtr –
-
inline virtual SP::SiconosVector velocity() const override¶
get velocity vector (pointer link)
- Returns
pointer on a SiconosVector
-
virtual void setVelocity(const SiconosVector &newValue) override¶
set velocity vector (copy)
- Parameters
newValue –
-
virtual void setVelocityPtr(SP::SiconosVector newPtr) override¶
set velocity vector (pointer link)
- Parameters
newPtr –
-
inline virtual SP::SiconosVector velocity0() const override¶
get initial velocity (pointer)
- Returns
pointer on a SiconosVector
-
virtual void setVelocity0(const SiconosVector &newValue) override¶
set initial velocity (copy)
- Parameters
newValue –
-
virtual void setVelocity0Ptr(SP::SiconosVector newPtr) override¶
set initial velocity (pointer link)
- Parameters
newPtr –
-
inline virtual SP::SiconosVector acceleration() const override¶
get acceleration (pointer link)
- Returns
pointer on a SiconosVector
-
inline SP::SiconosVector fInt() const¶
get $F_{int}$ (pointer link)
- Returns
pointer on a plugged vector
-
inline void setFIntPtr(SP::SiconosVector newPtr)¶
set $F_{int}$ (pointer link)
- Parameters
newPtr – a SP to plugged vector
-
inline SP::SiconosVector fExt() const¶
get \( F_{ext} \) , (pointer link)
- Returns
pointer on a plugged vector
-
inline void setFExtPtr(SP::SiconosVector newPtr)¶
set \( F_{ext} \) , (pointer link)
- Parameters
newPtr – a SP to a Simple vector
-
inline SP::SiconosVector fGyr() const¶
get \( F_{gyr} \) , (pointer link)
- Returns
pointer on a plugged vector
-
inline void setFGyrPtr(SP::SiconosVector newPtr)¶
set \( F_{gyr} \) , (pointer link)
- Parameters
newPtr – a SP to plugged vector
-
inline SP::SiconosMatrix jacobianFIntq() const¶
get \( \nabla_qF_{int} \) , (pointer link)
- Returns
pointer on a SiconosMatrix
-
inline SP::SiconosMatrix jacobianFIntqDot() const¶
get \( \nabla_{\dot q}F_{int} \) , (pointer link)
- Returns
pointer on a SiconosMatrix
-
inline void setJacobianFIntqPtr(SP::SiconosMatrix newPtr)¶
set \( \nabla_{q}F_{int} \) , (pointer link)
- Parameters
newPtr – a pointer to a SiconosMatrix
-
inline void setJacobianFIntqDotPtr(SP::SiconosMatrix newPtr)¶
set \( \nabla_{\dot q}F_{int} \) , (pointer link)
- Parameters
newPtr – a pointer to a SiconosMatrix
-
inline SP::SiconosMatrix jacobianFGyrq() const¶
get \( \nabla_{q}F_{gyr} \) , (pointer link)
- Returns
pointer on a SiconosMatrix
-
inline SP::SiconosMatrix jacobianFGyrqDot() const¶
get \( \nabla_{\dot q}F_{gyr} \) , (pointer link)
- Returns
pointer on a SiconosMatrix
-
inline void setJacobianFGyrqPtr(SP::SiconosMatrix newPtr)¶
get \( \nabla_{q}F_{gyr} \) , (pointer link)
- Parameters
newPtr – a SP SiconosMatrix
-
inline void setJacobianFGyrqDotPtr(SP::SiconosMatrix newPtr)¶
get \( \nabla_{\dot q}F_{gyr} \) , (pointer link)
- Parameters
newPtr – a SP SiconosMatrix
-
inline virtual SP::SiconosVector forces() const override¶
get \( F(v,q,t,z) \) (pointer link)
- Returns
pointer on a SiconosVector
-
inline virtual SP::SiconosMatrix jacobianqForces() const override¶
get \( \nabla_qF(v,q,t,z) \) (pointer link)
- Returns
pointer on a SiconosMatrix
-
inline virtual SP::SiconosMatrix jacobianvForces() const override¶
get \( \nabla_{\dot q}F(v,q,t,z) \) (pointer link)
- Returns
pointer on a SiconosMatrix
-
inline virtual const SiconosMemory &qMemory() override¶
get all the values of the state vector q stored in memory.
note: not const due to SchatzmanPaoliOSI::initializeWorkVectorsForDS
- Returns
a memory
-
inline virtual const SiconosMemory &velocityMemory() override¶
get all the values of the state vector velocity stored in memory.
note: not const due to SchatzmanPaoliOSI::initializeWorkVectorsForDS
- Returns
a memory
-
inline const SiconosMemory &pMemory(unsigned int level)¶
get all the values of the state vector p stored in memory
- Parameters
level –
- Returns
a memory
-
inline virtual const SiconosMemory &forcesMemory() override¶
get forces in memory buff
- Returns
pointer on a SiconosMemory
-
virtual void initMemory(unsigned int size) override¶
initialize the SiconosMemory objects with a positive size.
- Parameters
size – the size of the SiconosMemory. must be >= 0
-
virtual void swapInMemory() override¶
push the current values of x, q and r in the stored previous values xMemory, qMemory, rMemory,
- Todo:
Modify the function swapIn Memory with the new Object Memory
-
inline void setComputeMassFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute the mass
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
inline void setComputeMassFunction(FPtr7 fct)¶
set a specified function to compute Mass
- Parameters
fct – a pointer on the plugin function
-
void setComputeFIntFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute FInt
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
void setComputeFIntFunction(FPtr6 fct)¶
set a specified function to compute fInt
- Parameters
fct – a pointer on the plugin function
-
inline void setComputeFExtFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute Fext
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
inline void setComputeFExtFunction(VectorFunctionOfTime fct)¶
set a specified function to compute fExt
- Parameters
fct – a pointer on the plugin function
-
void setComputeFGyrFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute the inertia
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
void setComputeFGyrFunction(FPtr5 fct)¶
set a specified function to compute FGyr
- Parameters
fct – a pointer on the plugin function
-
void setComputeJacobianFIntqFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute the jacobian w.r.t q of the internal forces
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
void setComputeJacobianFIntqDotFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute the jacobian of the internal forces w.r.t.
q
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
void setComputeJacobianFIntqFunction(FPtr6 fct)¶
set a specified function to compute jacobian following q of the FInt
- Parameters
fct – a pointer on the plugin function
-
void setComputeJacobianFIntqDotFunction(FPtr6 fct)¶
set a specified function to compute jacobian following qDot of the FInt
- Parameters
fct – a pointer on the plugin function
-
void setComputeJacobianFGyrqFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute the jacobian w.r.t q of the the external forces
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
void setComputeJacobianFGyrqDotFunction(const std::string &pluginPath, const std::string &functionName)¶
allow to set a specified function to compute the jacobian w.r.t qDot of the the external strength
- Parameters
pluginPath – std::string : the complete path to the plugin
functionName – std::string : the name of the function to use in this plugin
-
void setComputeJacobianFGyrqFunction(FPtr5 fct)¶
set a specified function to compute the jacobian following q of FGyr
- Parameters
fct – a pointer on the plugin function
-
void setComputeJacobianFGyrqDotFunction(FPtr5 fct)¶
set a specified function to compute the jacobian following qDot of FGyr
- Parameters
fct – a pointer on the plugin function
-
virtual void computeMass() override¶
default function to compute the mass
-
virtual void computeMass(SP::SiconosVector position) override¶
function to compute the mass
- Parameters
position – value used to evaluate the mass matrix
-
virtual void computeFInt(double time)¶
default function to compute the internal strengths
- Parameters
time – the current time
-
virtual void computeFInt(double time, SP::SiconosVector position, SP::SiconosVector velocity)¶
function to compute the internal strengths with some specific values for position and velocity (ie not those of the current state).
- Parameters
time – the current time,
position – value used to evaluate the internal forces
velocity – value used to evaluate the internal forces
-
virtual void computeFExt(double time)¶
default function to compute the external strengths
- Parameters
time – the current time
-
virtual void computeFGyr()¶
default function to compute the inertia
-
virtual void computeFGyr(SP::SiconosVector position, SP::SiconosVector velocity)¶
function to compute the inertia with some specific values for q and velocity (ie not those of the current state).
- Parameters
position – value used to evaluate the inertia forces
velocity – value used to evaluate the inertia forces
-
virtual void computeJacobianFIntq(double time)¶
To compute the jacobian w.r.t q of the internal forces.
- Parameters
time – the current time
-
virtual void computeJacobianFIntqDot(double time)¶
To compute the jacobian w.r.t qDot of the internal forces.
- Parameters
time – the current time
-
virtual void computeJacobianFIntq(double time, SP::SiconosVector position, SP::SiconosVector velocity)¶
To compute the jacobian w.r.t q of the internal forces.
- Parameters
time – the current time
position – value used to evaluate the jacobian
velocity – value used to evaluate the jacobian
-
virtual void computeJacobianFIntqDot(double time, SP::SiconosVector position, SP::SiconosVector velocity)¶
To compute the jacobian w.r.t.
qDot of the internal forces
- Parameters
time – the current time
position – value used to evaluate the jacobian
velocity – value used to evaluate the jacobian
-
virtual void computeJacobianFGyrq()¶
function to compute the jacobian w.r.t.
q of the inertia forces
-
virtual void computeJacobianFGyrqDot()¶
function to compute the jacobian w.r.t.
qDot of the inertia forces
-
virtual void computeJacobianFGyrq(SP::SiconosVector position, SP::SiconosVector velocity)¶
function to compute the jacobian w.r.t.
q of the inertia forces
- Parameters
position – value used to evaluate the jacobian
velocity – value used to evaluate the jacobian
-
virtual void computeJacobianFGyrqDot(SP::SiconosVector position, SP::SiconosVector velocity)¶
function to compute the jacobian w.r.t.
qDot of the inertia forces
- Parameters
position – value used to evaluate the jacobian
velocity – value used to evaluate the jacobian
-
inline virtual void updatePlugins(double time) override¶
default function to update the plugins functions using a new time:
- Parameters
time – the current time
-
double computeKineticEnergy()¶
To compute the kinetic energy.
-
virtual void display(bool brief = true) const override¶
print the data of the dynamical system on the standard output
-
void computePostImpactVelocity()¶
Computes post-impact velocity, using pre-impact velocity and impulse (p) value.
Used in EventDriven (LsodarOSI->updateState)
-
void init_generalized_coordinates(unsigned int level)¶
Allocate memory for q[level], level > 1 Useful for some integrators that need q[2] or other coordinates vectors.
- Parameters
level – the required level
-
virtual void init_inverse_mass() override¶
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() override¶
Update the content of the lu factorization of the mass of the system, if required.
-
virtual void init_forces() override¶
Allocate memory for forces and its jacobian.