Class NewtonEulerDS#
Defined in Program listing for file kernel/src/modelingTools/NewtonEulerDS.hpp
-
class NewtonEulerDS : public SecondOrderDS#
NewtonEuler non linear dynamical systems.
The equations of motion in the Newton-Euler formalism can be stated as
\[\begin{split} \left\{\begin{array}{rcl} M \dot v + F_{int}(q,v, \Omega, t)&=& F_{ext}(t), \\ I \dot \Omega + \Omega \wedge I\Omega + M_{int}(q,v, \Omega, t) &=& M_{ext}(t), \\ \dot q &=& T(q) [ v, \Omega] \\ \dot R &=& R \tilde \Omega,\quad R^{-1}=R^T,\quad \det(R)=1 . \end{array}\right. \end{split}\]with
\( x_G,v_G \) position and velocity of the center of mass expressed in a inertial frame of reference (world frame)
\( \Omega \) angular velocity vector expressed in the body-fixed frame (frame attached to the object)
\( R \) rotation matrix form the inertial frame to the body-fixed frame \( R^{-1}=R^T, \det(R)=1 \), i.e \( R\in SO^+(3) \)
\( M=m\,I_{3\times 3} \) diagonal mass matrix with \( m \in \mathbb{R} \) the scalar mass
\( I \) constant inertia matrix
\( F_{ext} \) and \( M_{ext} \) are the external applied forces and moment
In the current implementation, \( R \) is parametrized by a unit quaternion.
Public Functions
-
NewtonEulerDS(SP::SiconosVector position, SP::SiconosVector twist, double mass, SP::SiconosMatrix inertia)#
constructor from a minimum set of data
- Parameters:
position – initial coordinates of this DynamicalSystem
twist – initial twist of this DynamicalSystem
mass – the mass
inertia – the inertia matrix
-
virtual ~NewtonEulerDS()#
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 SP::SiconosVector forces() const override#
get forces
- Returns:
pointer on a SiconosVector
-
inline virtual SP::SiconosMatrix jacobianqForces() const override#
get JacobianqForces
- Returns:
pointer on a SiconosMatrix
-
inline virtual SP::SiconosMatrix jacobianvForces() const override#
get JacobianvForces
- Returns:
pointer on a SiconosMatrix
-
inline virtual unsigned int getqDim() const#
Returns dimension of vector q.
-
inline virtual SP::SiconosVector q() const override#
get q (pointer link)
- 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 SP::SiconosVector twist() const#
get twist
- Returns:
pointer on a SiconosVector
-
inline virtual SP::SiconosVector velocity() const override#
get twist
- Returns:
pointer on a SiconosVector this accessor is left to get a uniform access to velocity. This should be removed with MechanicalDS class
-
inline virtual SP::SiconosVector velocity0() const override#
get initial velocity (pointer)
- Returns:
pointer on a SiconosVector
-
virtual void setVelocity(const SiconosVector &newValue) override#
set velocity (copy)
- Parameters:
newValue –
-
virtual void setVelocityPtr(SP::SiconosVector newPtr) override#
set velocity (pointer link)
- Parameters:
newPtr –
-
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 virtual void computeMass() override#
default function to compute the mass
-
inline virtual void computeMass(SP::SiconosVector position) override#
function to compute the mass
- Parameters:
position – value used to evaluate the mass matrix
-
SP::SiconosVector linearVelocity(bool absoluteRef) const#
Get the linear velocity in the absolute (inertial) or relative (body) frame of reference.
- Parameters:
absoluteRef – If true, velocity is returned in the inertial frame, otherwise velocity is returned in the body frame.
- Returns:
A SiconosVector of size 3 containing the linear velocity of this dynamical system.
-
void linearVelocity(bool absoluteRef, SiconosVector &v) const#
Fill a SiconosVector with the linear velocity in the absolute (inertial) or relative (body) frame of reference.
- Parameters:
absoluteRef – If true, velocity is returned in the inertial frame, otherwise velocity is returned in the body frame.
v – A SiconosVector of size 3 to receive the linear velocity.
-
SP::SiconosVector angularVelocity(bool absoluteRef) const#
Get the angular velocity in the absolute (inertial) or relative (body) frame of reference.
- Parameters:
absoluteRef – If true, velocity is returned in the inertial frame, otherwise velocity is returned in the body frame.
- Returns:
A SiconosVector of size 3 containing the angular velocity of this dynamical system.
-
void angularVelocity(bool absoluteRef, SiconosVector &w) const#
Fill a SiconosVector with the angular velocity in the absolute (inertial) or relative (body) frame of reference.
- Parameters:
absoluteRef – If true, velocity is returned in the inertial frame, otherwise velocity is returned in the body frame.
w – A SiconosVector of size 3 to receive the angular velocity.
-
inline double scalarMass() const#
get mass value
- Returns:
a double
-
inline void setScalarMass(double mass)#
Modify the scalar mass.
-
inline void setInertia(SP::SiconosMatrix newInertia)#
Modify the inertia matrix (pointer link)
- Parameters:
newInertia – the new inertia matrix
-
void setInertia(double ix, double iy, double iz)#
Modify the inertia matrix.
- Parameters:
ix – x component
iy – y component
iz – z component
-
void updateMassMatrix()#
to be called after scalar mass or inertia matrix have changed
-
inline void setFExtPtr(SP::SiconosVector newPtr)#
set fExt to pointer newPtr
- Parameters:
newPtr – a SP to a Simple vector
-
inline void setMExtPtr(SP::SiconosVector newPtr)#
set mExt to pointer newPtr
- Parameters:
newPtr – a SP to a Simple vector
-
inline virtual const SiconosMemory &qMemory() override#
get all the values of the state vector q stored in memory
- Returns:
a memory
-
inline const SiconosMemory &twistMemory()#
get all the values of the state vector twist stored in memory
- Returns:
a memory
-
inline virtual const SiconosMemory &velocityMemory() override#
get all the values of the state vector twist stored in memory
- Returns:
a memory
-
virtual void initMemory(unsigned int steps) override#
initialize the SiconosMemory objects with a positive size.
- Parameters:
steps – the size of the SiconosMemory (i)
-
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 virtual const SiconosMemory &forcesMemory() override#
get forces in memory buff
- Returns:
pointer on a SiconosMemory
-
double computeKineticEnergy()#
To compute the kinetic energy.
-
virtual void display(bool brief = true) const override#
print the data to the screen
-
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.
-
inline void setComputeFExtFunction(const std::string &pluginPath, const std::string &functionName)#
allow to set a specified function to compute _fExt
- Parameters:
pluginPath – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
inline void setComputeMExtFunction(const std::string &pluginPath, const std::string &functionName)#
allow to set a specified function to compute _mExt
- Parameters:
pluginPath – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
inline void setComputeFExtFunction(FExt_NE fct)#
set a specified function to compute _fExt
- Parameters:
fct – a pointer on the plugin function
-
inline void setComputeMExtFunction(FExt_NE fct)#
set a specified function to compute _mExt
- Parameters:
fct – a pointer on the plugin function
-
inline void setComputeFIntFunction(const std::string &pluginPath, const std::string &functionName)#
allow to set a specified function to compute _fInt
- Parameters:
pluginPath – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
inline void setComputeMIntFunction(const std::string &pluginPath, const std::string &functionName)#
allow to set a specified function to compute _mInt
- Parameters:
pluginPath – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
inline void setComputeFIntFunction(FInt_NE fct)#
set a specified function to compute _fInt
- Parameters:
fct – a pointer on the plugin function
-
inline void setComputeMIntFunction(FInt_NE fct)#
set a specified function to compute _mInt
- 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 setComputeJacobianFIntvFunction(const std::string &pluginPath, const std::string &functionName)#
allow to set a specified function to compute the jacobian following v of the internal forces w.r.t.
- Parameters:
pluginPath – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
void setComputeJacobianFIntqFunction(FInt_NE fct)#
set a specified function to compute jacobian following q of the FInt
- Parameters:
fct – a pointer on the plugin function
-
void setComputeJacobianFIntvFunction(FInt_NE fct)#
set a specified function to compute jacobian following v of the FInt
- Parameters:
fct – a pointer on the plugin function
-
void setComputeJacobianMIntqFunction(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 – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
void setComputeJacobianMIntvFunction(const std::string &pluginPath, const std::string &functionName)#
allow to set a specified function to compute the jacobian following v of the internal forces w.r.t.
- Parameters:
pluginPath – the complete path to the plugin
functionName – the name of the function to use in this plugin
-
void setComputeJacobianMIntqFunction(FInt_NE fct)#
set a specified function to compute jacobian following q of the FInt
- Parameters:
fct – a pointer on the plugin function
-
void setComputeJacobianMIntvFunction(FInt_NE fct)#
set a specified function to compute jacobian following v of the FInt
- Parameters:
fct – a pointer on the plugin function
-
virtual void computeFExt(double time)#
function to compute the external forces
- Parameters:
time – the current time
-
virtual void computeFExt(double time, SP::SiconosVector fExt)#
default function to compute the external forces
- Parameters:
time – the current time
fExt – the computed external force (in-out param)
-
virtual void computeMExt(double time, SP::SiconosVector mExt)#
function to compute the external moments The external moments are expressed by default in the body frame, since the Euler equation for Omega is written in the body—fixed frame.
Nevertheless, if _isMextExpressedInInertialFrame) is set to true, we assume that the external moment is given in the inertial frame and we perform the rotation afterwards
- Parameters:
time – the current time
mExt – the computed external moment (in-out param)
-
void addExtForceAtPos(SP::SiconosVector force, bool forceAbsRef, SP::SiconosVector pos = SP::SiconosVector(), bool posAbsRef = false)#
Adds a force/torque impulse to a body’s FExt and MExt vectors in either absolute (inertial) or relative (body) frame.
Modifies contents of _fExt and _mExt! Therefore these must have been set as constant vectors using setFExtPtr and setMExtPtr prior to calling this function. Adjustments to _mExt will take into account the value of _isMextExpressedInInertialFrame.
- Parameters:
force – A force vector to be added.
forceAbsRef – If true, force is in inertial frame, otherwise it is in body frame.
pos – A position at which force should be applied. If nullptr, the center of mass is assumed.
posAbsRef – If true, pos is in inertial frame, otherwise it is in body frame.
-
void computeFInt(double time, SP::SiconosVector q, SP::SiconosVector v)#
default function to compute the internal forces
- Parameters:
time – the current time function to compute the internal forces
time – the current time
q –
v –
-
virtual void computeFInt(double time, SP::SiconosVector q, SP::SiconosVector v, SP::SiconosVector fInt)#
default function to compute the internal forces
- Parameters:
time – the current time
q –
v –
fInt – the computed internal force vector
-
void computeMInt(double time, SP::SiconosVector q, SP::SiconosVector v)#
default function to compute the internal moments
- Parameters:
time – the current time
q –
v –
-
virtual void computeMInt(double time, SP::SiconosVector q, SP::SiconosVector v, SP::SiconosVector mInt)#
default function to compute the internal moments
- Parameters:
time – the current time
q –
v –
mInt – the computed internal moment vector
-
inline virtual void updatePlugins(double time) override#
default function to update the plugins functions using a new time:
- Parameters:
time – the current time
-
virtual void init_forces() override#
Allocate memory for forces and its jacobian.
-
virtual void computeForces(double time)#
Default function to compute forces.
- Parameters:
time – double, the current time
-
virtual void computeForces(double time, SP::SiconosVector q, SP::SiconosVector twist) override#
function to compute forces with some specific values for q and twist (ie not those of the current state).
- Parameters:
time – double : the current time
q – SP::SiconosVector: pointers on q
twist – SP::SiconosVector: pointers on twist
-
virtual void computeJacobianqForces(double time) override#
Default function to compute the jacobian w.r.t.
q of forces
- Parameters:
time – double, the current time
-
virtual void computeJacobianvForces(double time) override#
Default function to compute the jacobian w.r.t.
v of forces
- Parameters:
time – double, the current time
-
virtual void computeMGyr(SP::SiconosVector twist)#
function to compute gyroscopic forces with some specific values for q and twist (ie not those of the current state).
- Parameters:
twist – SP::SiconosVector: pointers on twist vector
-
virtual void computeMGyr(SP::SiconosVector twist, SP::SiconosVector mGyr)#
function to compute gyroscopic forces with some specific values for q and twist (ie not those of the current state).
- Parameters:
twist – pointer to twist vector
mGyr – pointer to gyroscopic forces
-
virtual void computeJacobianMGyrtwist(double time)#
Default function to compute the jacobian following q of mGyr.
- Parameters:
time – the current time
-
virtual void computeJacobianMGyrtwistByFD(double time, SP::SiconosVector q, SP::SiconosVector twist)#
Default function to compute the jacobian following q of mGyr by forward finite difference.
- Parameters:
time – the current time
q – current state
twist – pointer to twist vector
-
void computeJacobianFIntq(double time)#
To compute the jacobian w.r.t q of the internal forces.
- Parameters:
time – double : the current time
-
void computeJacobianFIntv(double time)#
To compute the jacobian w.r.t v of the internal forces.
- Parameters:
time – double : the current time
-
virtual void computeJacobianFIntq(double time, SP::SiconosVector position, SP::SiconosVector twist)#
To compute the jacobian w.r.t q of the internal forces.
- Parameters:
time – double
position – SP::SiconosVector
twist – SP::SiconosVector
-
void computeJacobianFIntqByFD(double time, SP::SiconosVector position, SP::SiconosVector twist)#
To compute the jacobian w.r.t q of the internal forces by forward finite difference.
- Parameters:
time – double
position – SP::SiconosVector
twist – SP::SiconosVector
-
virtual void computeJacobianFIntv(double time, SP::SiconosVector position, SP::SiconosVector twist)#
To compute the jacobian w.r.t.
v of the internal forces
- Parameters:
time – double: the current time
position – SP::SiconosVector
twist – SP::SiconosVector
-
void computeJacobianFIntvByFD(double time, SP::SiconosVector position, SP::SiconosVector twist)#
To compute the jacobian w.r.t v of the internal forces by forward finite difference.
- Parameters:
time – double
position – SP::SiconosVector
twist – SP::SiconosVector
-
virtual void computeJacobianMIntq(double time)#
To compute the jacobian w.r.t q of the internal forces.
- Parameters:
time – double : the current time
-
virtual void computeJacobianMIntv(double time)#
To compute the jacobian w.r.t v of the internal forces.
- Parameters:
time – double : the current time
-
virtual void computeJacobianMIntq(double time, SP::SiconosVector position, SP::SiconosVector twist)#
To compute the jacobian w.r.t q of the internal forces.
- Parameters:
time – double : the current time,
position – SP::SiconosVector
twist – SP::SiconosVector
-
void computeJacobianMIntqByFD(double time, SP::SiconosVector position, SP::SiconosVector twist)#
To compute the jacobian w.r.t q of the internal moments by forward finite difference.
- Parameters:
time – double
position – SP::SiconosVector
twist – SP::SiconosVector