timings-constrained-dynamics-derivatives.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2019-2025 LAAS-CNRS INRIA
3 //
4 //
5 #include "model-fixture.hpp"
6 
15 
16 #include <benchmark/benchmark.h>
17 
18 #include <iostream>
19 
21 {
22  void SetUp(benchmark::State & st)
23  {
25 
26  const std::string RF = "RLEG_LINK6";
27  const auto RF_id = model.frames[model.getFrameId(RF)].parentJoint;
28  const std::string LF = "LLEG_LINK6";
29  const auto LF_id = model.frames[model.getFrameId(LF)].parentJoint;
30 
31  ci_RF_6D = std::make_unique<pinocchio::RigidConstraintModel>(
33  ci_LF_6D = std::make_unique<pinocchio::RigidConstraintModel>(
35 
37 
38  contact_models_6D.clear();
39  contact_models_6D.push_back(*ci_RF_6D);
40 
41  contact_data_6D.clear();
42  contact_data_6D.push_back(pinocchio::RigidConstraintData(*ci_RF_6D));
43 
45 
46  contact_models_6D6D.clear();
47  contact_models_6D6D.push_back(*ci_RF_6D);
48  contact_models_6D6D.push_back(*ci_LF_6D);
49 
50  contact_data_6D6D.clear();
51  contact_data_6D6D.push_back(pinocchio::RigidConstraintData(*ci_RF_6D));
52  contact_data_6D6D.push_back(pinocchio::RigidConstraintData(*ci_LF_6D));
53 
55  }
56 
57  void TearDown(benchmark::State & st)
58  {
60  }
61 
62  std::unique_ptr<pinocchio::RigidConstraintModel> ci_RF_6D;
63  std::unique_ptr<pinocchio::RigidConstraintModel> ci_LF_6D;
64 
71 
75 };
76 
77 static void CustomArguments(benchmark::internal::Benchmark * b)
78 {
79  b->MinWarmUpTime(3.);
80 }
81 
82 // COMPUTE_ABA_DERIVATIVES
83 
85  const pinocchio::Model & model,
87  const Eigen::VectorXd & q,
88  const Eigen::VectorXd & v,
89  const Eigen::VectorXd & tau)
90 {
92 }
93 BENCHMARK_DEFINE_F(ContactFixture, COMPUTE_ABA_DERIVATIVES)(benchmark::State & st)
94 {
95  for (auto _ : st)
96  {
98  }
99 }
100 BENCHMARK_REGISTER_F(ContactFixture, COMPUTE_ABA_DERIVATIVES)->Apply(CustomArguments);
101 
102 // CONSTRAINT_DYNAMICS_DERIVATIVES_EMPTY
103 
105  const pinocchio::Model & model,
109 {
111 }
112 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINT_DYNAMICS_DERIVATIVES_EMPTY)(benchmark::State & st)
113 {
114  pinocchio::initConstraintDynamics(model, data, contact_models_empty);
115  pinocchio::constraintDynamics(model, data, q, v, tau, contact_models_empty, contact_data_empty);
116  for (auto _ : st)
117  {
118  constraintDynamicsDerivativesCall(model, data, contact_models_empty, contact_data_empty);
119  }
120 }
121 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINT_DYNAMICS_DERIVATIVES_EMPTY)->Apply(CustomArguments);
122 
123 // CONSTRAINT_DYNAMICS_DERIVATIVES_6D
124 
125 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINT_DYNAMICS_DERIVATIVES_6D)(benchmark::State & st)
126 {
127  pinocchio::initConstraintDynamics(model, data, contact_models_6D);
128  pinocchio::constraintDynamics(model, data, q, v, tau, contact_models_6D, contact_data_6D);
129  for (auto _ : st)
130  {
131  constraintDynamicsDerivativesCall(model, data, contact_models_6D, contact_data_6D);
132  }
133 }
134 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINT_DYNAMICS_DERIVATIVES_6D)->Apply(CustomArguments);
135 
136 // CONSTRAINT_DYNAMICS_DERIVATIVES_6D6D
137 
138 BENCHMARK_DEFINE_F(ContactFixture, CONSTRAINT_DYNAMICS_DERIVATIVES_6D6D)(benchmark::State & st)
139 {
140  pinocchio::initConstraintDynamics(model, data, contact_models_6D6D);
141  pinocchio::constraintDynamics(model, data, q, v, tau, contact_models_6D6D, contact_data_6D6D);
142  for (auto _ : st)
143  {
144  constraintDynamicsDerivativesCall(model, data, contact_models_6D6D, contact_data_6D6D);
145  }
146 }
147 BENCHMARK_REGISTER_F(ContactFixture, CONSTRAINT_DYNAMICS_DERIVATIVES_6D6D)->Apply(CustomArguments);
148 
pinocchio::ModelTpl::frames
FrameVector frames
Vector of operational frames registered on the model.
Definition: multibody/model.hpp:193
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
CustomArguments
static void CustomArguments(benchmark::internal::Benchmark *b)
Definition: timings-constrained-dynamics-derivatives.cpp:77
pinocchio::computeConstraintDynamicsDerivatives
void computeConstraintDynamicsDerivatives(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const std::vector< RigidConstraintModelTpl< Scalar, Options >, ConstraintModelAllocator > &contact_models, std::vector< RigidConstraintDataTpl< Scalar, Options >, ConstraintDataAllocator > &contact_data, const ProximalSettingsTpl< Scalar > &settings, const Eigen::MatrixBase< MatrixType1 > &ddq_partial_dq, const Eigen::MatrixBase< MatrixType2 > &ddq_partial_dv, const Eigen::MatrixBase< MatrixType3 > &ddq_partial_dtau, const Eigen::MatrixBase< MatrixType4 > &lambda_partial_dq, const Eigen::MatrixBase< MatrixType5 > &lambda_partial_dv, const Eigen::MatrixBase< MatrixType6 > &lambda_partial_dtau)
ModelFixture::model
pinocchio::Model model
Definition: model-fixture.hpp:79
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
contact-cholesky.contact_data
list contact_data
Definition: contact-cholesky.py:33
setup.data
data
Definition: cmake/cython/setup.in.py:48
BENCHMARK_DEFINE_F
BENCHMARK_DEFINE_F(ContactFixture, COMPUTE_ABA_DERIVATIVES)(benchmark
Definition: timings-constrained-dynamics-derivatives.cpp:93
ContactFixture::TearDown
void TearDown(benchmark::State &st)
Definition: timings-constrained-dynamics-derivatives.cpp:57
constrained-dynamics-derivatives.hpp
pinocchio::CONTACT_6D
@ CONTACT_6D
Point contact model.
Definition: algorithm/contact-info.hpp:22
aba-derivatives.hpp
b
Vec3f b
ModelFixture::SetUp
void SetUp(benchmark::State &)
Definition: model-fixture.hpp:39
anymal-simulation.model
model
Definition: anymal-simulation.py:8
forward-dynamics-derivatives.tau
tau
Definition: forward-dynamics-derivatives.py:25
joint-configuration.hpp
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::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::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.
computeABADerivativesCall
static PINOCCHIO_DONT_INLINE void computeABADerivativesCall(const pinocchio::Model &model, pinocchio::Data &data, const Eigen::VectorXd &q, const Eigen::VectorXd &v, const Eigen::VectorXd &tau)
Definition: timings-constrained-dynamics-derivatives.cpp:84
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
cholesky.hpp
pinocchio::computeABADerivatives
void computeABADerivatives(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< MatrixType1 > &aba_partial_dq, const Eigen::MatrixBase< MatrixType2 > &aba_partial_dv, const Eigen::MatrixBase< MatrixType3 > &aba_partial_dtau)
The derivatives of the Articulated-Body algorithm.
BENCHMARK_REGISTER_F
BENCHMARK_REGISTER_F(ContactFixture, COMPUTE_ABA_DERIVATIVES) -> Apply(CustomArguments)
pinocchio::ContactCholeskyDecomposition
ContactCholeskyDecompositionTpl< context::Scalar, context::Options > ContactCholeskyDecomposition
Definition: algorithm/fwd.hpp:17
ContactFixture::contact_chol_6D6D
pinocchio::ContactCholeskyDecomposition contact_chol_6D6D
Definition: timings-constrained-dynamics-derivatives.cpp:74
contact-cholesky.contact_models
list contact_models
Definition: contact-cholesky.py:22
constraintDynamicsDerivativesCall
static PINOCCHIO_DONT_INLINE void constraintDynamicsDerivativesCall(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-constrained-dynamics-derivatives.cpp:104
pinocchio::ContactCholeskyDecompositionTpl< context::Scalar, context::Options >
contact-dynamics.hpp
sample-models.hpp
PINOCCHIO_BENCHMARK_MAIN
PINOCCHIO_BENCHMARK_MAIN()
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
ContactFixture::SetUp
void SetUp(benchmark::State &st)
Definition: timings-constrained-dynamics-derivatives.cpp:22
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...
pinocchio::LOCAL
@ LOCAL
Definition: multibody/fwd.hpp:50


pinocchio
Author(s):
autogenerated on Thu Apr 10 2025 02:42:22