mbsim  4.0.0
MBSim Kernel
dynamic_system.h
1/* Copyright (C) 2004-2014 MBSim Development Team
2 * This library is free software; you can redistribute it and/or
3 * modify it under the terms of the GNU Lesser General Public
4 * License as published by the Free Software Foundation; either
5 * version 2.1 of the License, or (at your option) any later version.
6 *
7 * This library is distributed in the hope that it will be useful,
8 * but WITHOUT ANY WARRANTY; without even the implied warranty of
9 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
10 * Lesser General Public License for more details.
11 *
12 * You should have received a copy of the GNU Lesser General Public
13 * License along with this library; if not, write to the Free Software
14 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
15 *
16 * Contact: martin.o.foerg@googlemail.com
17 */
18
19#ifndef _DYNAMIC_SYSTEM_H_
20#define _DYNAMIC_SYSTEM_H_
21
22#include "mbsim/element.h"
23#include "mbsim/mbsim_event.h"
24
25namespace H5 {
26 class Group;
27}
28
29namespace MBSim {
30 class Frame;
31 class FixedRelativeFrame;
32 class Contour;
33 class RigidContour;
34 class Object;
35 class Link;
36 class Constraint;
37 class ModellingInterface;
38 class Contact;
39 class Observer;
40
58 class DynamicSystem : public Element {
59 public:
63 DynamicSystem(const std::string &name);
64
68 ~DynamicSystem() override;
69
70 /* INTERFACE FOR DERIVED CLASSES */
71 virtual void updateT();
72 virtual void updateh(int k=0);
73 virtual void updateM();
74 virtual void updateLLM();
75 virtual void updatedq();
76 virtual void updatedu();
77 virtual void updatedx();
78 virtual void updateqd();
79 virtual void updateud();
80 virtual void updatexd();
81 virtual void updatezd();
82 virtual void updatewb();
83 virtual void updateW(int j=0);
84 virtual void updateV(int j=0);
85 virtual void updater(int j=0);
86 virtual void updateJrla(int j=0);
87 virtual void updateg();
88 virtual void updategd();
89 virtual void updateStopVector();
90 virtual void updateStopVectorParameters();
91 virtual void updateLinkStatus();
92 virtual void updateLinkStatusReg();
93 virtual void updateWInverseKinetics();
94 virtual void updatebInverseKinetics();
95 virtual void sethSize(int hSize_, int j=0);
96 virtual int gethSize(int i=0) const { return hSize[i]; }
97 virtual int getqSize() const { return qSize; }
98 virtual int getuSize(int i=0) const { return uSize[i]; }
99 virtual void calcqSize();
100 virtual void calcuSize(int j=0);
101 virtual int getuInd(int i=0) { return uInd[i]; }
102 virtual void setqInd(int qInd_);
103 virtual void setuInd(int uInd_, int j=0);
104 virtual void sethInd(int hInd_, int j=0);
105 virtual void setxInd(int xInd_);
106 virtual const fmatvec::Vec& getq() const { return q; };
107 virtual fmatvec::Vec& getq() { return q; };
108 virtual const fmatvec::Vec& getu() const { return u; };
109 virtual fmatvec::Vec& getu() { return u; };
110 void setq(const fmatvec::Vec& q_) { q = q_; }
111 void setu(const fmatvec::Vec& u_) { u = u_; }
112 void setx(const fmatvec::Vec& x_) { x = x_; }
113 void setjsv(const fmatvec::VecInt& jsv_) { jsv = jsv_; }
114 void setInternalState(const fmatvec::Vec& internalState) { curis = internalState; nextis = internalState; }
115 virtual H5::GroupBase *getPlotGroup() { return plotGroup; }
116 std::shared_ptr<OpenMBV::Group> getOpenMBVGrp() override { return openMBVGrp; }
117 std::shared_ptr<OpenMBV::Group> getFramesOpenMBVGrp() override { return framesOpenMBVGrp; }
118 std::shared_ptr<OpenMBV::Group> getContoursOpenMBVGrp() override { return contoursOpenMBVGrp; }
119 std::shared_ptr<OpenMBV::Group> getGroupsOpenMBVGrp() override { return groupsOpenMBVGrp; }
120 std::shared_ptr<OpenMBV::Group> getObjectsOpenMBVGrp() override { return objectsOpenMBVGrp; }
121 std::shared_ptr<OpenMBV::Group> getLinksOpenMBVGrp() override { return linksOpenMBVGrp; }
122 std::shared_ptr<OpenMBV::Group> getConstraintsOpenMBVGrp() override { return constraintsOpenMBVGrp; }
123 std::shared_ptr<OpenMBV::Group> getObserversOpenMBVGrp() override { return observersOpenMBVGrp; }
124
125 virtual void calcxSize();
126 const fmatvec::Vec& getx() const { return x; };
127 fmatvec::Vec& getx() { return x; };
128 int getxSize() const { return xSize; }
129 void updatexRef(fmatvec::Vec &xParent);
130 void updatexdRef(fmatvec::Vec &xdParent);
131 void updatedxRef(fmatvec::Vec &dxParent);
132 void init(InitStage stage, const InitConfigSet &config) override;
133 virtual void initz();
134 virtual void writez(H5::GroupBase *parent);
135 virtual void readz0(H5::GroupBase *parent);
136 /*****************************************************/
137
138 /* INHERITED INTERFACE OF ELEMENT */
141 void plot() override;
142 void plotAtSpecialEvent() override;
143 /*****************************************************/
144
149 virtual int solveConstraintsFixpointSingle();
150
156 virtual int solveImpactsFixpointSingle();
157
162 virtual int solveConstraintsGaussSeidel();
163
169 virtual int solveImpactsGaussSeidel();
170
175 virtual int solveConstraintsRootFinding();
176
182 virtual int solveImpactsRootFinding();
183
187 virtual int jacobianConstraints();
188
192 virtual int jacobianImpacts();
193
197 virtual void checkConstraintsForTermination();
198
202 virtual void checkImpactsForTermination();
203
207 virtual void updaterFactors();
208
214 virtual Frame* getFrame(const std::string &name, bool check=true) const;
215
221 virtual Contour* getContour(const std::string &name, bool check=true) const;
222 /*****************************************************/
223
224 /* GETTER / SETTER */
225 void setFrameOfReference(Frame *frame) { R = frame; };
226 const Frame* getFrameOfReference() const { return R; };
227
228 const fmatvec::Mat& getT(bool check=true) const;
229 const fmatvec::Vec& geth(int i=0, bool check=true) const;
230 const fmatvec::SymMat& getM(bool check=true) const;
231 const fmatvec::SymMat& getLLM(bool check=true) const;
232 const fmatvec::Vec& getdq(bool check=true) const;
233 const fmatvec::Vec& getdu(bool check=true) const;
234 const fmatvec::Vec& getdx(bool check=true) const;
235 const fmatvec::Vec& getqd(bool check=true) const;
236 const fmatvec::Vec& getud(bool check=true) const;
237 const fmatvec::Vec& getxd(bool check=true) const;
238 const fmatvec::Mat& getW(int i=0, bool check=true) const;
239 const fmatvec::Mat& getV(int i=0, bool check=true) const;
240 const fmatvec::Vec& getwb(bool check=true) const;
241 const fmatvec::Vec& getla(bool check=true) const;
242 const fmatvec::Vec& getLa(bool check=true) const;
243 const fmatvec::Vec& getg(bool check=true) const;
244 const fmatvec::Vec& getgd(bool check=true) const;
245 fmatvec::Vec& getla(bool check=true);
246 fmatvec::Vec& getLa(bool check=true);
247 void setla(const fmatvec::Vec &la_) { la = la_; }
248 void setLa(const fmatvec::Vec &La_) { La = La_; }
249 fmatvec::VecInt& getjsv() { return jsv; }
250 const fmatvec::VecInt& getjsv() const { return jsv; }
251 fmatvec::Mat& getW(int i=0, bool check=true);
252 fmatvec::SymMat& getLLM(bool check=true);
253 fmatvec::Vec& getdq(bool check=true);
254 fmatvec::Vec& getdu(bool check=true);
255 fmatvec::Vec& getdx(bool check=true);
256
257 fmatvec::VecInt& getLinkStatus() { return LinkStatus; }
258 fmatvec::VecInt& getLinkStatusReg() { return LinkStatusReg; }
259 const fmatvec::VecInt& getLinkStatus() const { return LinkStatus; }
260 const fmatvec::VecInt& getLinkStatusReg() const { return LinkStatusReg; }
261
262 const fmatvec::Mat& evalT();
263 const fmatvec::Vec& evalh(int i=0);
264 const fmatvec::SymMat& evalM();
265 const fmatvec::SymMat& evalLLM();
266 const fmatvec::Mat& evalW(int i=0);
267 const fmatvec::Mat& evalV(int i=0);
268 const fmatvec::Vec& evalwb();
269 const fmatvec::Vec& evalr(int i=0);
270 const fmatvec::Mat& evalJrla(int i=0);
271 const fmatvec::Vec& evalrdt();
272 const fmatvec::Vec& evalg();
273 const fmatvec::Vec& evalgd();
274 const fmatvec::Mat& evalWInverseKinetics();
275 const fmatvec::Mat& evalbInverseKinetics();
276
277 void setqd(const fmatvec::Vec& qd_) { qd = qd_; }
278 void setud(const fmatvec::Vec& ud_) { ud = ud_; }
279 void setxd(const fmatvec::Vec& xd_) { xd = xd_; }
280
281 int getxInd() { return xInd; }
282 int getlaInd() const { return laInd; }
283
284 int gethInd(int i=0) { return hInd[i]; }
285 void setlaInd(int ind) { laInd = ind; }
286 void setisInd(int ind) { isInd = ind; }
287 void setgInd(int ind) { gInd = ind; }
288 void setgdInd(int ind) { gdInd = ind; }
289 void setrFactorInd(int ind) { rFactorInd = ind; }
290 virtual void setsvInd(int svInd_);
291 void setLinkStatusInd(int LinkStatusInd_) {LinkStatusInd = LinkStatusInd_;};
292 void setLinkStatusRegInd(int LinkStatusRegInd_) {LinkStatusRegInd = LinkStatusRegInd_;};
293
294 void setqSize(int qSize_) { qSize = qSize_; }
295 void setuSize(int uSize_, int i=0) { uSize[i] = uSize_; }
296 void setxSize(int xSize_) { xSize = xSize_; }
297
298 int getisSize() const { return isSize; }
299 int getlaSize() const { return laSize; }
300 int getgSize() const { return gSize; }
301 int getgdSize() const { return gdSize; }
302 int getrFactorSize() const { return rFactorSize; }
303 int getsvSize() const { return svSize; }
304 int getLinkStatusSize() const { return LinkStatusSize; }
305 int getLinkStatusRegSize() const { return LinkStatusRegSize; }
306 /*****************************************************/
307
308 const std::vector<Object*>& getObjects() const { return object; }
309 const std::vector<Link*>& getLinks() const { return link; }
310 const std::vector<DynamicSystem*>& getDynamicSystems() const { return dynamicsystem; }
311 const std::vector<Frame*>& getFrames() const { return frame; }
312 const std::vector<Contour*>& getContours() const { return contour; }
313 const std::vector<Link*>& getSetValuedLinks() const { return linkSetValued; }
314
319 void updateqRef(fmatvec::Vec &qParent);
320
325 void updateqdRef(fmatvec::Vec &qdParent);
326
327 void updatedqRef(fmatvec::Vec &dqParent);
328
333 void updateuRef(fmatvec::Vec &uParent);
334
339 void updateuallRef(fmatvec::Vec &uParent);
340
345 void updateudRef(fmatvec::Vec &udParent);
346
347 void updateduRef(fmatvec::Vec &duParent);
348
353 void updateudallRef(fmatvec::Vec &udParent);
354
362 void updatehRef(fmatvec::Vec &hParent, int j=0);
363
368 void updaterRef(fmatvec::Vec &rParent, int j=0);
369
370 void updateJrlaRef(fmatvec::Mat &rParent, int j=0);
371
376 void updaterdtRef(fmatvec::Vec &rdtParent);
377
382 void updateTRef(fmatvec::Mat &TParent);
383
389 void updateMRef(fmatvec::SymMat &MParent);
390
396 void updateLLMRef(fmatvec::SymMat &LLMParent);
397
401 virtual void updateInternalStateRef(fmatvec::Vec &curisParent, fmatvec::Vec &nextisParent);
402
407 virtual void updategRef(fmatvec::Vec &gParent);
408
413 virtual void updategdRef(fmatvec::Vec &gdParent);
414
419 void updatelaRef(fmatvec::Vec &laParent);
420
425 void updateLaRef(fmatvec::Vec &LaParent);
426
427 void updatelaInverseKineticsRef(fmatvec::Vec &laParent);
428 void updatebInverseKineticsRef(fmatvec::Mat &bParent);
429
434 virtual void updatewbRef(fmatvec::Vec &wbParent);
435
441 virtual void updateWRef(fmatvec::Mat &WParent, int j=0);
442
448 void updateWInverseKineticsRef(fmatvec::Mat &WParent);
449
455 virtual void updateVRef(fmatvec::Mat &VParent, int j=0);
456
461 void updatesvRef(fmatvec::Vec& svParent);
462
467 void updatejsvRef(fmatvec::VecInt &jsvParent);
468
473 void updateLinkStatusRef(fmatvec::VecInt &LinkStatusParent);
474
479 void updateLinkStatusRegRef(fmatvec::VecInt &LinkStatusRegParent);
480
485 void updateresRef(fmatvec::Vec &resParent);
486
491 void updaterFactorRef(fmatvec::Vec &rFactorParent);
492
493 void clearElementLists();
494
499 void buildListOfDynamicSystems(std::vector<DynamicSystem*> &sys);
500
505 void buildListOfObjects(std::vector<Object*> &obj);
506
511 void buildListOfLinks(std::vector<Link*> &lnk);
512
517 void buildListOfConstraints(std::vector<Constraint*> &crt);
518
523 void buildListOfFrames(std::vector<Frame*> &frm);
524
529 void buildListOfContours(std::vector<Contour*> &cnt);
530
535 void buildListOfModels(std::vector<ModellingInterface*> &modelList);
536
541 void buildListOfInverseKineticsLinks(std::vector<Link*> &iklnk);
542
547 void buildListOfObservers(std::vector<Observer*> &obsrv);
548
553
557 void setUpLinks();
558
562 bool gActiveChanged();
563
567 bool gActiveChangedReg();
568
572 bool detectImpact();
573
577 void calcsvSize();
578
582 void calclaSize(int j);
583
587 void calcLinkStatusSize();
588
593
598
603
604 void calcisSize();
605
611 virtual void aboutToUpdateInternalState();
612
613 virtual void postprocessing();
614
618 void calcgSize(int j);
619
630 void calcgdSize(int j);
631
635 void calcrFactorSize(int j);
636
641
645 void setUpActiveLinks();
646
661 void checkActive(int j);
662
666 virtual void setGeneralizedRelativePositionTolerance(double tol);
667
671 virtual void setGeneralizedRelativeVelocityTolerance(double tol);
672
676 virtual void setGeneralizedRelativeAccelerationTolerance(double tol);
677
681 virtual void setGeneralizedForceTolerance(double tol);
682
686 virtual void setGeneralizedImpulseTolerance(double tol);
687
690 virtual void setGeneralizedRelativePositionCorrectionValue(double corr);
691
695 virtual void setGeneralizedRelativeVelocityCorrectionValue(double corr);
696
700 void setrMax(double rMax);
701
706 int frameIndex(const Frame *frame_) const;
707
708 void addFrame(FixedRelativeFrame *frame);
709
710 void addContour(RigidContour *contour);
711
715 void addGroup(DynamicSystem *sys);
716
722 DynamicSystem* getGroup(const std::string &name,bool check=true) const;
723
727 void addObject(Object *obj);
728
734 Object* getObject(const std::string &name,bool check=true) const;
735
739 void addLink(Link *lnk);
740
744 void addConstraint(Constraint *crt);
745
749 void addInverseKineticsLink(Link *lnk);
750
751 Observer* getObserver(const std::string &name,bool check=true) const;
752 void addObserver(Observer *ele);
753
759 Link* getLink(const std::string &name,bool check=true) const;
760
766 Constraint* getConstraint(const std::string &name,bool check=true) const;
767
771 void addModel(ModellingInterface *model_);
772
778 ModellingInterface* getModel(const std::string &name, bool check=true) const;
779
781 Frame *getFrameI() { return I; }
782
783 Element *getChildByContainerAndName(const std::string &container, const std::string &name) const override;
784
785 virtual void updatecorr(int j);
786 void updatecorrRef(fmatvec::Vec &ref);
787 void calccorrSize(int j);
788
789 void checkRoot();
790
791 void resetUpToDate() override;
792
793 void updateStateTable();
794
795 H5::GroupBase *getFramesPlotGroup() override { return framesPlotGroup; }
796 H5::GroupBase *getContoursPlotGroup() override { return contoursPlotGroup; }
797 H5::GroupBase *getGroupsPlotGroup() override { return groupsPlotGroup; }
798 H5::GroupBase *getObjectsPlotGroup() override { return objectsPlotGroup; }
799 H5::GroupBase *getLinksPlotGroup() override { return linksPlotGroup; }
800 H5::GroupBase *getConstraintsPlotGroup() override { return constraintsPlotGroup; }
801 H5::GroupBase *getObserversPlotGroup() override { return observersPlotGroup; }
802
803 private:
804 friend class DynamicSystemSolver;
805 void addFrame(Frame *frame_);
806 void addContour(Contour *contour_);
807
808 protected:
813
817 std::vector<Object*> object;
818 std::vector<Object*> objectWithNonConstantMassMatrix;
819 std::vector<Link*> link;
820 std::vector<Link*> linkSingleValued;
821 std::vector<Link*> linkSetValued;
822 std::vector<Link*> linkSetValuedActive;
823 std::vector<Link*> linkWithStopVector;
824 std::vector<ModellingInterface*> model;
825 std::vector<DynamicSystem*> dynamicsystem;
826 std::vector<Link*> inverseKineticsLink;
827 std::vector<Observer*> observer;
828 std::vector<Constraint*> constraint;
829
833 fmatvec::Mat T;
834
838 fmatvec::SymMat M;
839
843 fmatvec::SymMat LLM;
844
848 fmatvec::Vec q, qd, dq;
849
853 fmatvec::Vec u, ud, du;
854
858 fmatvec::Vec x, xd, dx;
859
863 fmatvec::Vec h[2], r[2], rdt;
864
868 fmatvec::Mat Jrla[2];
869
873 fmatvec::Mat W[2], V[2];
874
878 fmatvec::Vec la, La;
879
880 fmatvec::Vec curis, nextis;
881
885 fmatvec::Vec g, gd;
886
890 fmatvec::Vec wb;
891
895 fmatvec::Vec res;
896
900 fmatvec::Vec rFactor;
901
905 fmatvec::Vec sv;
906
910 fmatvec::VecInt jsv;
911
915 fmatvec::VecInt LinkStatus;
916
920 fmatvec::VecInt LinkStatusReg;
921
925 int qSize, qInd;
926
930 int uSize[2], uInd[2];
931
935 int xSize, xInd;
936
940 int hSize[2], hInd[2];
941
942 int isSize, isInd;
943
947 int gSize, gInd;
948
952 int gdSize, gdInd;
953
957 int laSize, laInd;
958
962 int rFactorSize, rFactorInd;
963
967 int svSize, svInd;
968
972 int LinkStatusSize, LinkStatusInd;
973
977 int LinkStatusRegSize, LinkStatusRegInd;
978
982 std::vector<Frame*> frame;
983 std::vector<Contour*> contour;
984
985 std::shared_ptr<OpenMBV::Group> openMBVGrp;
986 std::shared_ptr<OpenMBV::Group> framesOpenMBVGrp;
987 std::shared_ptr<OpenMBV::Group> contoursOpenMBVGrp;
988 std::shared_ptr<OpenMBV::Group> groupsOpenMBVGrp;
989 std::shared_ptr<OpenMBV::Group> objectsOpenMBVGrp;
990 std::shared_ptr<OpenMBV::Group> linksOpenMBVGrp;
991 std::shared_ptr<OpenMBV::Group> constraintsOpenMBVGrp;
992 std::shared_ptr<OpenMBV::Group> observersOpenMBVGrp;
993 std::shared_ptr<H5::File> hdf5File;
994
997
1001 int laInverseKineticsSize, bInverseKineticsSize;
1002
1003 fmatvec::Mat WInverseKinetics, bInverseKinetics;
1004 fmatvec::Vec laInverseKinetics;
1005
1006 int corrSize, corrInd;
1007 fmatvec::Vec corr;
1008
1009 std::string saved_frameOfReference;
1010
1011 H5::GroupBase *framesPlotGroup;
1012 H5::GroupBase *contoursPlotGroup;
1013 H5::GroupBase *groupsPlotGroup;
1014 H5::GroupBase *objectsPlotGroup;
1015 H5::GroupBase *linksPlotGroup;
1016 H5::GroupBase *constraintsPlotGroup;
1017 H5::GroupBase *observersPlotGroup;
1018 };
1019}
1020
1021#endif
1022
Class for constraints between generalized coordinates of objects.
Definition: constraint.h:30
basic class for contour definition for rigid (which do not know about their shape) and flexible (they...
Definition: contour.h:40
solver interface for modelling and simulation of dynamic systems
Definition: dynamic_system_solver.h:61
dynamic system as topmost hierarchical level
Definition: dynamic_system.h:58
virtual void updatewbRef(fmatvec::Vec &wbParent)
references to TODO of dynamic system parent
Definition: dynamic_system.cc:854
DynamicSystem * getGroup(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1455
Frame * R
parent frame
Definition: dynamic_system.h:812
void updaterRef(fmatvec::Vec &rParent, int j=0)
references to nonsmooth right hand side of dynamic system parent
Definition: dynamic_system.cc:739
void init(InitStage stage, const InitConfigSet &config) override
plots time series header
Definition: dynamic_system.cc:407
fmatvec::SymMat M
mass matrix
Definition: dynamic_system.h:838
int svSize
size and local start index of stop vector relative to parent
Definition: dynamic_system.h:967
void buildListOfObservers(std::vector< Observer * > &obsrv)
build flat list of observers
Definition: dynamic_system.cc:1069
void calclaInverseKineticsSize()
calculates size of contact force parameters
Definition: dynamic_system.cc:1231
void updateWInverseKineticsRef(fmatvec::Mat &WParent)
references to contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:868
void updateqdRef(fmatvec::Vec &qdParent)
references to differentiated positions of dynamic system parent
Definition: dynamic_system.cc:621
DynamicSystem(const std::string &name)
constructor
Definition: dynamic_system.cc:45
ModellingInterface * getModel(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1551
virtual Frame * getFrame(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:583
void addObject(Object *obj)
Definition: dynamic_system.cc:1576
fmatvec::VecInt jsv
boolean evaluation of stop vector concerning roots
Definition: dynamic_system.h:910
std::vector< Object * > object
container for possible ingredients
Definition: dynamic_system.h:817
void updatehRef(fmatvec::Vec &hParent, int j=0)
references to smooth right hand side of dynamic system parent
Definition: dynamic_system.cc:726
void addConstraint(Constraint *crt)
Definition: dynamic_system.cc:1495
int uSize[2]
size and local start index of velocities relative to parent
Definition: dynamic_system.h:930
virtual void setGeneralizedForceTolerance(double tol)
Definition: dynamic_system.cc:1382
fmatvec::Vec sv
stop vector (root functions for event driven integration
Definition: dynamic_system.h:905
void updatesvRef(fmatvec::Vec &svParent)
references to stopvector (rootfunction for event driven integrator) of dynamic system parent
Definition: dynamic_system.cc:889
fmatvec::Vec u
velocities, differentiated velocities, initial velocities
Definition: dynamic_system.h:853
int frameIndex(const Frame *frame_) const
Definition: dynamic_system.cc:1447
void updateMRef(fmatvec::SymMat &MParent)
references to mass matrix of dynamic system parent
Definition: dynamic_system.cc:782
void calclaSize(int j)
calculates size of contact force parameters
Definition: dynamic_system.cc:1221
int laSize
size and local start index of contact force parameters relative to parent
Definition: dynamic_system.h:957
Frame * I
Definition: dynamic_system.h:996
virtual void setGeneralizedRelativePositionCorrectionValue(double corr)
Definition: dynamic_system.cc:1396
void plot() override
plots time dependent data
Definition: dynamic_system.cc:369
int hSize[2]
size and local start index of order smooth right hand side relative to parent
Definition: dynamic_system.h:940
fmatvec::VecInt LinkStatusReg
status of single-valued links
Definition: dynamic_system.h:920
fmatvec::Mat Jrla[2]
Jacobian dr/dla.
Definition: dynamic_system.h:868
void plotAtSpecialEvent() override
plots time dependent data at special events
Definition: dynamic_system.cc:388
void updatejsvRef(fmatvec::VecInt &jsvParent)
references to boolean evaluation of stopvector concerning roots of dynamic system parent
Definition: dynamic_system.cc:899
void updateLinkStatusRegRef(fmatvec::VecInt &LinkStatusRegParent)
references to status vector of single valued links
Definition: dynamic_system.cc:933
virtual void checkImpactsForTermination()
validate force laws concerning given tolerances on velocity level
Definition: dynamic_system.cc:571
void buildListOfContours(std::vector< Contour * > &cnt)
build flat list of contours
Definition: dynamic_system.cc:1042
std::vector< Frame * > frame
vector of frames and contours
Definition: dynamic_system.h:982
virtual void setGeneralizedImpulseTolerance(double tol)
Definition: dynamic_system.cc:1389
virtual void updateInternalStateRef(fmatvec::Vec &curisParent, fmatvec::Vec &nextisParent)
references to internal state of dynamic system parent
Definition: dynamic_system.cc:802
int qSize
size and local start index of positions relative to parent
Definition: dynamic_system.h:925
Frame * getFrameI()
Definition: dynamic_system.h:781
fmatvec::Vec g
relative distances and velocities
Definition: dynamic_system.h:885
fmatvec::Vec x
order one parameters, differentiated order one parameters, initial order one parameters
Definition: dynamic_system.h:858
void buildListOfInverseKineticsLinks(std::vector< Link * > &iklnk)
build flat list of inverse kinetics links
Definition: dynamic_system.cc:1060
void updateudallRef(fmatvec::Vec &udParent)
references to velocities of dynamic system parent
Definition: dynamic_system.cc:679
fmatvec::Vec la
contact force parameters
Definition: dynamic_system.h:878
void buildListOfLinks(std::vector< Link * > &lnk)
build flat list of links
Definition: dynamic_system.cc:1015
int LinkStatusSize
size and local start index of set-valued link status vector relative to parent
Definition: dynamic_system.h:972
virtual void updaterFactors()
update relaxation factors for contact equations
Definition: dynamic_system.cc:577
Constraint * getConstraint(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1526
virtual void updateVRef(fmatvec::Mat &VParent, int j=0)
references to condensed contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:882
void updatelaRef(fmatvec::Vec &laParent)
references to contact forces of dynamic system parent
Definition: dynamic_system.cc:833
virtual void setGeneralizedRelativeAccelerationTolerance(double tol)
Definition: dynamic_system.cc:1375
int rFactorSize
size and local start index of rfactors relative to parent
Definition: dynamic_system.h:962
~DynamicSystem() override
destructor
Definition: dynamic_system.cc:61
void calcbInverseKineticsSize()
calculates size of contact force parameters
Definition: dynamic_system.cc:1241
void updateLinkStatusRef(fmatvec::VecInt &LinkStatusParent)
references to status vector of set valued links with piecewise link equations (which piece is valid)
Definition: dynamic_system.cc:923
void updateresRef(fmatvec::Vec &resParent)
references to residuum of contact equations of dynamic system parent
Definition: dynamic_system.cc:909
void setUpObjectsWithNonConstantMassMatrix()
rearrange vector of active setvalued links
Definition: dynamic_system.cc:1337
void setUpInverseKinetics()
analyse constraints of dynamic systems for usage in inverse kinetics
Definition: dynamic_system.cc:1078
Element * getChildByContainerAndName(const std::string &container, const std::string &name) const override
Get the Element named name in the container named container.
Definition: dynamic_system.cc:1587
virtual void setGeneralizedRelativeVelocityCorrectionValue(double corr)
Definition: dynamic_system.cc:1403
bool detectImpact()
Definition: dynamic_system.cc:1130
virtual int jacobianConstraints()
compute JACOBIAN of contact equations
Definition: dynamic_system.cc:549
virtual int solveImpactsGaussSeidel()
solve impact equations with Gauss-Seidel scheme on velocity level
Definition: dynamic_system.cc:525
void buildListOfObjects(std::vector< Object * > &obj)
build flat list of objects
Definition: dynamic_system.cc:1006
virtual int solveConstraintsFixpointSingle()
solve contact equations with single step fixed point scheme
Definition: dynamic_system.cc:501
fmatvec::Vec h[2]
smooth, smooth with respect to objects, smooth with respect to links and nonsmooth
Definition: dynamic_system.h:863
void calcLinkStatusSize()
calculates size of set-valued link status vector
Definition: dynamic_system.cc:1178
virtual int jacobianImpacts()
compute JACOBIAN of contact equations on velocity level
Definition: dynamic_system.cc:557
void buildListOfModels(std::vector< ModellingInterface * > &modelList)
build flat list of models
Definition: dynamic_system.cc:1051
fmatvec::Vec q
positions, differentiated positions, initial positions
Definition: dynamic_system.h:848
fmatvec::Vec rFactor
rfactors for relaxation nonlinear contact equations
Definition: dynamic_system.h:900
void updateuRef(fmatvec::Vec &uParent)
references to velocities of dynamic system parent
Definition: dynamic_system.cc:641
fmatvec::Mat T
linear relation matrix of position and velocity parameters
Definition: dynamic_system.h:833
fmatvec::Vec res
residuum of nonlinear contact equations for Newton scheme
Definition: dynamic_system.h:895
virtual void updateWRef(fmatvec::Mat &WParent, int j=0)
references to contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:861
void setUpActiveLinks()
rearrange vector of active setvalued links
Definition: dynamic_system.cc:1345
void updaterdtRef(fmatvec::Vec &rdtParent)
references to nonsmooth right hand side of dynamic system parent
Definition: dynamic_system.cc:762
void setrMax(double rMax)
Definition: dynamic_system.cc:1410
int gSize
size and local start index of relative distances relative to parent
Definition: dynamic_system.h:947
virtual int solveImpactsRootFinding()
solve impact equations with Newton scheme on velocity level
Definition: dynamic_system.cc:541
void buildListOfDynamicSystems(std::vector< DynamicSystem * > &sys)
build flat list of dynamic systems
Definition: dynamic_system.cc:997
void updateTRef(fmatvec::Mat &TParent)
references to linear transformation matrix between differentiated positions and velocities of dynamic...
Definition: dynamic_system.cc:772
virtual Contour * getContour(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:597
int gdSize
size and local start index of relative velocities relative to parent
Definition: dynamic_system.h:952
void buildListOfConstraints(std::vector< Constraint * > &crt)
build flat list of all constraints
Definition: dynamic_system.cc:1024
void addGroup(DynamicSystem *sys)
Definition: dynamic_system.cc:1565
void addInverseKineticsLink(Link *lnk)
Definition: dynamic_system.cc:1507
bool gActiveChanged()
Definition: dynamic_system.cc:1110
virtual int solveImpactsFixpointSingle()
solve impact equations with single step fixed point scheme on velocity level
Definition: dynamic_system.cc:509
void calcgSize(int j)
calculates size of relative distances
Definition: dynamic_system.cc:1307
fmatvec::Vec wb
TODO.
Definition: dynamic_system.h:890
void updateudRef(fmatvec::Vec &udParent)
references to differentiated velocities of dynamic system parent
Definition: dynamic_system.cc:659
virtual int solveConstraintsRootFinding()
solve contact equations with Newton scheme
Definition: dynamic_system.cc:533
int xSize
size and local start index of order one parameters relative to parent
Definition: dynamic_system.h:935
virtual int solveConstraintsGaussSeidel()
solve contact equations with Gauss-Seidel scheme
Definition: dynamic_system.cc:517
int laInverseKineticsSize
size of contact force parameters of special links relative to parent
Definition: dynamic_system.h:1001
void updateqRef(fmatvec::Vec &qParent)
references to positions of dynamic system parent
Definition: dynamic_system.cc:611
virtual void setGeneralizedRelativeVelocityTolerance(double tol)
Definition: dynamic_system.cc:1368
fmatvec::VecInt LinkStatus
status of set-valued links
Definition: dynamic_system.h:915
void addLink(Link *lnk)
Definition: dynamic_system.cc:1483
Link * getLink(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1512
void calcgdSize(int j)
calculates size of gap velocities
Definition: dynamic_system.cc:1317
virtual void updategRef(fmatvec::Vec &gParent)
references to relative distances of dynamic system parent
Definition: dynamic_system.cc:819
void updateuallRef(fmatvec::Vec &uParent)
references to velocities of dynamic system parent
Definition: dynamic_system.cc:651
void addModel(ModellingInterface *model_)
Definition: dynamic_system.cc:1540
void setUpLinks()
distribute links to set- and single valued container
Definition: dynamic_system.cc:1089
fmatvec::SymMat LLM
Cholesky decomposition of mass matrix.
Definition: dynamic_system.h:843
Object * getObject(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1469
void updateLLMRef(fmatvec::SymMat &LLMParent)
references to Cholesky decomposition of dynamic system parent
Definition: dynamic_system.cc:792
void checkActive(int j)
check if set-valued contacts are active and set corresponding attributes
Definition: dynamic_system.cc:1355
void updateLaRef(fmatvec::Vec &LaParent)
references to contact impulses of dynamic system parent
Definition: dynamic_system.cc:840
virtual void setGeneralizedRelativePositionTolerance(double tol)
Definition: dynamic_system.cc:1361
void buildListOfFrames(std::vector< Frame * > &frm)
build flat list of frames
Definition: dynamic_system.cc:1033
virtual void aboutToUpdateInternalState()
Definition: dynamic_system.cc:1279
bool gActiveChangedReg()
Definition: dynamic_system.cc:1120
int LinkStatusRegSize
size and local start index of single-valued link status vector relative to parent
Definition: dynamic_system.h:977
void calcLinkStatusRegSize()
calculates size of single-valued link status vector
Definition: dynamic_system.cc:1193
void calcsvSize()
calculates size of stop vector
Definition: dynamic_system.cc:1208
void calcrFactorSize(int j)
calculates size of relaxation factors for contact equations
Definition: dynamic_system.cc:1327
virtual void updategdRef(fmatvec::Vec &gdParent)
references to relative velocities of dynamic system parent
Definition: dynamic_system.cc:826
void updaterFactorRef(fmatvec::Vec &rFactorParent)
references to relaxation factors for contact equations of dynamic system parent
Definition: dynamic_system.cc:916
void setDynamicSystemSolver(DynamicSystemSolver *sys) override
Definition: dynamic_system.cc:342
virtual void checkConstraintsForTermination()
validate force laws concerning given tolerances
Definition: dynamic_system.cc:565
basic class of MBSim mainly for plotting
Definition: element.h:56
H5::GroupBase * plotGroup
associated plot group
Definition: element.h:297
InitStage
The stages of the initialization.
Definition: element.h:62
std::string name
name of element
Definition: element.h:260
cartesian frame on rigid bodies
Definition: fixed_relative_frame.h:31
cartesian frame on bodies used for application of e.g. links and loads
Definition: frame.h:39
Interface for models of arbitrary domains, e.g. electrical components.
Definition: modelling_interface.h:35
class for all objects having own dynamics and mass
Definition: object.h:42
Definition: observer.h:26
basic class for rigid contours
Definition: rigid_contour.h:37
namespace MBSim
Definition: bilateral_constraint.cc:30