Loading [MathJax]/extensions/tex2jax.js
mbsim  4.0.0
MBSim Kernel
All Classes Namespaces Functions Variables Typedefs Enumerations Enumerator Pages
dynamic_system_solver.h
1/* Copyright (C) 2004-2018 MBSim Development Team
2 *
3 * This library is free software; you can redistribute it and/or
4 * modify it under the terms of the GNU Lesser General Public
5 * License as published by the Free Software Foundation; either
6 * version 2.1 of the License, or (at your option) any later version.
7 *
8 * This library is distributed in the hope that it will be useful,
9 * but WITHOUT ANY WARRANTY; without even the implied warranty of
10 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
11 * Lesser General Public License for more details.
12 *
13 * You should have received a copy of the GNU Lesser General Public
14 * License along with this library; if not, write to the Free Software
15 * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA
16 *
17 * Contact: martin.o.foerg@googlemail.com
18 */
19
20#ifndef _DYNAMIC_SYSTEM_SOLVER_H_
21#define _DYNAMIC_SYSTEM_SOLVER_H_
22
23#include "mbsim/group.h"
24#include "fmatvec/sparse_matrix.h"
25#include "mbsim/functions/function.h"
26#include "mbsim/environment.h"
27
28#include <atomic>
29
30namespace MBSim {
31
32 class Graph;
33 class MBSimEnvironment;
34 class MultiDimNewtonMethod;
35 class ConstraintResiduum;
36 class ConstraintJacobian;
37
38 struct StateTable {
39 std::string name;
40 char label;
41 int number;
42 StateTable() = default;
43 StateTable(const std::string &name_, char label_, int number_) : name(name_), label(label_), number(number_) { }
44 };
45
61 class DynamicSystemSolver : public Group {
62
63 class Residuum : public MBSim::Function<fmatvec::Vec(fmatvec::Vec)> {
64 public:
65 Residuum(MBSim::DynamicSystemSolver *sys_) : sys(sys_) { }
66 fmatvec::Vec operator()(const fmatvec::Vec &z);
67 private:
69 };
70
71 public:
72
76 enum Solver { fixedpoint, GaussSeidel, direct, rootfinding, unknownSolver, directNonlinear };
77
82 DynamicSystemSolver(const std::string &name="");
83
87 ~DynamicSystemSolver() override;
88
127 void initialize();
128
129 /*
130 * If true (the default) then
131 * the simulation output files (h5 files) are deleted/truncated/regenerated.
132 * If false, then these files left are untouched from a previous run.
133 * This is useful to regenerate the
134 * e.g. OpenMBV XML files without doing a new integration.
135 */
136 void setTruncateSimulationFiles(bool trunc) { truncateSimulationFiles=trunc; }
137
138 /* INHERITED INTERFACE OF GROUP */
139 void init(InitStage stage, const InitConfigSet &config) override;
140 using Group::plot;
141 /***************************************************/
142
143 /* INHERITED INTERFACE OF DYNAMICSYSTEM */
144 int solveConstraintsFixpointSingle() override;
145 int solveImpactsFixpointSingle() override;
146 int solveConstraintsGaussSeidel() override;
147 int solveImpactsGaussSeidel() override;
148 int solveConstraintsRootFinding() override;
149 int solveImpactsRootFinding() override;
150 void checkConstraintsForTermination() override;
151 void checkImpactsForTermination() override;
152 /***************************************************/
153
154 /* INHERITED INTERFACE OF OBJECTINTERFACE */
155 void updateT() override;
156 void updateh(int i=0) override;
157 void updateM() override;
158 void updateLLM() override;
159 void updatezd() override;
160 /***************************************************/
161
162 /* INHERITED INTERFACE OF LINKINTERFACE */
163
174 virtual void updater(int j=0);
175 virtual void updateJrla(int j=0);
176 virtual void updaterdt();
177 void updatewb() override;
178 void updateg() override;
179 void updategd() override;
180 void updateW(int j=0) override;
181 void updateV(int j=0) override;
182 virtual void updatebc();
183 virtual void updatebi();
184 virtual void updatela();
185 virtual void updateLa();
186 /***************************************************/
187
188 /* INHERITED INTERFACE OF ELEMENT */
189 /***************************************************/
190
191 /* GETTER / SETTER */
192
193 void setSmoothSolver(Solver solver_) { smoothSolver = solver_; }
194 void setConstraintSolver(Solver solver_) { contactSolver = solver_; }
195 void setImpactSolver(Solver solver_) { impactSolver = solver_; }
196 const Solver& getSmoothSolver() { return smoothSolver; }
197 const Solver& getConstraintSolver() { return contactSolver; }
198 const Solver& getImpactSolver() { return impactSolver; }
199 void setTermination(bool term_) { term = term_; }
200 void setMaximumNumberOfIterations(int iter) { maxIter = iter; }
201 void setHighNumberOfIterations(int iter) { highIter = iter; }
202 void setNumericalJacobian(bool numJac_) { numJac = numJac_; }
203 void setMaximumDampingSteps(int maxDSteps) { maxDampingSteps = maxDSteps; }
204 void setLevenbergMarquardtParamater(double lmParm_) { lmParm = lmParm_; }
205
206 void setUseOldla(bool flag) { useOldla = flag; }
207 void setDecreaseLevels(const fmatvec::VecInt &decreaseLevels_) { decreaseLevels = decreaseLevels_; }
208 void setCheckTermLevels(const fmatvec::VecInt &checkTermLevels_) { checkTermLevels = checkTermLevels_; }
209 void setCheckGSize(bool checkGSize_) { checkGSize = checkGSize_; }
210 void setLimitGSize(int limitGSize_) { limitGSize = limitGSize_; checkGSize = false; }
211
212 double& getTime() { return t; }
213 const double& getTime() const { return t; }
214 void setTime(double t_) { t = t_; }
215
216 double getStepSize() const { return dt; }
217 void setStepSize(double dt_) { dt = dt_; }
218
219 int getzSize() const { return zSize; }
220
221 fmatvec::Vec& getState() { return z; }
222 const fmatvec::Vec& getState() const { return z; }
223 void setState(const fmatvec::Vec &z_) { z = z_; }
224
225 const fmatvec::Vec& getzd(bool check=true) const { assert((not check) or (not updzd)); return zd; }
226 void setzd(const fmatvec::Vec &zd_) { zd = zd_; }
227
228 const fmatvec::SqrMat& getG(bool check=true) const { assert((not check) or (not updG)); return G; }
229 const fmatvec::SparseMat& getGs(bool check=true) const { assert((not check) or (not updG)); return Gs; }
230 const fmatvec::Vec& getbc(bool check=true) const { assert((not check) or (not updbc)); return bc; }
231 const fmatvec::Vec& getbi(bool check=true) const { assert((not check) or (not updbi)); return bi; }
232 const fmatvec::SqrMat& getJprox() const { return Jprox; }
233 fmatvec::SqrMat& getG(bool check=true) { assert((not check) or (not updG)); return G; }
234 fmatvec::SparseMat& getGs(bool check=true) { assert((not check) or (not updG)); return Gs; }
235 fmatvec::Vec& getbc(bool check=true) { assert((not check) or (not updbc)); return bc; }
236 fmatvec::Vec& getbi(bool check=true) { assert((not check) or (not updbi)); return bi; }
237 fmatvec::SqrMat& getJprox() { return Jprox; }
238 const fmatvec::Vec& getr(int j=0, bool check=true) const { assert((not check) or (not updr[j])); return r[j]; }
239 fmatvec::Vec& getr(int j=0, bool check=true) { assert((not check) or (not updr[j])); return r[j]; }
240
241 // Jacobian of dr/dla
242 const fmatvec::Mat& getJrla(int j=0, bool check=true) const { assert((not check) or (not updJrla[j])); return Jrla[j]; }
243 fmatvec::Mat& getJrla(int j=0, bool check=true) { assert((not check) or (not updJrla[j])); return Jrla[j]; }
244
245 const fmatvec::Vec& evaldq() { if(upddq) updatedq(); return dq; }
246 const fmatvec::Vec& evaldu() { if(upddu) updatedu(); return du; }
247 const fmatvec::Vec& evaldx() { if(upddx) updatedx(); return dx; }
248 const fmatvec::Vec& evalzd();
249 const fmatvec::SqrMat& evalG() { if(updG) updateG(); return G; }
250 const fmatvec::SparseMat& evalGs() { if(updG) updateG(); return Gs; }
251 const fmatvec::Vec& evalbc() { if(updbc) updatebc(); return bc; }
252 const fmatvec::Vec& evalbi() { if(updbi) updatebi(); return bi; }
253 const fmatvec::Vec& evalsv();
254 const fmatvec::Vec& evalz0();
255 const fmatvec::Vec& evalla() { if(updla) updatela(); return la; }
256 const fmatvec::Vec& evalLa() { if(updLa) updateLa(); return La; }
257
258 fmatvec::Vec& getzParent() { return zParent; }
259 fmatvec::Vec& getzdParent() { return zdParent; }
260 fmatvec::Vec& getlaParent() { return laParent; }
261 fmatvec::Vec& getLaParent() { return LaParent; }
262// const fmatvec::Vec& getzParent() const { return zParent; }
263// const fmatvec::Vec& getzdParent() const { return zdParent; }
264// const fmatvec::Mat& getWParent(int i=0) const { return WParent[i]; }
265 fmatvec::Mat& getWParent(int i=0) { return WParent[i]; }
266 fmatvec::Mat& getVParent(int i=0) { return VParent[i]; }
267// fmatvec::Vec& getlaParent() { return laParent; }
268// fmatvec::Vec& getLaParent() { return LaParent; }
269// const fmatvec::Vec& getgParent() const { return gParent; }
270 fmatvec::Vec& getgParent() { return gParent; }
271// const fmatvec::Vec& getgdParent() const { return gdParent; }
272 fmatvec::Vec& getgdParent() { return gdParent; }
273 fmatvec::Vec& getresParent() { return resParent; }
274 fmatvec::Vec& getrFactorParent() { return rFactorParent; }
275
276 void resizezParent(int nz) { zParent.resize(nz); }
277 void resizezdParent(int nz) { zdParent.resize(nz); }
278
279 DynamicSystemSolver* getDynamicSystemSolver() { return this; }
280 int getMaxIter() {return maxIter;}
281 int getHighIter() {return highIter;}
282 int getIterC() {return iterc;}
283 int getIterI() {return iteri;}
284 /***************************************************/
285
290
291 int (DynamicSystemSolver::*solveSmooth_)();
292
298
305
312
319
327
328 int solveImpactsNonlinearEquations();
329
334 virtual void updateG();
335
339 void decreaserFactors();
340
347 virtual const fmatvec::Vec& shift(std::optional<std::reference_wrapper<bool>> &&velProjWasCalled={},
348 std::optional<std::reference_wrapper<bool>> &&posProjWasCalled={});
349
354 void getLinkStatus(fmatvec::VecInt &LinkStatusExt);
355
360 void getLinkStatusReg(fmatvec::VecInt &LinkStatusRegExt);
361
365 bool positionDriftCompensationNeeded(double gmax);
366
370 bool velocityDriftCompensationNeeded(double gdmax);
371
375 void projectGeneralizedPositions(int mode, bool fullUpdate=false);
376
380 void projectGeneralizedVelocities(int mode);
381
386 void savela();
387
392 void initla();
393
398 void saveLa();
399
404 void initLa();
405
410 void addElement(Element *element_);
411
417 Element* getElement(const std::string &name);
418
422 std::string getSolverInfo();
423
428 void setStopIfNoConvergence(bool flag, bool dropInfo = false) { stopIfNoConvergence = flag; dropContactInfo=dropInfo; }
429
430 bool getStopIfNoConvergence() { return stopIfNoConvergence; }
431
435 void dropContactMatrices();
436
437 // MBSim signal handler
438#ifndef SWIG
440 public:
443 private:
444 using SigHandle = void (*)(int);
445 #ifndef _WIN32
446 SigHandle oldSigHup;
447 #endif
448 SigHandle oldSigInt;
449 SigHandle oldSigTerm;
450 };
451#endif
452
453 void throwIfExitRequested() {
454 if(exitRequest) {
455 if(!exitRequestPrinted) {
456 exitRequestPrinted = true;
457 msg(fmatvec::Atom::Error)<<"User requested a exit (throw exception now)."<<std::endl;
458 }
459 throw std::runtime_error("Exception due to user requested exit.");
460 }
461 }
462
463 bool exitRequested() {
464 if(exitRequest)
465 if(!exitRequestPrinted) {
466 exitRequestPrinted = true;
467 msg(fmatvec::Atom::Error)<<"User requested a exit (caller will handle this request now)."<<std::endl;
468 }
469 return exitRequest;
470 }
471
472 // TODO just for testing
473 void setPartialEventDrivenSolver(bool peds_) { peds = peds_; }
474
480 void writez(std::string fileName, bool formatH5=true);
481
487 void writeStateTable(std::string fileName);
488
493 void readz0(std::string fileName);
494
495 void initializeUsingXML(xercesc::DOMElement *element) override;
496
501 void setProjectionTolerance(double tol) { tolProj = tol; }
502
507 void setLocalSolverTolerance(double tol) { tolLocalSolver = tol; }
508 double getLocalSolverTolerance() { return tolLocalSolver; }
509
514 void setDynamicSystemSolverTolerance(double tol) { tolDSS = tol; }
515 double getDynamicSystemSolverTolerance() { return tolDSS; }
516
521 void updatezRef(fmatvec::Vec &ext);
522
527 void updatezdRef(fmatvec::Vec &ext);
528
529 void setAlwaysConsiderContact(bool alwaysConsiderContact_) { alwaysConsiderContact = alwaysConsiderContact_; }
530
531 void setInverseKinetics(bool inverseKinetics_) { inverseKinetics = inverseKinetics_; }
532 bool getInverseKinetics() const { return inverseKinetics; }
533
534 void setInitialProjection(bool initialProjection_) { initialProjection = initialProjection_; }
535 bool getInitialProjection() const { return initialProjection; }
536
537 void setDetermineEquilibriumState(bool determineEquilibriumState_) { determineEquilibriumState = determineEquilibriumState_; }
538 bool getDetermineEquilibriumState() const { return determineEquilibriumState; }
539
540 void setUseConstraintSolverForPlot(bool useConstraintSolverForPlot_) { useConstraintSolverForPlot = useConstraintSolverForPlot_; }
541 bool getUseConstraintSolverForPlot() const { return useConstraintSolverForPlot; }
542
543 fmatvec::Mat dhdq(int lb=0, int ub=0);
544 fmatvec::Mat dhdu(int lb=0, int ub=0);
545 fmatvec::Mat dhdx();
546 fmatvec::Vec dhdt();
547
548 void setRootID(int ID) {rootID = ID;}
549 int getRootID() const {return rootID;}
550
551 void resetUpToDate() override;
552
553 bool getUpdateT() { return updT; }
554 bool getUpdateM() { return updM; }
555 bool getUpdateLLM() { return updLLM; }
556 bool getUpdateh(int j) { return updh[j]; }
557 bool getUpdater(int j) { return updr[j]; }
558 bool getUpdateJrla(int j) { return updJrla[j]; }
559 bool getUpdaterdt() { return updrdt; }
560 bool getUpdateW(int j) { return updW[j]; }
561 bool getUpdateV(int j) { return updV[j]; }
562 bool getUpdatewb() { return updwb; }
563 bool getUpdateg() { return updg; }
564 bool getUpdategd() { return updgd; }
565 bool getUpdatela() { return updla; }
566 bool getUpdateLa() { return updLa; }
567 bool getUpdatezd() { return updzd; }
568 bool getUpdatedq() { return upddq; }
569 bool getUpdatedu() { return upddu; }
570 bool getUpdatedx() { return upddx; }
571 void setUpdatela(bool updla_) { updla = updla_; }
572 void setUpdateLa(bool updLa_) { updLa = updLa_; }
573 void setUpdateG(bool updG_) { updG = updG_; }
574 void setUpdatebi(bool updbi_) { updbi = updbi_; }
575 void setUpdatebc(bool updbc_) { updbc = updbc_; }
576 void setUpdatezd(bool updzd_) { updzd = updzd_; }
577 void setUpdateW(bool updW_, int i=0) { updW[i] = updW_; }
578
583 void updategRef(fmatvec::Vec &ref) override { Group::updategRef(ref); updg = true; }
584
589 void updategdRef(fmatvec::Vec &ref) override { Group::updategdRef(ref); updgd = true; }
590
595 void updatewbRef(fmatvec::Vec &ref) override { Group::updatewbRef(ref); updwb = true; }
596
602 void updateWRef(fmatvec::Mat &ref, int i=0) override { Group::updateWRef(ref,i); updW[i] = true; }
603
609 void updateVRef(fmatvec::Mat &ref, int i=0) override { Group::updateVRef(ref,i); updV[i] = true; }
610
614 virtual void updatelaInverseKinetics();
615
616 void updatedq() override;
617 void updatedu() override;
618 void updatedx() override;
619
620 void updateStopVector() override;
621
622 void updateInternalState();
623
624 void plot() override;
625
626 void addEnvironment(Environment* env);
627
632 template<class Env>
633 Env* getEnvironment();
634
638
639 std::vector<StateTable>& getStateTable() { return tabz; }
640
641 void setCompressionLevel(int level) { compressionLevel=level; }
642 void setChunkSize(int size) { chunkSize=size; }
643 void setCacheSize(int size) { cacheSize=size; }
644
645 void setqdequ(bool qdequ_) { qdequ = qdequ_; }
646 bool getqdequ() { return qdequ; }
647
648 protected:
652 double t;
653
657 double dt;
658
662 int zSize;
663
667 fmatvec::Vec z;
668
672 fmatvec::Vec zd;
673
677 fmatvec::SymMat MParent;
678
682 fmatvec::SymMat LLMParent;
683
687 fmatvec::Mat TParent;
688
692 fmatvec::Mat WParent[2];
693
697 fmatvec::Mat VParent[2];
698
702 fmatvec::Vec wbParent;
703
707 fmatvec::Vec laParent, LaParent;
708
712 fmatvec::Vec rFactorParent;
713
717 fmatvec::Vec sParent;
718
722 fmatvec::Vec resParent;
723
724 fmatvec::Vec curisParent, nextisParent;
725
729 fmatvec::Vec gParent;
730
734 fmatvec::Vec gdParent;
735
739 fmatvec::Vec zParent;
740
744 fmatvec::Vec zdParent;
745
746 fmatvec::Vec dxParent, dqParent, duParent;
747
751 fmatvec::Vec hParent[2];
752
756 fmatvec::Vec rParent[2], rdtParent;
757
758 fmatvec::Mat JrlaParent[2];
759
763 fmatvec::Vec svParent;
764
768 fmatvec::VecInt jsvParent;
769
773 fmatvec::VecInt LinkStatusParent;
774
779 fmatvec::VecInt LinkStatusRegParent;
780
784 fmatvec::SparseMat Gs;
785
789 fmatvec::SqrMat Jprox;
790
794 fmatvec::SqrMat G;
795
799 fmatvec::Vec bc, bi;
800
804 bool term;
805
809 int maxIter, highIter, maxDampingSteps, iterc, iteri;
810
814 double lmParm;
815
819 Solver smoothSolver, contactSolver, impactSolver;
820
825
830
835
839 bool numJac;
840
844 fmatvec::VecInt decreaseLevels;
845
849 fmatvec::VecInt checkTermLevels;
850
855
860
864 bool peds;
865
869 double tolProj;
870
874 double tolLocalSolver { 1e-10 };
875
879 double tolDSS { 1e-9 };
880
884 fmatvec::Vec laInverseKineticsParent;
885 fmatvec::Mat bInverseKineticsParent;
886
890 fmatvec::Mat WInverseKineticsParent;
891
892 bool alwaysConsiderContact;
893 bool inverseKinetics;
894 bool initialProjection;
895 bool determineEquilibriumState;
896 bool useConstraintSolverForPlot;
897
898 fmatvec::Vec corrParent;
899
900 int rootID;
901
902 double gTol, gdTol, gddTol, laTol, LaTol;
903
904 bool updT, updh[2], updr[2], updJrla[2], updrdt, updM, updLLM, updW[2], updV[2], updwb, updg, updgd, updG, updbc, updbi, updsv, updzd, updla, updLa, upddq, upddu, upddx;
905
906 bool useSmoothSolver;
907
908 bool qdequ;
909
910 std::vector<StateTable> tabz;
911
912 int compressionLevel { H5::File::getDefaultCompression() };
913 int chunkSize { H5::File::getDefaultChunkSize() };
914 int cacheSize { H5::File::getDefaultCacheSize() };
915
916 private:
920 static void sigInterruptHandler(int);
921
925 void constructor();
926
930 static std::atomic<bool> exitRequest;
931 static_assert(decltype(exitRequest)::is_always_lock_free);
932 bool exitRequestPrinted;
933
937 bool READZ0;
938
947 void addToGraph(Graph* graph, fmatvec::SqrMat &A, int i, std::vector<Element*> &objList);
948
949 bool truncateSimulationFiles;
950
951 // Holds the dynamic systems before the "reorganize hierarchy" takes place.
952 // This is required since all elements of all other containers from DynamicSystem are readded to DynamicSystemSolver,
953 // except this container (which is a "pure" container = no calculation is done in DynamicSystem)
954 // However, we must hold this container until the dtor of DynamicSystemSolver is called to avoid the deallocation of other
955 // elements hold by DynamicSystem elements (especially (currently only) hdf5File)
956 std::vector<DynamicSystem*> dynamicsystemPreReorganize;
957
958 double facSizeGs;
959
960 std::vector<std::unique_ptr<Environment>> environments;
961 MBSimEnvironment *mbsimEnvironment = nullptr;
962
963 bool firstPlot { true };
964
965 std::unique_ptr<MultiDimNewtonMethod> nonlinearConstraintNewtonSolver;
966 std::unique_ptr<ConstraintResiduum> constraintResiduum;
967 std::unique_ptr<ConstraintJacobian> constraintJacobian;
968 };
969
970 template<class Env>
972 // get the Environment of type Env
973 auto &reqType=typeid(Env);
974 for(auto &e : environments) {
975 auto &e_=*e;
976 if(reqType==typeid(e_))
977 return static_cast<Env*>(e.get());
978 }
979 auto newEnv=new Env;
980 addEnvironment(newEnv);
981 return newEnv;
982 }
983
984}
985
986#endif /* _DYNAMIC_SYSTEM_SOLVER_H_ */
Definition: dynamic_system_solver.h:63
Definition: dynamic_system_solver.h:439
solver interface for modelling and simulation of dynamic systems
Definition: dynamic_system_solver.h:61
void checkConstraintsForTermination() override
validate force laws concerning given tolerances
Definition: dynamic_system_solver.cc:724
fmatvec::Vec sParent
TODO.
Definition: dynamic_system_solver.h:717
int(DynamicSystemSolver::* solveImpacts_)()
function pointer for election of prox-solver for impact equations on velocity level
Definition: dynamic_system_solver.h:304
void plot() override
plots time dependent data
Definition: dynamic_system_solver.cc:1748
bool term
boolean to check for termination of contact equations solution
Definition: dynamic_system_solver.h:804
int solveConstraintsFixpointSingle() override
solve contact equations with single step fixed point scheme
Definition: dynamic_system_solver.cc:492
void computeInitialCondition()
compute initial condition for links for event driven integrator
Definition: dynamic_system_solver.cc:899
fmatvec::Vec hParent[2]
smooth, smooth with respect to objects, smooth with respect to links right hand side
Definition: dynamic_system_solver.h:751
void setProjectionTolerance(double tol)
set tolerance for projection of generalized position
Definition: dynamic_system_solver.h:501
fmatvec::Vec zd
derivative of state vector
Definition: dynamic_system_solver.h:672
bool useOldla
flag if contac force parameter of last time step should be used
Definition: dynamic_system_solver.h:834
fmatvec::Mat WParent[2]
contact force directions
Definition: dynamic_system_solver.h:692
fmatvec::VecInt LinkStatusParent
status vector of set valued links with piecewise link equation (which piece is valid)
Definition: dynamic_system_solver.h:773
fmatvec::Vec bc
TODO.
Definition: dynamic_system_solver.h:799
int solveConstraintsLinearEquations()
solution of contact equations with direct linear solver using minimum of least squares
Definition: dynamic_system_solver.cc:949
void setLocalSolverTolerance(double tol)
set tolerance for local none-linear solver (solvers on element level), like the Newton-Solver in Join...
Definition: dynamic_system_solver.h:507
fmatvec::SqrMat Jprox
JACOBIAN of contact equations for Newton scheme.
Definition: dynamic_system_solver.h:789
void initialize()
Initialize the system.
Definition: dynamic_system_solver.cc:104
void updateWRef(fmatvec::Mat &ref, int i=0) override
references to contact force direction matrix of dynamic system parent
Definition: dynamic_system_solver.h:602
fmatvec::Mat VParent[2]
condensed contact force directions
Definition: dynamic_system_solver.h:697
double lmParm
Levenberg-Marquard parameter.
Definition: dynamic_system_solver.h:814
fmatvec::Vec gParent
relative distances
Definition: dynamic_system_solver.h:729
virtual void updatelaInverseKinetics()
update inverse kinetics constraint forces
Definition: dynamic_system_solver.cc:1668
fmatvec::Vec z
state vector
Definition: dynamic_system_solver.h:667
bool peds
TODO, flag for occuring impact and sticking in event driven solver.
Definition: dynamic_system_solver.h:864
static void sigInterruptHandler(int)
handler for user interrupt signal
Definition: dynamic_system_solver.cc:1359
void writeStateTable(std::string fileName)
writes state table to a file
Definition: dynamic_system_solver.cc:1385
bool numJac
flag if Jacobian for Newton scheme should be calculated numerically
Definition: dynamic_system_solver.h:839
void saveLa()
save contact impulses for use as starting value in next time step
Definition: dynamic_system_solver.cc:1256
fmatvec::Vec wbParent
TODO.
Definition: dynamic_system_solver.h:702
double t
time
Definition: dynamic_system_solver.h:652
Solver smoothSolver
solver for contact equations and impact equations
Definition: dynamic_system_solver.h:819
int solveConstraintsRootFinding() override
solve contact equations with Newton scheme
Definition: dynamic_system_solver.cc:592
double tolProj
Tolerance for projection of generalized position.
Definition: dynamic_system_solver.h:869
fmatvec::Vec resParent
residuum of contact equations
Definition: dynamic_system_solver.h:722
int(DynamicSystemSolver::* solveConstraints_)()
Definition: dynamic_system_solver.h:297
bool dropContactInfo
flag if contact matrices for debugging should be dropped in no-convergence case
Definition: dynamic_system_solver.h:829
fmatvec::VecInt decreaseLevels
decreasing relaxation factors is done in levels containing the number of contact iterations as condit...
Definition: dynamic_system_solver.h:844
fmatvec::VecInt checkTermLevels
TODO.
Definition: dynamic_system_solver.h:849
bool READZ0
is a state read from a file
Definition: dynamic_system_solver.h:937
int solveConstraintsGaussSeidel() override
solve contact equations with Gauss-Seidel scheme
Definition: dynamic_system_solver.cc:552
fmatvec::SqrMat G
mass action matrix
Definition: dynamic_system_solver.h:794
std::string getSolverInfo()
Definition: dynamic_system_solver.cc:1305
int solveConstraintsNonlinearEquations()
solution of contact equations with direct nonlinear newton solver using minimum of least squares in e...
Definition: dynamic_system_solver.cc:968
fmatvec::VecInt LinkStatusRegParent
status vector of single valued links
Definition: dynamic_system_solver.h:779
void updategRef(fmatvec::Vec &ref) override
references to relative distances of dynamic system parent
Definition: dynamic_system_solver.h:583
void updatezRef(fmatvec::Vec &ext)
references to external state
Definition: dynamic_system_solver.cc:1404
void setStopIfNoConvergence(bool flag, bool dropInfo=false)
Definition: dynamic_system_solver.h:428
int solveImpactsLinearEquations()
solution of contact equations with Cholesky decomposition on velocity level
Definition: dynamic_system_solver.cc:978
void writez(std::string fileName, bool formatH5=true)
writes state to a file
Definition: dynamic_system_solver.cc:1363
void addToGraph(Graph *graph, fmatvec::SqrMat &A, int i, std::vector< Element * > &objList)
adds list of objects to tree
Definition: dynamic_system_solver.cc:1551
int solveImpactsRootFinding() override
solve impact equations with Newton scheme on velocity level
Definition: dynamic_system_solver.cc:659
static std::atomic< bool > exitRequest
boolean signal evaluation for end integration set by user
Definition: dynamic_system_solver.h:930
~DynamicSystemSolver() override
destructor
Definition: dynamic_system_solver.cc:96
double tolLocalSolver
Tolerance for local none-linear solvers (solvers on element level)
Definition: dynamic_system_solver.h:874
void updatezdRef(fmatvec::Vec &ext)
references to differentiated external state
Definition: dynamic_system_solver.cc:1416
void addElement(Element *element_)
Definition: dynamic_system_solver.cc:1266
void checkImpactsForTermination() override
validate force laws concerning given tolerances on velocity level
Definition: dynamic_system_solver.cc:735
bool positionDriftCompensationNeeded(double gmax)
check if drift compensation on position level is needed
Definition: dynamic_system_solver.cc:1110
fmatvec::Vec gdParent
relative velocities
Definition: dynamic_system_solver.h:734
DynamicSystemSolver(const std::string &name="")
constructor
Definition: dynamic_system_solver.cc:85
MBSimEnvironment * getMBSimEnvironment()
Definition: dynamic_system_solver.cc:1791
Env * getEnvironment()
Definition: dynamic_system_solver.h:971
fmatvec::Vec zParent
state
Definition: dynamic_system_solver.h:739
bool stopIfNoConvergence
flag if the contact equations should be stopped if there is no convergence
Definition: dynamic_system_solver.h:824
fmatvec::Vec laParent
contact force parameters
Definition: dynamic_system_solver.h:707
void projectGeneralizedPositions(int mode, bool fullUpdate=false)
drift projection for positions
Definition: dynamic_system_solver.cc:1136
fmatvec::SymMat LLMParent
Cholesky decomposition of mass matrix.
Definition: dynamic_system_solver.h:682
fmatvec::Vec rFactorParent
relaxation parameters for contact equations
Definition: dynamic_system_solver.h:712
void init(InitStage stage, const InitConfigSet &config) override
plots time series header
Definition: dynamic_system_solver.cc:120
Element * getElement(const std::string &name)
Definition: dynamic_system_solver.cc:1278
fmatvec::Vec svParent
stopvector (rootfunctions for event driven integration
Definition: dynamic_system_solver.h:763
void constructor()
set plot feature default values
Definition: dynamic_system_solver.cc:1430
double dt
step size
Definition: dynamic_system_solver.h:657
int zSize
size of state vector
Definition: dynamic_system_solver.h:662
fmatvec::Mat TParent
matrix of linear relation between differentiated positions and velocities
Definition: dynamic_system_solver.h:687
void updategdRef(fmatvec::Vec &ref) override
references to relative velocities of dynamic system parent
Definition: dynamic_system_solver.h:589
virtual const fmatvec::Vec & shift(std::optional< std::reference_wrapper< bool > > &&velProjWasCalled={}, std::optional< std::reference_wrapper< bool > > &&posProjWasCalled={})
update for event driven integrator for event
Definition: dynamic_system_solver.cc:1564
virtual void updateG()
updates mass action matrix
Definition: dynamic_system_solver.cc:987
void updateVRef(fmatvec::Mat &ref, int i=0) override
references to condensed contact force direction matrix of dynamic system parent
Definition: dynamic_system_solver.h:609
fmatvec::Vec zdParent
differentiated state
Definition: dynamic_system_solver.h:744
void setDynamicSystemSolverTolerance(double tol)
set tolerance for global none-linear solver (solvers on DynamicSystemSolver level),...
Definition: dynamic_system_solver.h:514
void initla()
load contact forces for use as starting value
Definition: dynamic_system_solver.cc:1251
void readz0(std::string fileName)
reads state from a file
Definition: dynamic_system_solver.cc:1396
fmatvec::SparseMat Gs
sparse mass action matrix
Definition: dynamic_system_solver.h:784
void updatewbRef(fmatvec::Vec &ref) override
references to TODO of dynamic system parent
Definition: dynamic_system_solver.h:595
fmatvec::SymMat MParent
mass matrix
Definition: dynamic_system_solver.h:677
void dropContactMatrices()
Definition: dynamic_system_solver.cc:1322
bool velocityDriftCompensationNeeded(double gdmax)
check if drift compensation on velocity level is needed
Definition: dynamic_system_solver.cc:1123
fmatvec::VecInt jsvParent
boolean evaluation of stopvector
Definition: dynamic_system_solver.h:768
int solveImpactsFixpointSingle() override
solve impact equations with single step fixed point scheme on velocity level
Definition: dynamic_system_solver.cc:522
fmatvec::Vec rParent[2]
nonsmooth right hand side
Definition: dynamic_system_solver.h:756
Solver
solver for contact equations
Definition: dynamic_system_solver.h:76
double tolDSS
Tolerance for global none-linear solvers (solvers on DynamicSystemSolver level)
Definition: dynamic_system_solver.h:879
int maxIter
maximum number of contact iterations, high number of contact iterations for warnings,...
Definition: dynamic_system_solver.h:809
int solveImpactsGaussSeidel() override
solve impact equations with Gauss-Seidel scheme on velocity level
Definition: dynamic_system_solver.cc:572
void projectGeneralizedVelocities(int mode)
drift projection for positions
Definition: dynamic_system_solver.cc:1200
bool checkGSize
boolean if force action matrix should be resized in each step
Definition: dynamic_system_solver.h:854
void initLa()
load contact impulses for use as starting value
Definition: dynamic_system_solver.cc:1261
virtual void updater(int j=0)
update smooth link force law
Definition: dynamic_system_solver.cc:848
void savela()
save contact forces for use as starting value in next time step
Definition: dynamic_system_solver.cc:1246
int limitGSize
TODO.
Definition: dynamic_system_solver.h:859
void decreaserFactors()
decrease relaxation factors if mass action matrix is not diagonal dominant
Definition: dynamic_system_solver.cc:1088
virtual void updatewbRef(fmatvec::Vec &wbParent)
references to TODO of dynamic system parent
Definition: dynamic_system.cc:854
void plot() override
plots time dependent data
Definition: dynamic_system.cc:369
fmatvec::Mat Jrla[2]
Jacobian dr/dla.
Definition: dynamic_system.h:868
fmatvec::Vec la
contact force parameters
Definition: dynamic_system.h:878
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
virtual void updateWRef(fmatvec::Mat &WParent, int j=0)
references to contact force direction matrix of dynamic system parent
Definition: dynamic_system.cc:861
virtual void updategRef(fmatvec::Vec &gParent)
references to relative distances of dynamic system parent
Definition: dynamic_system.cc:819
virtual void updategdRef(fmatvec::Vec &gdParent)
references to relative velocities of dynamic system parent
Definition: dynamic_system.cc:826
Element(const std::string &name)
constructor
Definition: element.cc:63
InitStage
The stages of the initialization.
Definition: element.h:62
std::string name
name of element
Definition: element.h:260
Definition: environment.h:37
Definition: function.h:53
class for tree-structured mechanical systems with recursive and flat memory mechanism
Definition: graph.h:36
group ingredients do not depend on each other
Definition: group.h:35
Definition: environment.h:59
namespace MBSim
Definition: bilateral_constraint.cc:30
Definition: dynamic_system_solver.h:38