unittest/crba.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2015-2023 CNRS INRIA
3 //
4 
5 /*
6  * Test the CRBA algorithm. The code validates both the computation times and
7  * the value by comparing the results of the CRBA with the reconstruction of
8  * the mass matrix using the RNEA.
9  * For a strong timing benchmark, see benchmark/timings.
10  *
11  */
12 
13 #ifndef EIGEN_RUNTIME_NO_MALLOC
14  #define EIGEN_RUNTIME_NO_MALLOC
15 #endif
16 
27 
28 #include <iostream>
29 
30 #include <boost/test/unit_test.hpp>
31 #include <boost/utility/binary.hpp>
32 
33 template<typename JointModel>
34 static void addJointAndBody(
37  const std::string & parent_name,
38  const std::string & name,
40  bool setRandomLimits = true)
41 {
42  using namespace pinocchio;
43  typedef typename JointModel::ConfigVector_t CV;
44  typedef typename JointModel::TangentVector_t TV;
45 
47 
48  if (setRandomLimits)
49  idx = model.addJoint(
50  model.getJointId(parent_name), joint, SE3::Random(), name + "_joint",
51  TV::Random() + TV::Constant(1), TV::Random() + TV::Constant(1),
52  CV::Random() - CV::Constant(1), CV::Random() + CV::Constant(1));
53  else
54  idx = model.addJoint(model.getJointId(parent_name), joint, placement, name + "_joint");
55 
56  model.addJointFrame(idx);
57 
58  model.appendBodyToJoint(idx, Inertia::Random(), SE3::Identity());
59  model.addBodyFrame(name + "_body", idx);
60 }
61 
62 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
63 
65 {
69 
70 #ifdef NDEBUG
71  #ifdef _INTENSE_TESTING_
72  const size_t NBT = 1000 * 1000;
73  #else
74  const size_t NBT = 10;
75  #endif
76 
77  Eigen::VectorXd q = Eigen::VectorXd::Zero(model.nq);
78 
80  timer.tic();
81  SMOOTH(NBT)
82  {
84  }
85  timer.toc(std::cout, NBT);
86 
87 #else
88  const size_t NBT = 1;
89  using namespace Eigen;
90  using namespace pinocchio;
91 
92  Eigen::MatrixXd M(model.nv, model.nv);
93  Eigen::VectorXd q = Eigen::VectorXd::Ones(model.nq);
94  q.segment<4>(3).normalize();
95  Eigen::VectorXd v = Eigen::VectorXd::Zero(model.nv);
96  Eigen::VectorXd a = Eigen::VectorXd::Zero(model.nv);
97  data.M.fill(0);
99  data.M.triangularView<Eigen::StrictlyLower>() =
100  data.M.transpose().triangularView<Eigen::StrictlyLower>();
101 
102  /* Joint inertia from iterative crba. */
103  const Eigen::VectorXd bias = rnea(model, data, q, v, a);
104  for (int i = 0; i < model.nv; ++i)
105  {
106  M.col(i) = rnea(model, data, q, v, Eigen::VectorXd::Unit(model.nv, i)) - bias;
107  }
108 
109  // std::cout << "Mcrb = [ " << data.M << " ];" << std::endl;
110  // std::cout << "Mrne = [ " << M << " ]; " << std::endl;
111  BOOST_CHECK(M.isApprox(data.M, 1e-12));
112 
113  std::cout << "(the time score in debug mode is not relevant) ";
114 
115  q = Eigen::VectorXd::Zero(model.nq);
116 
118  timer.tic();
119  SMOOTH(NBT)
120  {
122  }
123  timer.toc(std::cout, NBT);
124 
125 #endif // ifndef NDEBUG
126 }
127 
128 BOOST_AUTO_TEST_CASE(test_minimal_crba)
129 {
132  pinocchio::Data data(model), data_ref(model);
133 
134  model.lowerPositionLimit.head<7>().fill(-1.);
135  model.upperPositionLimit.head<7>().fill(1.);
136 
137  Eigen::VectorXd q =
138  randomConfiguration(model, model.lowerPositionLimit, model.upperPositionLimit);
139  Eigen::VectorXd v(Eigen::VectorXd::Random(model.nv));
140 
142  data_ref.M.triangularView<Eigen::StrictlyLower>() =
143  data_ref.M.transpose().triangularView<Eigen::StrictlyLower>();
144 
146  data.M.triangularView<Eigen::StrictlyLower>() =
147  data.M.transpose().triangularView<Eigen::StrictlyLower>();
148 
149  BOOST_CHECK(data.M.isApprox(data_ref.M));
150 
151  ccrba(model, data_ref, q, v);
152  computeJointJacobians(model, data_ref, q);
153  BOOST_CHECK(data.Ag.isApprox(data_ref.Ag));
154  BOOST_CHECK(data.J.isApprox(data_ref.J));
155 
156  // Check double call
157  {
158  pinocchio::Data data2(model);
161 
162  BOOST_CHECK(data2.Ycrb[0] == data.Ycrb[0]);
163  }
164 }
165 
166 BOOST_AUTO_TEST_CASE(test_roto_inertia_effects)
167 {
168  pinocchio::Model model, model_ref;
170  model_ref = model;
171 
172  BOOST_CHECK(model == model_ref);
173 
174  pinocchio::Data data(model), data_ref(model_ref);
175 
176  model.armature = Eigen::VectorXd::Random(model.nv) + Eigen::VectorXd::Constant(model.nv, 1.);
177 
178  Eigen::VectorXd q = randomConfiguration(model);
179  pinocchio::crba(model_ref, data_ref, q, pinocchio::Convention::WORLD);
180  data_ref.M.triangularView<Eigen::StrictlyLower>() =
181  data_ref.M.transpose().triangularView<Eigen::StrictlyLower>();
182  data_ref.M.diagonal() += model.armature;
183 
185  data.M.triangularView<Eigen::StrictlyLower>() =
186  data.M.transpose().triangularView<Eigen::StrictlyLower>();
187 
188  BOOST_CHECK(data.M.isApprox(data_ref.M));
189 }
190 
191 #ifndef NDEBUG
192 
193 BOOST_AUTO_TEST_CASE(test_crba_malloc)
194 {
195  using namespace pinocchio;
198 
199  model.addJoint(
200  size_t(model.njoints - 1), pinocchio::JointModelRevoluteUnaligned(SE3::Vector3::UnitX()),
201  SE3::Random(), "revolute_unaligned");
203 
204  const Eigen::VectorXd q = pinocchio::neutral(model);
205  Eigen::internal::set_is_malloc_allowed(false);
206  crba(model, data, q);
207  Eigen::internal::set_is_malloc_allowed(true);
208 }
209 
210 #endif
211 
212 BOOST_AUTO_TEST_SUITE_END()
pinocchio::DataTpl::M
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
Definition: multibody/data.hpp:199
pinocchio::WORLD
@ WORLD
Definition: multibody/fwd.hpp:48
display-shapes-meshcat.placement
placement
Definition: display-shapes-meshcat.py:23
PinocchioTicToc::toc
double toc()
Definition: timer.hpp:88
Eigen
pinocchio::DataTpl
Definition: context/generic.hpp:25
pinocchio::JointModelBase
Definition: joint-model-base.hpp:75
pinocchio::SE3Tpl< context::Scalar, context::Options >
pinocchio::Convention::WORLD
@ WORLD
inverse-kinematics.i
int i
Definition: inverse-kinematics.py:20
pinocchio::name
std::string name(const LieGroupGenericTpl< LieGroupCollection > &lg)
Visit a LieGroupVariant to get the name of it.
model.hpp
setup.data
data
Definition: cmake/cython/setup.in.py:48
pinocchio::buildModels::humanoidRandom
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
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::computeJointJacobians
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & computeJointJacobians(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Computes the full model Jacobian, i.e. the stack of all motion subspace expressed in the world frame....
pinocchio::randomConfiguration
void randomConfiguration(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorIn1 > &lowerLimits, const Eigen::MatrixBase< ConfigVectorIn2 > &upperLimits, const Eigen::MatrixBase< ReturnType > &qout)
Generate a configuration vector uniformly sampled among provided limits.
Definition: joint-configuration.hpp:325
rnea.hpp
timer.hpp
anymal-simulation.model
model
Definition: anymal-simulation.py:12
pinocchio::ModelTpl::JointIndex
pinocchio::JointIndex JointIndex
Definition: multibody/model.hpp:68
center-of-mass.hpp
SMOOTH
#define SMOOTH(s)
Definition: timer.hpp:38
joint-configuration.hpp
pinocchio::placement
const MotionDense< Motion2 > const SE3Tpl< SE3Scalar, SE3Options > & placement
Definition: spatial/classic-acceleration.hpp:122
autodiff-rnea.v
v
Definition: autodiff-rnea.py:15
BOOST_AUTO_TEST_CASE
BOOST_AUTO_TEST_CASE(test_crba)
Definition: unittest/crba.cpp:64
data.hpp
pinocchio::ccrba
const DataTpl< Scalar, Options, JointCollectionTpl >::Matrix6x & ccrba(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v)
Computes the Centroidal Momentum Matrix, the Composite Ridig Body Inertia as well as the centroidal m...
pinocchio::InertiaTpl< context::Scalar, context::Options >::Random
static InertiaTpl Random()
Definition: spatial/inertia.hpp:356
M
M
pinocchio::DataTpl::Ag
Matrix6x Ag
Centroidal Momentum Matrix.
Definition: multibody/data.hpp:284
setup.name
name
Definition: cmake/cython/setup.in.py:179
q
q
a
Vec3f a
addJointAndBody
static void addJointAndBody(pinocchio::Model &model, const pinocchio::JointModelBase< JointModel > &joint, const std::string &parent_name, const std::string &name, const pinocchio::SE3 placement=pinocchio::SE3::Random(), bool setRandomLimits=true)
Definition: unittest/crba.cpp:34
PinocchioTicToc::US
@ US
Definition: timer.hpp:53
pinocchio::bias
MotionTpl< Scalar, Options > bias(const JointDataTpl< Scalar, Options, JointCollectionTpl > &jdata)
Visit a JointDataTpl through JointBiasVisitor to get the joint bias as a dense motion.
pinocchio::DataTpl::J
Matrix6x J
Jacobian of joint placements.
Definition: multibody/data.hpp:360
pinocchio::JointModelRevoluteUnalignedTpl
Definition: multibody/joint/fwd.hpp:38
fill
fill
pinocchio::SE3Tpl< context::Scalar, context::Options >::Random
static SE3Tpl Random()
Definition: spatial/se3-tpl.hpp:154
centroidal.hpp
pinocchio::rnea
const DataTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType & rnea(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 > &a)
The Recursive Newton-Euler algorithm. It computes the inverse dynamics, aka the joint torques accordi...
pinocchio::Convention::LOCAL
@ LOCAL
pinocchio::SE3Tpl< context::Scalar, context::Options >::Identity
static SE3Tpl Identity()
Definition: spatial/se3-tpl.hpp:136
pinocchio::neutral
void neutral(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ReturnType > &qout)
Return the neutral configuration element related to the model configuration space.
Definition: joint-configuration.hpp:375
sample-models.hpp
jacobian.hpp
pinocchio::ModelTpl
Definition: context/generic.hpp:20
crba.hpp
PinocchioTicToc
Definition: timer.hpp:47
pinocchio::normalize
void normalize(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &qout)
Normalize a configuration vector.
Definition: joint-configuration.hpp:914
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1116
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:27
BOOST_CHECK
#define BOOST_CHECK(check)
Definition: overview-urdf.cpp:34
PinocchioTicToc::tic
void tic()
Definition: timer.hpp:82


pinocchio
Author(s):
autogenerated on Sat Jun 22 2024 02:41:46