timings-contact-dynamics.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2019-2025 INRIA
3 //
4 
5 #include "model-fixture.hpp"
6 
21 
22 #include <benchmark/benchmark.h>
23 
24 #include <iostream>
25 
27 {
28  void SetUp(benchmark::State & st)
29  {
31 
32  const std::string RF = "RLEG_LINK6";
33  const auto RF_id = model.frames[model.getFrameId(RF)].parentJoint;
34  const std::string LF = "LLEG_LINK6";
35  const auto LF_id = model.frames[model.getFrameId(LF)].parentJoint;
36 
37  ci_RF_6D = std::make_unique<pinocchio::RigidConstraintModel>(
39  ci_LF_6D = std::make_unique<pinocchio::RigidConstraintModel>(
41 
43 
44  contact_models_6D.clear();
45  contact_models_6D.push_back(*ci_RF_6D);
46 
47  contact_data_6D.clear();
48  contact_data_6D.push_back(pinocchio::RigidConstraintData(*ci_RF_6D));
49 
51 
52  contact_models_6D6D.clear();
53  contact_models_6D6D.push_back(*ci_RF_6D);
54  contact_models_6D6D.push_back(*ci_LF_6D);
55 
56  contact_data_6D6D.clear();
57  contact_data_6D6D.push_back(pinocchio::RigidConstraintData(*ci_RF_6D));
58  contact_data_6D6D.push_back(pinocchio::RigidConstraintData(*ci_LF_6D));
59 
61 
63  prox_settings.mu = 1e8;
64 
65  num_constraints = 12;
66 
68  + Eigen::MatrixXd::Identity(num_constraints, num_constraints);
69  }
70 
71  void TearDown(benchmark::State & st)
72  {
74  }
75 
76  std::unique_ptr<pinocchio::RigidConstraintModel> ci_RF_6D;
77  std::unique_ptr<pinocchio::RigidConstraintModel> ci_LF_6D;
78 
85 
89 
91 
92  int num_constraints = 12;
93 
94  Eigen::MatrixXd col_major_square_matrices;
95 };
96 
97 static void CustomArguments(benchmark::internal::Benchmark * b)
98 {
99  b->MinWarmUpTime(3.);
100 }
101 
102 // ABA_WORLD
103 
105  const pinocchio::Model & model,
107  const Eigen::VectorXd & q,
108  const Eigen::VectorXd & v,
109  const Eigen::VectorXd & tau)
110 {
112 }
113 BENCHMARK_DEFINE_F(ContactFixture, ABA_WORLD)(benchmark::State & st)
114 {
115  for (auto _ : st)
116  {
117  abaWorldCall(model, data, q, v, tau);
118  }
119 }
121 
122 // CHOLESKY_DECOMPOSE
123 
124 PINOCCHIO_DONT_INLINE static void
126 {
128 }
129 BENCHMARK_DEFINE_F(ContactFixture, CHOLESKY_DECOMPOSE)(benchmark::State & st)
130 {
132  for (auto _ : st)
133  {
135  }
136 }
137 BENCHMARK_REGISTER_F(ContactFixture, CHOLESKY_DECOMPOSE)->Apply(CustomArguments);
138 
139 // CONTACT_ABA_EMPTY
140 
142  const pinocchio::Model & model,
144  const Eigen::VectorXd & q,
145  const Eigen::VectorXd & v,
146  const Eigen::VectorXd & tau,
149 {
151 }
152 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_ABA_EMPTY)(benchmark::State & st)
153 {
154  for (auto _ : st)
155  {
156  contactABACall(model, data, q, v, tau, contact_models_empty, contact_data_empty);
157  }
158 }
159 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_ABA_EMPTY)->Apply(CustomArguments);
160 
161 // PV_EMPTY
162 
164  const pinocchio::Model & model,
166  const Eigen::VectorXd & q,
167  const Eigen::VectorXd & v,
168  const Eigen::VectorXd & tau,
172 {
174 }
175 BENCHMARK_DEFINE_F(ContactFixture, PV_EMPTY)(benchmark::State & st)
176 {
177  pinocchio::initPvSolver(model, data, contact_models_empty);
178  for (auto _ : st)
179  {
180  pvCall(model, data, q, v, tau, contact_models_empty, contact_data_empty, prox_settings);
181  }
182 }
184 
185 // CONSTRAINED_ABA_EMPTY
186 
188  const pinocchio::Model & model,
190  const Eigen::VectorXd & q,
191  const Eigen::VectorXd & v,
192  const Eigen::VectorXd & tau,
196 {
198 }
199 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINED_ABA_EMPTY)(benchmark::State & st)
200 {
201  pinocchio::initPvSolver(model, data, contact_models_empty);
202  for (auto _ : st)
203  {
205  model, data, q, v, tau, contact_models_empty, contact_data_empty, prox_settings);
206  }
207 }
208 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINED_ABA_EMPTY)->Apply(CustomArguments);
209 
210 // CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_EMPTY
211 
214  const pinocchio::Model & model,
218 {
219  contact.compute(model, data, contact_models, contact_data);
220 }
221 
222 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_EMPTY)(
223  benchmark::State & st)
224 {
226  for (auto _ : st)
227  {
229  contact_chol_empty, model, data, contact_models_empty, contact_data_empty);
230  }
231 }
232 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_EMPTY)
233  ->Apply(CustomArguments);
234 
235 // CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_EMPTY
236 
238  pinocchio::ContactCholeskyDecomposition & contact, Eigen::MatrixXd & H_inverse)
239 {
240  contact.inverse(H_inverse);
241 }
242 
243 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_EMPTY)(
244  benchmark::State & st)
245 {
247  contact_chol_empty.compute(model, data, contact_models_empty, contact_data_empty);
248  Eigen::MatrixXd H_inverse(contact_chol_empty.size(), contact_chol_empty.size());
249  for (auto _ : st)
250  {
251  contactCholeskyDecompositionInverseCall(contact_chol_empty, H_inverse);
252  }
253 }
254 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_EMPTY)
255  ->Apply(CustomArguments);
256 
257 // CONSTRAINT_DYNAMICS_EMPTY
258 
260  const pinocchio::Model & model,
262  const Eigen::VectorXd & q,
263  const Eigen::VectorXd & v,
264  const Eigen::VectorXd & tau,
267 {
269 }
270 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINT_DYNAMICS_EMPTY)(benchmark::State & st)
271 {
272  pinocchio::initConstraintDynamics(model, data, contact_models_empty);
273  for (auto _ : st)
274  {
275  constraintDynamicsCall(model, data, q, v, tau, contact_models_empty, contact_data_empty);
276  }
277 }
278 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINT_DYNAMICS_EMPTY)->Apply(CustomArguments);
279 
280 // CONTACT_ABA_6D
281 
282 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_ABA_6D)(benchmark::State & st)
283 {
284  for (auto _ : st)
285  {
286  contactABACall(model, data, q, v, tau, contact_models_6D, contact_data_6D);
287  }
288 }
289 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_ABA_6D)->Apply(CustomArguments);
290 
291 // PV_6D
292 
293 BENCHMARK_DEFINE_F(ContactFixture, PV_6D)(benchmark::State & st)
294 {
295  pinocchio::initPvSolver(model, data, contact_models_6D);
296  for (auto _ : st)
297  {
298  pvCall(model, data, q, v, tau, contact_models_6D, contact_data_6D, prox_settings);
299  }
300 }
302 
303 // CONSTRAINED_ABA_6D
304 
305 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINED_ABA_6D)(benchmark::State & st)
306 {
307  pinocchio::initPvSolver(model, data, contact_models_6D);
308  for (auto _ : st)
309  {
310  constrainedABACall(model, data, q, v, tau, contact_models_6D, contact_data_6D, prox_settings);
311  }
312 }
313 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINED_ABA_6D)->Apply(CustomArguments);
314 
315 // CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_6D
316 
317 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_6D)(benchmark::State & st)
318 {
320  for (auto _ : st)
321  {
323  contact_chol_6D, model, data, contact_models_6D, contact_data_6D);
324  }
325 }
326 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_6D)
327  ->Apply(CustomArguments);
328 
329 // CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_6D
330 
331 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_6D)(benchmark::State & st)
332 {
334  contact_chol_6D.compute(model, data, contact_models_6D, contact_data_6D);
335  Eigen::MatrixXd H_inverse(contact_chol_6D.size(), contact_chol_6D.size());
336  for (auto _ : st)
337  {
338  contactCholeskyDecompositionInverseCall(contact_chol_6D, H_inverse);
339  }
340 }
341 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_6D)
342  ->Apply(CustomArguments);
343 
344 // CONSTRAINT_DYNAMICS_6D
345 
346 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINT_DYNAMICS_6D)(benchmark::State & st)
347 {
348  pinocchio::initConstraintDynamics(model, data, contact_models_6D);
349  for (auto _ : st)
350  {
351  constraintDynamicsCall(model, data, q, v, tau, contact_models_6D, contact_data_6D);
352  }
353 }
354 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINT_DYNAMICS_6D)->Apply(CustomArguments);
355 
356 // GET_KKT_CONTACT_DYNAMIC_MATRIX_INVERSE_6D
357 
359  const pinocchio::Model & model,
361  const Eigen::MatrixXd & J,
362  const Eigen::MatrixXd & MJtJ_inv)
363 {
366 }
367 
368 BENCHMARK_DEFINE_F(ContactFixture, GET_KKT_CONTACT_DYNAMIC_MATRIX_INVERSE_6D)(benchmark::State & st)
369 {
370  Eigen::MatrixXd J(contact_chol_6D.constraintDim(), model.nv);
371  J.setZero();
372 
373  Eigen::MatrixXd MJtJ_inv(
374  model.nv + contact_chol_6D.constraintDim(), model.nv + contact_chol_6D.constraintDim());
375  MJtJ_inv.setZero();
376 
377  Eigen::VectorXd gamma(contact_chol_6D.constraintDim());
378  gamma.setZero();
379 
382  model, data, ci_RF_6D->joint1_id, ci_RF_6D->reference_frame, J.middleRows<6>(0));
384 
385  for (auto _ : st)
386  {
388  }
389 }
390 BENCHMARK_REGISTER_F(ContactFixture, GET_KKT_CONTACT_DYNAMIC_MATRIX_INVERSE_6D)
391  ->Apply(CustomArguments);
392 
393 // CONTACT_ABA_6D6D
394 
395 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_ABA_6D6D)(benchmark::State & st)
396 {
397  for (auto _ : st)
398  {
399  contactABACall(model, data, q, v, tau, contact_models_6D6D, contact_data_6D6D);
400  }
401 }
402 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_ABA_6D6D)->Apply(CustomArguments);
403 
404 // PV_6D6D
405 
406 BENCHMARK_DEFINE_F(ContactFixture, PV_6D6D)(benchmark::State & st)
407 {
408  pinocchio::initPvSolver(model, data, contact_models_6D6D);
409  for (auto _ : st)
410  {
411  pvCall(model, data, q, v, tau, contact_models_6D6D, contact_data_6D6D, prox_settings);
412  }
413 }
415 
416 // CONSTRAINED_ABA_6D6D
417 
418 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINED_ABA_6D6D)(benchmark::State & st)
419 {
420  pinocchio::initPvSolver(model, data, contact_models_6D6D);
421  for (auto _ : st)
422  {
424  model, data, q, v, tau, contact_models_6D6D, contact_data_6D6D, prox_settings);
425  }
426 }
427 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINED_ABA_6D6D)->Apply(CustomArguments);
428 
429 // CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_6D6D
430 
431 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_6D6D)(
432  benchmark::State & st)
433 {
435  for (auto _ : st)
436  {
438  contact_chol_6D6D, model, data, contact_models_6D6D, contact_data_6D6D);
439  }
440 }
441 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_COMPUTE_6D6D)
442  ->Apply(CustomArguments);
443 
444 // CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_6D6D
445 
446 BENCHMARK_DEFINE_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_6D6D)(
447  benchmark::State & st)
448 {
450  contact_chol_6D6D.compute(model, data, contact_models_6D6D, contact_data_6D6D);
451  Eigen::MatrixXd H_inverse(contact_chol_6D6D.size(), contact_chol_6D6D.size());
452  for (auto _ : st)
453  {
454  contactCholeskyDecompositionInverseCall(contact_chol_6D6D, H_inverse);
455  }
456 }
457 BENCHMARK_REGISTER_F(ContactFixture, CONTACT_CHOLESKY_DECOMPOSITION_INVERSE_6D6D)
458  ->Apply(CustomArguments);
459 
460 // CONSTRAINT_DYNAMICS_6D6D
461 
462 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINT_DYNAMICS_6D6D)(benchmark::State & st)
463 {
464  pinocchio::initConstraintDynamics(model, data, contact_models_6D6D);
465  for (auto _ : st)
466  {
467  constraintDynamicsCall(model, data, q, v, tau, contact_models_6D6D, contact_data_6D6D);
468  }
469 }
470 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINT_DYNAMICS_6D6D)->Apply(CustomArguments);
471 
472 // GET_KKT_CONTACT_DYNAMIC_MATRIX_INVERSE_6D6D
473 
474 BENCHMARK_DEFINE_F(ContactFixture, GET_KKT_CONTACT_DYNAMIC_MATRIX_INVERSE_6D6D)(
475  benchmark::State & st)
476 {
477  Eigen::MatrixXd J(contact_chol_6D6D.constraintDim(), model.nv);
478  J.setZero();
479 
480  Eigen::MatrixXd MJtJ_inv(
481  model.nv + contact_chol_6D6D.constraintDim(), model.nv + contact_chol_6D6D.constraintDim());
482  MJtJ_inv.setZero();
483 
484  Eigen::VectorXd gamma(contact_chol_6D6D.constraintDim());
485  gamma.setZero();
486 
488  getJointJacobian(model, data, ci_RF_6D->joint1_id, ci_RF_6D->reference_frame, J.middleRows<6>(0));
489  getJointJacobian(model, data, ci_LF_6D->joint1_id, ci_LF_6D->reference_frame, J.middleRows<6>(6));
490  forwardDynamics(model, data, q, v, tau, J, gamma);
491 
492  for (auto _ : st)
493  {
495  }
496 }
497 BENCHMARK_REGISTER_F(ContactFixture, GET_KKT_CONTACT_DYNAMIC_MATRIX_INVERSE_6D6D)
498  ->Apply(CustomArguments);
499 
500 // FORWARD_DYNAMICS_6D6D
501 
503  const pinocchio::Model & model,
507  const Eigen::VectorXd & q,
508  const Eigen::VectorXd & v,
509  const Eigen::VectorXd & tau,
510  Eigen::MatrixXd & J,
511  const Eigen::VectorXd & gamma)
512 {
514  pinocchio::getJointJacobian(model, data, c1.joint1_id, c1.reference_frame, J.middleRows<6>(0));
515  pinocchio::getJointJacobian(model, data, c2.joint1_id, c2.reference_frame, J.middleRows<6>(6));
517 }
518 
519 BENCHMARK_DEFINE_F(ContactFixture, FORWARD_DYNAMICS_6D6D)(benchmark::State & st)
520 {
521  Eigen::MatrixXd J(contact_chol_6D6D.constraintDim(), model.nv);
522  J.setZero();
523 
524  Eigen::VectorXd gamma(contact_chol_6D6D.constraintDim());
525  gamma.setZero();
526 
527  initConstraintDynamics(model, data, contact_models_6D6D);
528  for (auto _ : st)
529  {
530  forwardDynamisCall(model, data, *ci_RF_6D, *ci_LF_6D, q, v, tau, J, gamma);
531  }
532 }
533 BENCHMARK_REGISTER_F(ContactFixture, FORWARD_DYNAMICS_6D6D)->Apply(CustomArguments);
534 
BENCHMARK_REGISTER_F
BENCHMARK_REGISTER_F(ContactFixture, ABA_WORLD) -> Apply(CustomArguments)
abaWorldCall
static PINOCCHIO_DONT_INLINE void abaWorldCall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau)
Definition: timings-contact-dynamics.cpp:104
ContactFixture::col_major_square_matrices
Eigen::MatrixXd col_major_square_matrices
Definition: timings-contact-dynamics.cpp:94
pinocchio::getKKTContactDynamicMatrixInverse
PINOCCHIO_DEPRECATED void getKKTContactDynamicMatrixInverse(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Eigen::MatrixBase< KKTMatrixType > &KKTMatrix_inv)
Computes the inverse of the KKT matrix for dynamics with contact constraints.
constrainedABACall
static PINOCCHIO_DONT_INLINE void constrainedABACall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau, const PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintModel) &contact_models, PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintData) &contact_data, pinocchio::ProximalSettings &prox_settings)
Definition: timings-contact-dynamics.cpp:187
pinocchio::ModelTpl::frames
FrameVector frames
Vector of operational frames registered on the model.
Definition: multibody/model.hpp:193
pinocchio::constrainedABA
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & constrainedABA(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ContactModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, ContactDataAllocator > &contact_datas, ProximalSettingsTpl< Scalar > &settings)
The constrained Articulated Body Algorithm (constrainedABA). It computes constrained forward dynamics...
pinocchio::DataTpl
Definition: context/generic.hpp:25
PINOCCHIO_DONT_INLINE
#define PINOCCHIO_DONT_INLINE
Function attribute to forbid inlining. This is a compiler hint that can be not respected.
Definition: include/pinocchio/macros.hpp:53
forwardDynamisCall
static PINOCCHIO_DONT_INLINE void forwardDynamisCall(const pinocchio::Model &model, pinocchio::Data &data, const pinocchio::RigidConstraintModel &c1, const pinocchio::RigidConstraintModel &c2, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau, Eigen::MatrixXd &J, const Eigen::VectorXd &gamma)
Definition: timings-contact-dynamics.cpp:502
pinocchio::ProximalSettingsTpl::mu
Scalar mu
Regularization parameter of the proximal algorithm.
Definition: algorithm/proximal.hpp:95
ModelFixture::model
pinocchio::Model model
Definition: model-fixture.hpp:79
choleskyDecomposeCall
static PINOCCHIO_DONT_INLINE void choleskyDecomposeCall(const pinocchio::Model &model, pinocchio::Data &data)
Definition: timings-contact-dynamics.cpp:125
ModelFixture::TearDown
void TearDown(benchmark::State &)
Definition: model-fixture.hpp:51
ModelFixture
Definition: model-fixture.hpp:37
ContactFixture::contact_chol_6D
pinocchio::ContactCholeskyDecomposition contact_chol_6D
Definition: timings-constrained-dynamics-derivatives.cpp:73
pinocchio::pv
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & pv(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ContactModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, ContactDataAllocator > &contact_datas, ProximalSettingsTpl< Scalar > &settings)
The Popov-Vereshchagin algorithm. It computes constrained forward dynamics, aka the joint acceleratio...
pinocchio::Convention::WORLD
@ WORLD
contact-cholesky.contact_data
list contact_data
Definition: contact-cholesky.py:33
kinematics.hpp
rnea-derivatives.hpp
setup.data
data
Definition: cmake/cython/setup.in.py:48
ContactFixture::TearDown
void TearDown(benchmark::State &st)
Definition: timings-contact-dynamics.cpp:71
inverse-kinematics-3d.J
J
Definition: inverse-kinematics-3d.py:28
pinocchio::crba
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & crba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Convention convention=Convention::LOCAL)
Computes the upper triangular part of the joint space inertia matrix M by using the Composite Rigid B...
pinocchio::aba
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & aba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const Convention convention=Convention::LOCAL)
The Articulated-Body algorithm. It computes the forward dynamics, aka the joint accelerations given t...
pinocchio::initPvSolver
void initPvSolver(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const std::vector< RigidConstraintModelTpl< Scalar, Options >, Allocator > &contact_models)
Init the data according to the contact information contained in contact_models.
pinocchio::CONTACT_6D
@ CONTACT_6D
Point contact model.
Definition: algorithm/contact-info.hpp:22
contactABACall
static PINOCCHIO_DONT_INLINE void contactABACall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau, const PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintModel) &contact_models, PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintData) &contact_data)
Definition: timings-contact-dynamics.cpp:141
rnea.hpp
aba-derivatives.hpp
b
Vec3f b
ModelFixture::SetUp
void SetUp(benchmark::State &)
Definition: model-fixture.hpp:39
aba.hpp
anymal-simulation.model
model
Definition: anymal-simulation.py:8
contactCholeskyDecompositionComputeCall
static PINOCCHIO_DONT_INLINE void contactCholeskyDecompositionComputeCall(pinocchio::ContactCholeskyDecomposition &contact, const pinocchio::Model &model, pinocchio::Data &data, const PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintModel) &contact_models, PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintData) &contact_data)
Definition: timings-contact-dynamics.cpp:212
forward-dynamics-derivatives.tau
tau
Definition: forward-dynamics-derivatives.py:25
joint-configuration.hpp
pinocchio::ProximalSettingsTpl
Structure containing all the settings parameters for the proximal algorithms.
Definition: algorithm/fwd.hpp:13
ContactFixture::contact_chol_empty
pinocchio::ContactCholeskyDecomposition contact_chol_empty
Definition: timings-constrained-dynamics-derivatives.cpp:72
pinocchio::RigidConstraintModelTpl
Definition: algorithm/constraints/fwd.hpp:14
ContactFixture::ci_RF_6D
std::unique_ptr< pinocchio::RigidConstraintModel > ci_RF_6D
Definition: timings-constrained-dynamics-derivatives.cpp:62
model-fixture.hpp
pinocchio::getJointJacobian
void getJointJacobian(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const DataTpl< Scalar, Options, JointCollectionTpl > &data, const JointIndex joint_id, const ReferenceFrame reference_frame, const Eigen::MatrixBase< Matrix6Like > &J)
Computes the Jacobian of a specific joint frame expressed in one of the pinocchio::ReferenceFrame opt...
pinocchio::ModelTpl::getFrameId
FrameIndex getFrameId(const std::string &name, const FrameType &type=(FrameType)(JOINT|FIXED_JOINT|BODY|OP_FRAME|SENSOR)) const
Returns the index of a frame given by its name.
urdf.hpp
pinocchio::forwardDynamics
const PINOCCHIO_DEPRECATED DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & forwardDynamics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const Eigen::MatrixBase< ConstraintMatrixType > &J, const Eigen::MatrixBase< DriftVectorType > &gamma, const Scalar inv_damping=0.)
Compute the forward dynamics with contact constraints. Internally, pinocchio::computeAllTerms is call...
pinocchio::RigidConstraintDataTpl
Definition: algorithm/constraints/fwd.hpp:16
cartpole.v
v
Definition: cartpole.py:153
pinocchio::initConstraintDynamics
void initConstraintDynamics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const std::vector< RigidConstraintModelTpl< Scalar, Options >, Allocator > &contact_models)
Init the forward dynamics data according to the contact information contained in contact_models.
inverse-dynamics._
_
Definition: inverse-dynamics.py:22
ContactFixture
Definition: timings-constrained-dynamics-derivatives.cpp:20
q
q
PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR
#define PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(T)
Definition: container/aligned-vector.hpp:13
getKKTContactDynamicMatrixInverseCall
static PINOCCHIO_DONT_INLINE void getKKTContactDynamicMatrixInverseCall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::MatrixXd &J, const Eigen::MatrixXd &MJtJ_inv)
Definition: timings-contact-dynamics.cpp:358
cholesky.hpp
pinocchio::cholesky::decompose
const DataTpl< Scalar, Options, JointCollectionTpl >::MatrixXs & decompose(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data)
Compute the Cholesky decomposition of the joint space inertia matrix M contained in data.
CustomArguments
static void CustomArguments(benchmark::internal::Benchmark *b)
Definition: timings-contact-dynamics.cpp:97
contactCholeskyDecompositionInverseCall
static PINOCCHIO_DONT_INLINE void contactCholeskyDecompositionInverseCall(pinocchio::ContactCholeskyDecomposition &contact, Eigen::MatrixXd &H_inverse)
Definition: timings-contact-dynamics.cpp:237
pinocchio::contactABA
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & contactABA(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, DataAllocator > &contact_data)
Definition: constrained-dynamics.hpp:172
ContactFixture::prox_settings
pinocchio::ProximalSettings prox_settings
Definition: timings-contact-dynamics.cpp:90
pinocchio::ContactCholeskyDecomposition
ContactCholeskyDecompositionTpl< context::Scalar, context::Options > ContactCholeskyDecomposition
Definition: algorithm/fwd.hpp:17
pinocchio::ProximalSettingsTpl::max_iter
int max_iter
Maximal number of iterations.
Definition: algorithm/proximal.hpp:98
PINOCCHIO_BENCHMARK_MAIN
PINOCCHIO_BENCHMARK_MAIN()
kinematics-derivatives.hpp
pvCall
static PINOCCHIO_DONT_INLINE void pvCall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau, const PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintModel) &contact_models, PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintData) &contact_data, pinocchio::ProximalSettings &prox_settings)
Definition: timings-contact-dynamics.cpp:163
ContactFixture::contact_chol_6D6D
pinocchio::ContactCholeskyDecomposition contact_chol_6D6D
Definition: timings-constrained-dynamics-derivatives.cpp:74
BENCHMARK_DEFINE_F
BENCHMARK_DEFINE_F(ContactFixture, ABA_WORLD)(benchmark
Definition: timings-contact-dynamics.cpp:113
contact-cholesky.contact_models
list contact_models
Definition: contact-cholesky.py:22
pinocchio::ContactCholeskyDecompositionTpl< context::Scalar, context::Options >
contact-dynamics.hpp
sample-models.hpp
c2
c2
pinocchio::ModelTpl< context::Scalar, context::Options >
ContactFixture::PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR
PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintModel) contact_models_empty
constrained-dynamics.hpp
simulation-closed-kinematic-chains.prox_settings
prox_settings
Definition: simulation-closed-kinematic-chains.py:158
ContactFixture::SetUp
void SetUp(benchmark::State &st)
Definition: timings-contact-dynamics.cpp:28
crba.hpp
pinocchio::computeAllTerms
void computeAllTerms(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes efficiently all the terms needed for dynamic simulation. It is equivalent to the call at the...
ContactFixture::ci_LF_6D
std::unique_ptr< pinocchio::RigidConstraintModel > ci_LF_6D
Definition: timings-constrained-dynamics-derivatives.cpp:63
pinocchio::constraintDynamics
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & constraintDynamics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType1 > &v, const Eigen::MatrixBase< TangentVectorType2 > &tau, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ConstraintModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, ConstraintDataAllocator > &contact_datas, ProximalSettingsTpl< Scalar > &settings)
Computes the forward dynamics with contact constraints according to a given list of contact informati...
ContactFixture::num_constraints
int num_constraints
Definition: timings-contact-dynamics.cpp:92
pv.hpp
constraintDynamicsCall
static PINOCCHIO_DONT_INLINE void constraintDynamicsCall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau, const PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintModel) &contact_models, PINOCCHIO_STD_VECTOR_WITH_EIGEN_ALLOCATOR(pinocchio::RigidConstraintData) &contact_data)
Definition: timings-contact-dynamics.cpp:259
pinocchio::LOCAL
@ LOCAL
Definition: multibody/fwd.hpp:50


pinocchio
Author(s):
autogenerated on Wed Apr 16 2025 02:41:51