All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Pages
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 
25 namespace H5 {
26  class Group;
27 }
28 
29 namespace MBSim {
30  class Frame;
31  class FixedRelativeFrame;
32  class Contour;
33  class Object;
34  class Link;
35  class ModellingInterface;
36  class Contact;
37  class InverseKineticsJoint;
38  class Observer;
39 
57  class DynamicSystem : public Element {
58  public:
62  DynamicSystem(const std::string &name);
63 
67  virtual ~DynamicSystem();
68 
69  /* INTERFACE FOR DERIVED CLASSES */
70  virtual void updateT(double t);
71  virtual void updateh(double t, int i=0);
72  virtual void updateStateDependentVariables(double t) = 0;
73  virtual void updateStateDerivativeDependentVariables(double t);
74  virtual void updateM(double t, int i=0);
75  virtual void updateJacobians(double t, int j=0) = 0;
76  virtual void updatedq(double t, double dt);
77  virtual void updateud(double t, int i=0) { THROW_MBSIMERROR("(DynamicSystem::updateud): Not implemented!"); }
78  virtual void updatezd(double t) = 0;
79  virtual void updatedu(double t, double dt) = 0;
80  virtual void updateqd(double t);
81  virtual void sethSize(int hSize_, int i=0);
82  virtual int gethSize(int i=0) const { return hSize[i]; }
83  virtual int getqSize() const { return qSize; }
84  virtual int getuSize(int i=0) const { return uSize[i]; }
85  virtual void calcqSize();
86  virtual void calcuSize(int j=0);
87  //virtual int getqInd(DynamicSystem* sys);
88  virtual int getuInd(int i=0) { return uInd[i]; }
89  //virtual int getuInd(DynamicSystem* sys, int i=0);
90  //virtual void setqInd(int qInd_) { qInd = qInd_; }
91  //virtual void setuInd(int uInd_, int i=0) { uInd[i] = uInd_; }
92  virtual void setqInd(int qInd_);
93  virtual void setuInd(int uInd_, int i=0);
94  virtual void sethInd(int hInd_, int i=0);
95  virtual void setxInd(int xInd_);
96  //virtual int gethInd(DynamicSystem* sys, int i=0);
97  virtual const fmatvec::Vec& getq() const { return q; };
98  virtual fmatvec::Vec& getq() { return q; };
99  virtual const fmatvec::Vec& getu() const { return u; };
100  virtual fmatvec::Vec& getu() { return u; };
101  virtual H5::GroupBase *getPlotGroup() { return plotGroup; }
102  virtual PlotFeatureStatus getPlotFeature(PlotFeature fp) { return Element::getPlotFeature(fp); };
103  virtual PlotFeatureStatus getPlotFeatureForChildren(PlotFeature fp) { return Element::getPlotFeatureForChildren(fp); };
104 #ifdef HAVE_OPENMBVCPPINTERFACE
105  virtual boost::shared_ptr<OpenMBV::Group> getOpenMBVGrp();
106 #endif
107 
108  virtual void updatewb(double t, int j=0);
109  virtual void updateW(double t, int j=0);
110  virtual void updateV(double t, int j=0);
111  virtual void updateg(double t);
112  virtual void updategd(double t);
113  virtual void updateStopVector(double t);
114  virtual void updateLinkStatus(double t);
115  virtual void updateLinkStatusReg(double t);
116 
117  virtual void updategInverseKinetics(double t);
118  virtual void updategdInverseKinetics(double t);
119  virtual void updateWInverseKinetics(double t, int j=0);
120  virtual void updatehInverseKinetics(double t, int j=0);
121  virtual void updateJacobiansInverseKinetics(double t, int j=0);
122  virtual void updatebInverseKinetics(double t);
123 
124  virtual void updatedx(double t, double dt);
125  virtual void updatexd(double t);
126  virtual void calcxSize();
127  const fmatvec::Vec& getx() const { return x; };
128  fmatvec::Vec& getx() { return x; };
129  int getxSize() const { return xSize; }
130  void updatexRef(const fmatvec::Vec &ref);
131  void updatexdRef(const fmatvec::Vec &ref);
132  virtual void init(InitStage stage);
133  virtual void initz();
134  virtual void writez(H5::GroupBase *group);
135  virtual void readz0(H5::GroupBase *parent);
136  /*****************************************************/
137 
138  /* INHERITED INTERFACE OF ELEMENT */
140  virtual std::string getType() const { return "DynamicSystem"; }
141  virtual void setDynamicSystemSolver(DynamicSystemSolver* sys);
142  virtual void plot(double t, double dt);
143  virtual void plotAtSpecialEvent(double t, double dt=1.);
144  virtual void closePlot();
145  /*****************************************************/
146 
147  /* INTERFACE FOR DERIVED CLASSES */
151  virtual void facLLM(int i=0) = 0;
152 
157  virtual int solveConstraintsFixpointSingle();
158 
164  virtual int solveImpactsFixpointSingle(double dt);
165 
170  virtual int solveConstraintsGaussSeidel();
171 
177  virtual int solveImpactsGaussSeidel(double dt);
178 
183  virtual int solveConstraintsRootFinding();
184 
190  virtual int solveImpactsRootFinding(double dt);
191 
195  virtual int jacobianConstraints();
196 
200  virtual int jacobianImpacts();
201 
205  virtual void checkConstraintsForTermination();
206 
210  virtual void checkImpactsForTermination(double dt);
211 
215  virtual void updaterFactors();
216 
222  virtual Frame* getFrame(const std::string &name, bool check=true) const;
223 
229  virtual Contour* getContour(const std::string &name, bool check=true) const;
230  /*****************************************************/
231 
232  /* GETTER / SETTER */
233  void setPosition(const fmatvec::Vec3& PrPF_) { PrPF = PrPF_; }
234  void setOrientation(const fmatvec::SqrMat3& APF_) { APF = APF_; }
235  void setFrameOfReference(Frame *frame) { R = frame; };
236  const fmatvec::Vec3& getPosition() const { return PrPF; }
237  const fmatvec::SqrMat3& getOrientation() const { return APF; }
238  const Frame* getFrameOfReference() const { return R; };
239 
240  const fmatvec::Vec& getxd() const { return xd; };
241  fmatvec::Vec& getxd() { return xd; };
242  const fmatvec::Vec& getx0() const { return x0; };
243  fmatvec::Vec& getx0() { return x0; };
244 
245  const fmatvec::Mat& getT() const { return T; };
246  const fmatvec::SymMat& getM(int i=0) const { return M[i]; };
247  const fmatvec::SymMat& getLLM(int i=0) const { return LLM[i]; };
248  fmatvec::SymMat& getLLM(int i=0) { return LLM[i]; };
249  const fmatvec::Vec& geth(int j=0) const { return h[j]; };
250  fmatvec::Vec& geth(int j=0) { return h[j]; };
251  const fmatvec::Vec& getf() const { return f; };
252  fmatvec::Vec& getf() { return f; };
253 
254  const fmatvec::Mat& getW(int i=0) const { return W[i]; }
255  fmatvec::Mat& getW(int i=0) { return W[i]; }
256  const fmatvec::Mat& getV(int i=0) const { return V[i]; }
257  fmatvec::Mat& getV(int i=0) { return V[i]; }
258 
259  const fmatvec::Vec& getla() const { return la; }
260  fmatvec::Vec& getla() { return la; }
261  const fmatvec::Vec& getg() const { return g; }
262  fmatvec::Vec& getg() { return g; }
263  const fmatvec::Vec& getgd() const { return gd; }
264  fmatvec::Vec& getgd() { return gd; }
265  const fmatvec::Vec& getrFactor() const { return rFactor; }
266  fmatvec::Vec& getrFactor() { return rFactor; }
267  fmatvec::Vec& getsv() { return sv; }
268  const fmatvec::Vec& getsv() const { return sv; }
269  fmatvec::VecInt& getjsv() { return jsv; }
270  const fmatvec::VecInt& getjsv() const { return jsv; }
271  fmatvec::VecInt& getLinkStatus() { return LinkStatus; }
272  fmatvec::VecInt& getLinkStatusReg() { return LinkStatusReg; }
273  const fmatvec::VecInt& getLinkStatus() const { return LinkStatus; }
274  const fmatvec::VecInt& getLinkStatusReg() const { return LinkStatusReg; }
275  const fmatvec::Vec& getres() const { return res; }
276  fmatvec::Vec& getres() { return res; }
277  const fmatvec::Vec& getcorr() const { return corr; };
278  fmatvec::Vec& getcorr() { return corr; };
279 
280  void setx(const fmatvec::Vec& x_) { x = x_; }
281  void setx0(const fmatvec::Vec &x0_) { x0 = x0_; }
282  void setx0(double x0_) { x0 = fmatvec::Vec(1,fmatvec::INIT,x0_); }
283 
284  int getxInd() { return xInd; }
285  int getlaInd() const { return laInd; }
286 
287  int gethInd(int i=0) { return hInd[i]; }
288  void setlaInd(int ind) { laInd = ind; }
289  void setgInd(int ind) { gInd = ind; }
290  void setgdInd(int ind) { gdInd = ind; }
291  void setrFactorInd(int ind) { rFactorInd = ind; }
292  virtual void setsvInd(int svInd_);
293  void setLinkStatusInd(int LinkStatusInd_) {LinkStatusInd = LinkStatusInd_;};
294  void setLinkStatusRegInd(int LinkStatusRegInd_) {LinkStatusRegInd = LinkStatusRegInd_;};
295 
296  int getzSize() const { return qSize + uSize[0] + xSize; }
297 
298  void setqSize(int qSize_) { qSize = qSize_; }
299  void setuSize(int uSize_, int i=0) { uSize[i] = uSize_; }
300  void setxSize(int xSize_) { xSize = xSize_; }
301 
302  int getlaSize() const { return laSize; }
303  int getgSize() const { return gSize; }
304  int getgdSize() const { return gdSize; }
305  int getrFactorSize() const { return rFactorSize; }
306  int getsvSize() const { return svSize; }
307  int getLinkStatusSize() const { return LinkStatusSize; }
308  int getLinkStatusRegSize() const { return LinkStatusRegSize; }
309  /*****************************************************/
310 
311  const std::vector<Object*>& getObjects() const { return object; }
312  const std::vector<Link*>& getLinks() const { return link; }
313  const std::vector<DynamicSystem*>& getDynamicSystems() const { return dynamicsystem; }
314  const std::vector<Frame*>& getFrames() const { return frame; }
315  const std::vector<Contour*>& getContours() const { return contour; }
316 
321  void updateqRef(const fmatvec::Vec &ref);
322 
327  void updateqdRef(const fmatvec::Vec &ref);
328 
333  void updateuRef(const fmatvec::Vec &ref);
334 
339  void updateuallRef(const fmatvec::Vec &ref);
340 
345  void updateudRef(const fmatvec::Vec &ref, int i=0);
346 
351  void updateudallRef(const fmatvec::Vec &ref, int i=0);
352 
360  void updatehRef(const fmatvec::Vec &hRef, int i=0);
361 
366  void updaterRef(const fmatvec::Vec &ref, int j=0);
367 
372  void updateTRef(const fmatvec::Mat &ref);
373 
379  void updateMRef(const fmatvec::SymMat &ref, int i=0);
380 
386  void updateLLMRef(const fmatvec::SymMat &ref, int i=0);
387 
392  void updategRef(const fmatvec::Vec &ref);
393 
398  void updategdRef(const fmatvec::Vec &ref);
399 
404  void updatelaRef(const fmatvec::Vec &ref);
405 
406  void updatelaInverseKineticsRef(const fmatvec::Vec &ref);
407  void updatebInverseKineticsRef(const fmatvec::Mat &ref);
408 
413  void updatewbRef(const fmatvec::Vec &ref);
414 
420  void updateWRef(const fmatvec::Mat &ref, int i=0);
421 
427  void updateWInverseKineticsRef(const fmatvec::Mat &ref, int i=0);
428 
434  void updateVRef(const fmatvec::Mat &ref, int i=0);
435 
440  void updatesvRef(const fmatvec::Vec& ref);
441 
446  void updatejsvRef(const fmatvec::VecInt &ref);
447 
452  void updateLinkStatusRef(const fmatvec::VecInt &LinkStatusParent);
453 
458  void updateLinkStatusRegRef(const fmatvec::VecInt &LinkStatusRegParent);
459 
464  void updateresRef(const fmatvec::Vec &ref);
465 
470  void updaterFactorRef(const fmatvec::Vec &ref);
471 
472  void clearElementLists();
473 
478  void buildListOfDynamicSystems(std::vector<DynamicSystem*> &sys);
479 
484  void buildListOfObjects(std::vector<Object*> &obj);
485 
490  void buildListOfLinks(std::vector<Link*> &lnk);
491 
496  void buildListOfSetValuedLinks(std::vector<Link*> &lnk);
497 
502  void buildListOfFrames(std::vector<Frame*> &frm);
503 
508  void buildListOfContours(std::vector<Contour*> &cnt);
509 
514  void buildListOfModels(std::vector<ModellingInterface*> &model);
515 
520  void buildListOfInverseKineticsLinks(std::vector<Link*> &lnk);
521 
526  void buildListOfObservers(std::vector<Observer*> &obsrv);
527 
531  void setUpInverseKinetics();
532 
536  void setUpLinks();
537 
541  bool gActiveChanged();
542 
546  bool gActiveChangedReg();
547 
551  bool detectImpact();
552 
556  void calcsvSize();
557 
561  void calclaSize(int j);
562 
566  void calcLinkStatusSize();
567 
571  void calcLinkStatusRegSize();
572 
577 
582 
586  void calcgSize(int j);
587 
598  void calcgdSize(int j);
599 
603  void calcrFactorSize(int j);
604 
608  void setUpActiveLinks();
609 
624  void checkActive(int i);
625 
629  void checkActiveReg(int i);
630 
634  virtual void setgTol(double tol);
635 
639  virtual void setgdTol(double tol);
640 
644  virtual void setgddTol(double tol);
645 
649  virtual void setlaTol(double tol);
650 
654  virtual void setLaTol(double tol);
655 
659  void setrMax(double rMax);
660 
661  void addFrame(FixedRelativeFrame *frame);
662 
667  int frameIndex(const Frame *frame_) const;
668 
672  void addGroup(DynamicSystem *dynamicsystem);
673 
679  DynamicSystem* getGroup(const std::string &name,bool check=true) const;
680 
684  void addObject(Object *object);
685 
691  Object* getObject(const std::string &name,bool check=true) const;
692 
696  void addLink(Link *link);
697 
701  void addInverseKineticsLink(Link *link);
702 
703  Observer* getObserver(const std::string &name,bool check=true) const;
704  void addObserver(Observer *element);
705 
711  Link* getLink(const std::string &name,bool check=true) const;
712 
716  void addModel(ModellingInterface *modell);
717 
723  ModellingInterface* getModel(const std::string &name, bool check=true) const;
724 
727 
728  virtual Element *getChildByContainerAndName(const std::string &container, const std::string &name) const;
729 
730  virtual void updatecorr(int j);
731  void updatecorrRef(const fmatvec::Vec &ref);
732  void calccorrSize(int j);
733 
734  void checkRoot();
735 
736  protected:
741 
745  fmatvec::Vec3 PrPF;
746 
750  fmatvec::SqrMat3 APF;
751 
755  std::vector<Object*> object;
756  std::vector<Link*> link;
757  std::vector<Link*> linkSingleValued;
758  std::vector<Link*> linkSetValued;
759  std::vector<Link*> linkSetValuedActive;
760  std::vector<ModellingInterface*> model;
761  std::vector<DynamicSystem*> dynamicsystem;
762  std::vector<Link*> inverseKineticsLink;
763  std::vector<Observer*> observer;
764  std::vector< std::vector<Link*> > linkOrdered;
765 
770 
775 
780 
784  fmatvec::Vec q, qd, q0;
785 
789  fmatvec::Vec u, ud[2], u0;
790 
794  fmatvec::Vec x, xd, x0;
795 
799  fmatvec::Vec h[2], r[2], f;
800 
804  fmatvec::Mat W[2], V[2];
805 
810 
815 
820 
825 
830 
835 
839  fmatvec::VecInt jsv;
840 
844  fmatvec::VecInt LinkStatus;
845 
849  fmatvec::VecInt LinkStatusReg;
850 
854  int qSize, qInd;
855 
859  int uSize[2], uInd[2];
860 
864  int xSize, xInd;
865 
869  int hSize[2], hInd[2];
870 
874  int gSize, gInd;
875 
879  int gdSize, gdInd;
880 
884  int laSize, laInd;
885 
889  int rFactorSize, rFactorInd;
890 
894  int svSize, svInd;
895 
899  int LinkStatusSize, LinkStatusInd;
900 
904  int LinkStatusRegSize, LinkStatusRegInd;
905 
909  std::vector<Frame*> frame;
910  std::vector<Contour*> contour;
911 
912 #ifdef HAVE_OPENMBVCPPINTERFACE
913  boost::shared_ptr<OpenMBV::Group> openMBVGrp;
914 #endif
915  boost::shared_ptr<H5::File> hdf5File;
916 
920  void addContour(Contour* contour);
921 
924 
928  int laInverseKineticsSize, bInverseKineticsSize;
929 
930  fmatvec::Mat WInverseKinetics[2], bInverseKinetics;
931  fmatvec::Vec laInverseKinetics;
932 
933  int corrSize, corrInd;
934  fmatvec::Vec corr;
935 
936  std::string saved_frameOfReference;
937  };
938 }
939 
940 #endif
941 
int rFactorSize
size and local start index of rfactors relative to parent
Definition: dynamic_system.h:889
bool gActiveChangedReg()
Definition: dynamic_system.cc:1083
virtual void setlaTol(double tol)
Definition: dynamic_system.cc:1292
fmatvec::Vec3 PrPF
relative translation with respect to parent frame
Definition: dynamic_system.h:745
fmatvec::Vec q
positions, differentiated positions, initial positions
Definition: dynamic_system.h:784
virtual int jacobianConstraints()
compute JACOBIAN of contact equations
Definition: dynamic_system.cc:581
void buildListOfModels(std::vector< ModellingInterface * > &model)
build flat list of models
Definition: dynamic_system.cc:1009
void updateLLMRef(const fmatvec::SymMat &ref, int i=0)
references to Cholesky decomposition of dynamic system parent
Definition: dynamic_system.cc:775
int laInverseKineticsSize
size of contact force parameters of special links relative to parent
Definition: dynamic_system.h:928
Link * getLink(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1386
fmatvec::SymMat M[2]
mass matrix
Definition: dynamic_system.h:774
int laSize
size and local start index of contact force parameters relative to parent
Definition: dynamic_system.h:884
Vector< Ref, double > Vec
virtual void init(InitStage stage)
plots time series header
Definition: dynamic_system.cc:418
solver interface for modelling and simulation of dynamic systems
Definition: dynamic_system_solver.h:49
void addModel(ModellingInterface *modell)
Definition: dynamic_system.cc:1400
void addGroup(DynamicSystem *dynamicsystem)
Definition: dynamic_system.cc:1423
void buildListOfObjects(std::vector< Object * > &obj)
build flat list of objects
Definition: dynamic_system.cc:965
void checkActiveReg(int i)
check if single-valued contacts are active
Definition: dynamic_system.cc:1266
void buildListOfContours(std::vector< Contour * > &cnt)
build flat list of contours
Definition: dynamic_system.cc:1000
void updateWRef(const fmatvec::Mat &ref, int i=0)
references to contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:823
fmatvec::Mat T
linear relation matrix of position and velocity parameters
Definition: dynamic_system.h:769
int frameIndex(const Frame *frame_) const
Definition: dynamic_system.cc:1331
fmatvec::VecInt LinkStatus
status of set-valued links
Definition: dynamic_system.h:844
void setUpLinks()
distribute links to set- and single valued container
Definition: dynamic_system.cc:1044
fmatvec::Vec u
velocities, differentiated velocities, initial velocities
Definition: dynamic_system.h:789
virtual void updaterFactors()
update relaxation factors for contact equations
Definition: dynamic_system.cc:609
void updateqdRef(const fmatvec::Vec &ref)
references to differentiated positions of dynamic system parent
Definition: dynamic_system.cc:653
int uSize[2]
size and local start index of velocities relative to parent
Definition: dynamic_system.h:859
fmatvec::Vec x
order one parameters, differentiated order one parameters, initial order one parameters ...
Definition: dynamic_system.h:794
cartesian frame on rigid bodies
Definition: frame.h:187
bool detectImpact()
Definition: dynamic_system.cc:1097
void updaterRef(const fmatvec::Vec &ref, int j=0)
references to nonsmooth right hand side of dynamic system parent
Definition: dynamic_system.cc:742
virtual int solveConstraintsRootFinding()
solve contact equations with Newton scheme
Definition: dynamic_system.cc:565
void updateudallRef(const fmatvec::Vec &ref, int i=0)
references to velocities of dynamic system parent
Definition: dynamic_system.cc:691
void updateresRef(const fmatvec::Vec &ref)
references to residuum of contact equations of dynamic system parent
Definition: dynamic_system.cc:871
int LinkStatusSize
size and local start index of set-valued link status vector relative to parent
Definition: dynamic_system.h:899
basic class for contour definition for rigid (which do not know about their shape) and flexible (they...
Definition: contour.h:51
virtual void closePlot()
closes plot file
Definition: dynamic_system.cc:402
fmatvec::SymMat LLM[2]
Cholesky decomposition of mass matrix.
Definition: dynamic_system.h:779
virtual void setgTol(double tol)
Definition: dynamic_system.cc:1271
basic class of MBSim mainly for plotting
Definition: element.h:58
void updateqRef(const fmatvec::Vec &ref)
references to positions of dynamic system parent
Definition: dynamic_system.cc:643
virtual void plotAtSpecialEvent(double t, double dt=1.)
plots time dependent data at special events
Definition: dynamic_system.cc:387
void updateWInverseKineticsRef(const fmatvec::Mat &ref, int i=0)
references to contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:830
int gSize
size and local start index of relative distances relative to parent
Definition: dynamic_system.h:874
virtual void checkImpactsForTermination(double dt)
validate force laws concerning given tolerances on velocity level
Definition: dynamic_system.cc:603
virtual void setLaTol(double tol)
Definition: dynamic_system.cc:1299
FixedRelativeFrame * I
Definition: dynamic_system.h:923
void updatesvRef(const fmatvec::Vec &ref)
references to stopvector (rootfunction for event driven integrator) of dynamic system parent ...
Definition: dynamic_system.cc:851
virtual Contour * getContour(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:629
int gdSize
size and local start index of relative velocities relative to parent
Definition: dynamic_system.h:879
virtual void plot(double t, double dt)
plots time dependent data
Definition: dynamic_system.cc:368
fmatvec::Vec h[2]
smooth, smooth with respect to objects, smooth with respect to links, nonsmooth and order one right h...
Definition: dynamic_system.h:799
virtual int solveImpactsFixpointSingle(double dt)
solve impact equations with single step fixed point scheme on velocity level
Definition: dynamic_system.cc:541
void calcrFactorSize(int j)
calculates size of relaxation factors for contact equations
Definition: dynamic_system.cc:1240
void calcbInverseKineticsSize()
calculates size of contact force parameters
Definition: dynamic_system.cc:1210
virtual int solveImpactsRootFinding(double dt)
solve impact equations with Newton scheme on velocity level
Definition: dynamic_system.cc:573
std::vector< Frame * > frame
vector of frames and contours
Definition: dynamic_system.h:909
void updateLinkStatusRegRef(const fmatvec::VecInt &LinkStatusRegParent)
references to status vector of single valued links
Definition: dynamic_system.cc:895
fmatvec::Vec res
residuum of nonlinear contact equations for Newton scheme
Definition: dynamic_system.h:824
void updateudRef(const fmatvec::Vec &ref, int i=0)
references to differentiated velocities of dynamic system parent
Definition: dynamic_system.cc:681
void updategdRef(const fmatvec::Vec &ref)
references to relative velocities of dynamic system parent
Definition: dynamic_system.cc:792
void buildListOfSetValuedLinks(std::vector< Link * > &lnk)
build flat list of all setvalued links
Definition: dynamic_system.cc:983
fmatvec::SqrMat3 APF
relative rotation with respect to parent frame
Definition: dynamic_system.h:750
void updaterFactorRef(const fmatvec::Vec &ref)
references to relaxation factors for contact equations of dynamic system parent
Definition: dynamic_system.cc:878
void buildListOfFrames(std::vector< Frame * > &frm)
build flat list of frames
Definition: dynamic_system.cc:991
dynamic system as topmost hierarchical level
Definition: dynamic_system.h:57
void updateuallRef(const fmatvec::Vec &ref)
references to velocities of dynamic system parent
Definition: dynamic_system.cc:673
void buildListOfObservers(std::vector< Observer * > &obsrv)
build flat list of observers
Definition: dynamic_system.cc:1027
virtual std::string getType() const
Definition: dynamic_system.h:140
Frame * R
parent frame
Definition: dynamic_system.h:740
void calclaInverseKineticsSize()
calculates size of contact force parameters
Definition: dynamic_system.cc:1200
void calcLinkStatusSize()
calculates size of set-valued link status vector
Definition: dynamic_system.cc:1147
void updategRef(const fmatvec::Vec &ref)
references to relative distances of dynamic system parent
Definition: dynamic_system.cc:785
fmatvec::VecInt LinkStatusReg
status of single-valued links
Definition: dynamic_system.h:849
void checkActive(int i)
check if set-valued contacts are active and set corresponding attributes
Definition: dynamic_system.cc:1260
InitStage
The stages of the initialization.
Definition: element.h:97
virtual Element * getChildByContainerAndName(const std::string &container, const std::string &name) const
Get the Element named name in the container named container.
Definition: dynamic_system.cc:1441
int qSize
size and local start index of positions relative to parent
Definition: dynamic_system.h:854
virtual int solveConstraintsGaussSeidel()
solve contact equations with Gauss-Seidel scheme
Definition: dynamic_system.cc:549
DynamicSystem(const std::string &name)
constructor
Definition: dynamic_system.cc:51
virtual int solveImpactsGaussSeidel(double dt)
solve impact equations with Gauss-Seidel scheme on velocity level
Definition: dynamic_system.cc:557
FixedRelativeFrame * getFrameI()
Definition: dynamic_system.h:726
void setUpActiveLinks()
rearrange vector of active setvalued links
Definition: dynamic_system.cc:1250
PlotFeatureStatus getPlotFeature(PlotFeature pf)
Definition: element.h:217
virtual void setgdTol(double tol)
Definition: dynamic_system.cc:1278
void buildListOfLinks(std::vector< Link * > &lnk)
build flat list of links
Definition: dynamic_system.cc:974
void updateTRef(const fmatvec::Mat &ref)
references to linear transformation matrix between differentiated positions and velocities of dynamic...
Definition: dynamic_system.cc:755
Object * getObject(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1353
void updateVRef(const fmatvec::Mat &ref, int i=0)
references to condensed contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:844
void buildListOfInverseKineticsLinks(std::vector< Link * > &lnk)
build flat list of inverse kinetics links
Definition: dynamic_system.cc:1018
void setrMax(double rMax)
Definition: dynamic_system.cc:1306
std::string name
name of element
Definition: element.h:290
virtual void facLLM(int i=0)=0
compute Cholesky decomposition of mass matrix TODO necessary?
int LinkStatusRegSize
size and local start index of single-valued link status vector relative to parent ...
Definition: dynamic_system.h:904
DynamicSystem * getGroup(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1339
void addLink(Link *link)
Definition: dynamic_system.cc:1367
void calcgdSize(int j)
calculates size of gap velocities
Definition: dynamic_system.cc:1230
ModellingInterface * getModel(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:1409
fmatvec::Vec wb
TODO.
Definition: dynamic_system.h:819
void updatejsvRef(const fmatvec::VecInt &ref)
references to boolean evaluation of stopvector concerning roots of dynamic system parent ...
Definition: dynamic_system.cc:861
void addContour(Contour *contour)
Definition: dynamic_system.cc:1322
virtual void checkConstraintsForTermination()
validate force laws concerning given tolerances
Definition: dynamic_system.cc:597
int svSize
size and local start index of stop vector relative to parent
Definition: dynamic_system.h:894
cartesian frame on bodies used for application of e.g. links and loads
Definition: frame.h:39
fmatvec::Vec rFactor
rfactors for relaxation nonlinear contact equations
Definition: dynamic_system.h:829
void calcgSize(int j)
calculates size of relative distances
Definition: dynamic_system.cc:1220
fmatvec::Vec sv
stop vector (root functions for event driven integration
Definition: dynamic_system.h:834
fmatvec::Vec g
relative distances and velocities
Definition: dynamic_system.h:814
void calcLinkStatusRegSize()
calculates size of single-valued link status vector
Definition: dynamic_system.cc:1162
bool gActiveChanged()
Definition: dynamic_system.cc:1069
void updatehRef(const fmatvec::Vec &hRef, int i=0)
references to smooth right hand side of dynamic system parent
Definition: dynamic_system.cc:725
void addInverseKineticsLink(Link *link)
Definition: dynamic_system.cc:1377
H5::GroupBase * plotGroup
associated plot group
Definition: element.h:324
void calcsvSize()
calculates size of stop vector
Definition: dynamic_system.cc:1177
virtual ~DynamicSystem()
destructor
Definition: dynamic_system.cc:67
std::vector< Object * > object
container for possible ingredients
Definition: dynamic_system.h:755
fmatvec::Vec la
contact force parameters
Definition: dynamic_system.h:809
virtual void setDynamicSystemSolver(DynamicSystemSolver *sys)
sets the used dynamics system solver to the element
Definition: dynamic_system.cc:356
void updateuRef(const fmatvec::Vec &ref)
references to velocities of dynamic system parent
Definition: dynamic_system.cc:663
void setUpInverseKinetics()
analyse constraints of dynamic systems for usage in inverse kinetics
Definition: dynamic_system.cc:1036
int hSize[2]
size and local start index of order smooth right hand side relative to parent
Definition: dynamic_system.h:869
virtual Frame * getFrame(const std::string &name, bool check=true) const
Definition: dynamic_system.cc:615
int xSize
size and local start index of order one parameters relative to parent
Definition: dynamic_system.h:864
void calclaSize(int j)
calculates size of contact force parameters
Definition: dynamic_system.cc:1190
fmatvec::VecInt jsv
boolean evaluation of stop vector concerning roots
Definition: dynamic_system.h:839
PlotFeatureStatus
Plot feature status.
Definition: element.h:61
virtual int solveConstraintsFixpointSingle()
solve contact equations with single step fixed point scheme
Definition: dynamic_system.cc:533
virtual int jacobianImpacts()
compute JACOBIAN of contact equations on velocity level
Definition: dynamic_system.cc:589
void updateLinkStatusRef(const fmatvec::VecInt &LinkStatusParent)
references to status vector of set valued links with piecewise link equations (which piece is valid) ...
Definition: dynamic_system.cc:885
void buildListOfDynamicSystems(std::vector< DynamicSystem * > &sys)
build flat list of dynamic systems
Definition: dynamic_system.cc:956
PlotFeature
Plot Features.
Definition: element.h:74
void updatewbRef(const fmatvec::Vec &ref)
references to TODO of dynamic system parent
Definition: dynamic_system.cc:816
void updateMRef(const fmatvec::SymMat &ref, int i=0)
references to mass matrix of dynamic system parent
Definition: dynamic_system.cc:765
virtual void setgddTol(double tol)
Definition: dynamic_system.cc:1285
PlotFeatureStatus getPlotFeatureForChildren(PlotFeature pf)
Definition: element.h:222
void addObject(Object *object)
Definition: dynamic_system.cc:1432
void updatelaRef(const fmatvec::Vec &ref)
references to contact force parameters of dynamic system parent
Definition: dynamic_system.cc:799

Impressum / Disclaimer / Datenschutz Generated by doxygen 1.8.5 Valid HTML