Project: | Siconos |
---|---|
Internal Release Number: | 1.0 |
Last update: | November 29, 2006 |
Related Documents: |
What for? | all classes/objects to define a simulation to solve the pre-defined NonSmoothDynamicalSystem. |
---|---|
Feature sets: | F-2001. and F-2.100, F-2.106; |
Depencies: | Siconos/Numerics,
Input-output (xml), Modeling |
Interfaces: | Siconos/Numerics |
Sources directory name: | Kernel/src/simulationTools |
So before starting any simulation process, you need to have a complete non smooth dynamical system (roughly some DynamicalSystems and some Interactions), attached to a Model. Then the key points of the Simulation are:
Above all these points, a global strategy must be defined. At the time, two possibilities exist in Siconos:
Finally to build a Simulation, four steps are required:
and so, the C++ code looks like:
// Prerequisite: // a Model*, say m, that holds a NonSmoothDynamicalSystem* named nsds. // // Step 1: choose the strategy, build the Simulation and link it to the Model TimeStepping * s = new TimeStepping(m); // or EventDriven * s = new EventDriven(m); // Step 2: the time discretisation, linked with the Simulation when built. TimeDiscretisation * td = new TimeDiscretisation(...,s) // Step 3: OneStepIntegrators; applied to DynamicalSystems that belong to nsds. Moreau * integrator = new Moreau(someDS,...,s) ; // or Lsodar * integrator = new Lsodar(someDS,s); // Step 4: OneStepNSProblems; applied to all the Interactions of the nsds. OneStepNSProblem * lcp = new LCP(s, ...)
Example:
Model * m = new Model(t0); // only initial time is provided Simulation * s = new TimeStepping(m); // s belongs to the Model m double h = ...; unsigned int nSteps = ...; TimeDiscretisation * td = new TimeDiscretisation(h, nSteps); // tk and final time will be automatically computed.
main functions:
Main variables: (a double) and
(a matrix of double). One value of each parameter is attached to each DynamicalSystem concerned by the OSI through a map (STL). For example:
// allDS is a set of DynamicalSystems : ds1, ds2 ... // s is a Simulation* Moreau * integrator = new Moreau(allDS, 0.5, s); // all values of theta are initialized to 0.5 double theta = integrator->getTheta(ds1); // returns the value of theta for the dynamical system ds1. SiconosMatrix * W = integrator->getW(ds1); // W is a SimpleMatrix *, that represents the W linked to ds1.
Construction:
// allDS is a set of DynamicalSystems and s a Simulation* Lsodar * integrator = new Lsodar(allDS,s);
For each of them, the required parameters are indicated between parenthesis. Their default values are:
Obviously, depending on the problem type, some solvers fit more than others.
where and
are the unknowns and
and
.
Thus, "formalizing" a LCP consists in writing M matrix and q vector (functions assembleM and computeQ ).
plus a Coulomb-Friction law.
classes: UnitaryRelation, UnitaryRelationsSet
The Simulation class has a specific member named indexSets, which is a vector (STL) of UnitaryRelationsSet.
In that case, if D is diagonal, the previous set is equivalent to n independent scalar equations, n being the size of the Interaction:
So a model with n Interactions of size one with the n relation will be equivalent.
But if D is not diagonal, all relations are to be treated all together, in one block. The same think occurs for a friction problem, where the minimum block size is 2 or 3 (for 2 or 3 dimensional friction) because we need to consider normal and tangential part(s) of y.
So to avoid confusion during Simulation process, we define UnitaryRelation object, which "cuts" the Interaction and its vectors into the smallest possible blocks. The size of these blocks corresponds to the non-smooth law size.
What is important is that during Simulation process we do not handles Interactions but UnitaryRelations. From the user point of view, that did not change anything and since Interactions and UnitaryRelations are linked you can goes on using Interactions in a transparent way. But if you intend to use or have a look into IndexSets, you should know what they handle, UnitaryRelations.
with UR=UnitaryRelation, and the i-1 derivative of the relation of ur.
Once again the proper definition of f and g depends on the system type, so for details see Event Driven algorithm for Lagrangian systems.
The possible types (derived classes) are:
The events manager is initialized with time-discretisation events, from the user time-discretisation. Then, each time a new event is detected (added by user or when a root is found during integration) it is scheduled in the manager.
The manager has also a processEvents function, which moves currentEvent to past events set, processes nextEvent and prepare the next step.
Other useful functions are:
EventDriven * s = new EventDriven(myModel); s->initialize(); // We get the events manager EventsManager * eventsManager = s->getEventsManagerPtr(); // while there are some events ... while(eventsManager->hasNextEvent()) { // integrate between current and next event s->advanceToEvent(); // solve the non-smooth dynamics, if necessary ... eventsManager->processEvents(); } // Or in one step: while(eventsManager->hasNextEvent()) { s->computeOneStep(); }
Because of the unilateral constraints, the evolution of the considered system may be non-smooth. Some jumps can occur in the velocity and the "acceleration" may not be defined everywhere. The generalized coordinates, assumed to be absolutely continuous are:
We will index with "+" and "-" right and left values of the variable at discontinuity.
The equations of motion are written in terms of a measure differential equation:
r being the generalized force due the unilateral constraints. Using the Lebesgue decomposition theorem and its variants, the differential measure dv and dr are decomposed in:
First term of the decomposition corresponds to the smooth part, with , the acceleration in the usual sense. The second term corresponds to the behavior at times of discontinuities, (
: Dirac), and the last term, a singular measure, will be neglected.
Thanks to these decompositions, the non-smooth Dynamics can be split into "impact equations", that will correspond to the non-smooth events, and some "smooth Dynamics". These equations are completed by the constraints, formulated at different kinematics levels, as shown in the following paragraphs.
is like an impulsion.
This equation will be solved at the time of impact together with an impact law. That is for a Newton impact law
This problem can be reduced on the local unknowns if the matrix
is assumed to be invertible, leading to the following Linear Complementarity Problem at time
of discontinuities of v:
Later this system will be identified as "LCP at velocity level".
where we assume that .
The following smooth systems are then to be solved:
To solve these systems, at each time, i.e. to known the configuration after each events and to integrate it numerically, it is useful to express the complementarity laws at different kinematics level. We also introduce the pre-defined index sets (about index sets, see Index Sets definition):
is the set of all the potential UnitaryRelations (UR).
(or if the UR is in
then contact occurs).
(or if the UR is in
, contact remains, no take off).
This results in the new writing of the "Bilateral Smooth Dynamics:
which can be reduced on variable and
, if M(q) is invertible, when
:
Later this system will be identified as "LCP at acceleration level".
knowing the value of and
at the beginning of the time step
:
This results in the computation of at this new point and to an update of the index sets
and
.
Next, term will be forget and considered as included in
.
And Lagrangian relations are (see Lagrangian Relation):
Q (resp. P) being a collection of all the q (resp. p) of the Dynamical Systems involved in the Interaction.
As we have seen in the previous section, the notion of kinematics level is really important. We introduce this in Siconos thanks to "[i]" notation. More precisely, for each Unitary Relation, we define y[i] as the derivative number i of variable y, according to time. In the same way, we denote the variable that is linked with y[i] through a Non-Smooth law (usually a complementarity). Finally to each
corresponds a p[i].
To make things clearer, let us rewrite the previous defined systems with Siconos notations:
with roots finding of:
Then, to build an EventDriven simulation, it is necessary to define two OneStepNSProblems, one at velocity and one at acceleration level. So here is a classical code for simulation construction:
EventDriven* s = new EventDriven(ball); // -- Time discretisation -- TimeDiscretisation * t = new TimeDiscretisation(timeStep,s); // -- OneStepIntegrators -- OneStepIntegrator * OSI = new Lsodar(setOfDS,s); // -- OneStepNsProblem -- OneStepNSProblem * impact = new LCP(s, "impact",solverName,101, 0.0001,"max",0.6); OneStepNSProblem * acceleration = new LCP(s, "acceleration",solverName,101, 0.0001,"max",0.6);
Finally, the algorithm described earlier is:
Corresponding code:
s->advanceToEvent() // This results in a call to Lsodar->integrate and to schedule of new non-smooth events if necessary
simulation->updateOutput(0, 1); simulation->updateIndexSets();
simulation->computeOneStepNSProblem("impact");
simulation->computeOneStepNSProblem("acceleration");
simulation->updateIndexSetsWithDoubleCondition();
simulation->nextStep();
The figure below represents the class architecture for TimeStepping simulation.
We will now describe discretisation process for two differents cases: first order linear systems and then Lagrangian systems, linear and time invariant.
with linear relations (LinearTIR):
The left-hand term is .
Right-hand terms are approximated this way :
is approximated using a
-method
since the second integral comes from independent sources, it can be evaluated with whatever quadrature method, for instance a -method
the third integral is approximated like in an implicit Euler integration
By replacing the accurate solution by the approximated value
, we get :
Assuming that is invertible, matrix
is defined as
. We get then :
An intermediate variable related to the smooth part of the system is defined as :
Thus the calculus of becomes :
with Lagrangian Linear Relations (LagrangianR, see Lagrangian Linear Relation):
Due to the non smooth character of the motion, the first term is integrated by an one order scheme( backward Euler-like) such that :
For simplicity sake, we note the approximation of q and :
For the other terms, a -method is used :
For the term which contains the reaction force, we state a new variable such that :
The displacement is integrated through the velocity with :
Finally, we obtain the time discretized equation of motion:
which can be written :
where
The free velocity correponds to the velocity of the system without any constraints.
In the Moreau's time--stepping, we use a reformulation of the unilateral constraints in terms of velocity:
which leads to the following discretisation :
where is a prediction of the position at time
, for instance,
.
If we want to introduce now the Newton impact law, we consider an equivalent velocity defined by
and we apply the constraints directly on this velocity :
This set of equations can be reduced to a "condensed" system in terms of and
:
is the set of all the potential UnitaryRelations (UR).
For second order systems:
.
Thus, the LCP is built only for unitary relations that belongs to , level=0 for first order and level=1 for second order systems.
Then, the steps of a Moreau's Time-Stepping simulation will be:
Knowing all values at the beginning of the time step ,
TimeStepping * s = new TimeStepping(myModel); TimeDiscretisation * t = new TimeDiscretisation(timeStep,s); s->initialize(); int k = t->getK(); // Current step int N = t->getNSteps(); // Number of time steps // --- Time loop --- while(k < N)// for each time step ... { // transfer of state i+1 into state i and time incrementation s->nextStep(); // get current time step k = t->getK(); // compute xFree, or qFree,vFree s->computeFreeStep(); // Update the index sets (ie compute prediction yp for Lagrangian Systems) s->updateIndexSets(); // Formalize and solve a LCP computeOneStepNSProblem("timeStepping"); // Update state, using last computed values s->update(level); // }
Note that all time-independent operators are computed during simulation initialisation.