unittest/model.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2016-2022 CNRS INRIA
3 //
4 
7 
15 
18 
19 #include <boost/test/unit_test.hpp>
20 #include <boost/utility/binary.hpp>
21 
22 using namespace pinocchio;
23 
24 BOOST_AUTO_TEST_SUITE(BOOST_TEST_MODULE)
25 
26 BOOST_AUTO_TEST_CASE(test_model_subtree)
27 {
28  Model model;
30 
31  Model::JointIndex idx_larm1 = model.getJointId("larm1_joint");
32  BOOST_CHECK(idx_larm1 < (Model::JointIndex)model.njoints);
33  Model::IndexVector & subtree = model.subtrees[idx_larm1];
34  BOOST_CHECK(subtree.size() == 6);
35 
36  for (size_t joint_id = 0; joint_id < model.joints.size(); ++joint_id)
37  {
38  const Model::IndexVector & children = model.children[joint_id];
39  for (size_t i = 0; i < children.size(); ++i)
40  {
41  BOOST_CHECK(model.parents[children[i]] == joint_id);
42  }
43  }
44 
45  for (size_t i = 1; i < subtree.size(); ++i)
46  {
47  BOOST_CHECK(model.parents[subtree[i]] == subtree[i - 1]);
48  }
49 
50  // Check that i starts subtree[i]
51  for (JointIndex joint_id = 1; joint_id < (JointIndex)model.njoints; ++joint_id)
52  {
53  BOOST_CHECK(model.subtrees[joint_id][0] == joint_id);
54  }
55 
56  // Check that subtree[0] contains all joint ids
57  for (JointIndex joint_id = 1; joint_id < (JointIndex)model.njoints; ++joint_id)
58  {
59  BOOST_CHECK(model.subtrees[0][joint_id - 1] == joint_id);
60  }
61 }
62 
63 BOOST_AUTO_TEST_CASE(test_model_get_frame_id)
64 {
65  Model model;
67 
68  for (FrameIndex i = 0; i < static_cast<FrameIndex>(model.nframes); i++)
69  {
70  BOOST_CHECK_EQUAL(i, model.getFrameId(model.frames[i].name));
71  }
72  BOOST_CHECK_EQUAL(model.nframes, model.getFrameId("NOT_A_FRAME"));
73 }
74 
75 BOOST_AUTO_TEST_CASE(test_model_support)
76 {
77  Model model;
79  const Model::IndexVector support0_ref(1, 0);
80  BOOST_CHECK(model.supports[0] == support0_ref);
81 
82  // Check that i ends supports[i]
83  for (JointIndex joint_id = 1; joint_id < (JointIndex)model.njoints; ++joint_id)
84  {
85  BOOST_CHECK(model.supports[joint_id].back() == joint_id);
86  Model::IndexVector & support = model.supports[joint_id];
87 
88  size_t current_id = support.size() - 2;
89  for (JointIndex parent_id = model.parents[joint_id]; parent_id > 0;
90  parent_id = model.parents[parent_id], current_id--)
91  {
92  BOOST_CHECK(parent_id == support[current_id]);
93  }
94  }
95 }
96 
99 BOOST_AUTO_TEST_CASE(test_model_mimic_joint_support)
100 {
101  Model model;
102  // j1 -- j2 -- j3
103  // -- j4 -- j5
104  // with j5 mimicking j3
105  // with j2 mimicking j1
106  model.addJoint(0, JointModelRX(), SE3::Identity(), "j1");
107  model.addJoint(
108  1, JointModelMimic(JointModelRX(), model.joints[1], 1., 0.), SE3::Identity(), "j2");
109  model.addJoint(2, JointModelRX(), SE3::Identity(), "j3");
110  model.addJoint(1, JointModelRX(), SE3::Identity(), "j4");
111  model.addJoint(
112  4, JointModelMimic(JointModelRX(), model.joints[3], 1., 0.), SE3::Identity(), "j5");
113 
114  BOOST_REQUIRE_EQUAL(model.mimic_joint_supports.size(), 6);
115  BOOST_CHECK((model.mimic_joint_supports[0] == Model::IndexVector({0})));
116  BOOST_CHECK((model.mimic_joint_supports[1] == Model::IndexVector({0})));
117  // Mimic joint should support itself
118  BOOST_CHECK((model.mimic_joint_supports[2] == Model::IndexVector({0, 2})));
119  // j3 is supported by j2
120  BOOST_CHECK((model.mimic_joint_supports[3] == Model::IndexVector({0, 2})));
121  BOOST_CHECK((model.mimic_joint_supports[4] == Model::IndexVector({0})));
122  // Mimic joint should support itself
123  BOOST_CHECK((model.mimic_joint_supports[5] == Model::IndexVector({0, 5})));
124 }
125 
126 BOOST_AUTO_TEST_CASE(test_model_subspace_dimensions)
127 {
128  Model model;
130 
131  // Check that i ends supports[i]
132  for (JointIndex joint_id = 1; joint_id < (JointIndex)model.njoints; ++joint_id)
133  {
134  const Model::JointModel & jmodel = model.joints[joint_id];
135 
136  BOOST_CHECK(model.nqs[joint_id] == jmodel.nq());
137  BOOST_CHECK(model.idx_qs[joint_id] == jmodel.idx_q());
138 
139  BOOST_CHECK(model.nvs[joint_id] == jmodel.nv());
140  BOOST_CHECK(model.idx_vs[joint_id] == jmodel.idx_v());
141 
142  BOOST_CHECK(model.nvExtendeds[joint_id] == jmodel.nvExtended());
143  BOOST_CHECK(model.idx_vExtendeds[joint_id] == jmodel.idx_vExtended());
144  }
145 }
146 
148 {
149  Model model;
151 
152  BOOST_CHECK(model == model);
153 }
154 
156 {
157  Model model;
159 
160  BOOST_CHECK(model.cast<double>() == model.cast<double>());
161  BOOST_CHECK(model.cast<double>().cast<long double>() == model.cast<long double>());
162 
163  typedef ModelTpl<long double> Modelld;
164 
165  Modelld model2(model);
166  BOOST_CHECK(model2 == model.cast<long double>());
167 }
168 
169 BOOST_AUTO_TEST_CASE(test_std_vector_of_Model)
170 {
171  Model model;
173 
175  for (size_t k = 0; k < 20; ++k)
176  {
177  models.push_back(Model());
178  buildModels::humanoid(models.back());
179 
180  BOOST_CHECK(model == models.back());
181  }
182 }
183 
184 #ifdef PINOCCHIO_WITH_HPP_FCL
185 struct AddPrefix
186 {
187  std::string p;
188  std::string operator()(const std::string & n)
189  {
190  return p + n;
191  }
192  Frame operator()(const Frame & _f)
193  {
194  Frame f(_f);
195  f.name = p + f.name;
196  return f;
197  }
198  AddPrefix(const char * _p)
199  : p(_p)
200  {
201  }
202 };
203 
205 {
207  GeometryModel geomManipulator, geomHumanoid;
208 
210  buildModels::manipulatorGeometries(manipulator, geomManipulator);
211  geomManipulator.addAllCollisionPairs();
212  // Add prefix to joint and frame names
213  AddPrefix addManipulatorPrefix("manipulator/");
214  std::transform(
215  ++manipulator.names.begin(), manipulator.names.end(), ++manipulator.names.begin(),
216  addManipulatorPrefix);
217  std::transform(
218  ++manipulator.frames.begin(), manipulator.frames.end(), ++manipulator.frames.begin(),
219  addManipulatorPrefix);
220 
221  BOOST_TEST_MESSAGE(manipulator);
222 
224  buildModels::humanoidGeometries(humanoid, geomHumanoid);
225  geomHumanoid.addAllCollisionPairs();
226  // Add prefix to joint and frame names
227  AddPrefix addHumanoidPrefix("humanoid/");
228  std::transform(
229  ++humanoid.names.begin(), humanoid.names.end(), ++humanoid.names.begin(), addHumanoidPrefix);
230  std::transform(
231  ++humanoid.frames.begin(), humanoid.frames.end(), ++humanoid.frames.begin(), addHumanoidPrefix);
232 
233  BOOST_TEST_MESSAGE(humanoid);
234 
235  typename Model::ConfigVectorType humanoid_config_vector(humanoid.nq);
236  typename Model::ConfigVectorType manipulator_config_vector(manipulator.nq);
237  humanoid_config_vector = randomConfiguration(humanoid);
238  manipulator_config_vector = randomConfiguration(manipulator);
239  humanoid.referenceConfigurations.insert(std::make_pair("common_key", humanoid_config_vector));
240  manipulator.referenceConfigurations.insert(
241  std::make_pair("common_key", manipulator_config_vector));
242 
243  humanoid_config_vector = randomConfiguration(humanoid);
244  manipulator_config_vector = randomConfiguration(manipulator);
245  humanoid.referenceConfigurations.insert(std::make_pair("humanoid_key", humanoid_config_vector));
246  manipulator.referenceConfigurations.insert(
247  std::make_pair("manipulator_key", manipulator_config_vector));
248 
249  // TODO fix inertia of the base
250  manipulator.inertias[0].setRandom();
251  SE3 aMb = SE3::Random();
252 
253  // First append a model to the universe frame.
254  Model model1;
255  GeometryModel geomModel1;
256  FrameIndex fid = 0;
257  appendModel(
258  humanoid, manipulator, geomHumanoid, geomManipulator, fid, SE3::Identity(), model1, geomModel1);
259  typedef typename Model::ConfigVectorMap ConfigVectorMap;
260 
261  typename Model::ConfigVectorType neutral_config_vector(model1.nq);
262  neutral(model1, neutral_config_vector);
263 
264  BOOST_CHECK(model1.referenceConfigurations.size() == 3);
265  for (typename ConfigVectorMap::const_iterator config_it = model1.referenceConfigurations.begin();
266  config_it != model1.referenceConfigurations.end(); ++config_it)
267  {
268  const std::string & config_name = config_it->first;
269  const typename Model::ConfigVectorType & config_vector = config_it->second;
270 
271  typename ConfigVectorMap::const_iterator humanoid_config =
272  humanoid.referenceConfigurations.find(config_name);
273  typename ConfigVectorMap::const_iterator manipulator_config =
274  manipulator.referenceConfigurations.find(config_name);
275  for (JointIndex joint_id = 1; joint_id < model1.joints.size(); ++joint_id)
276  {
277  const JointModel & joint_model1 = model1.joints[joint_id];
278  if (
279  humanoid_config != humanoid.referenceConfigurations.end()
280  && humanoid.existJointName(model1.names[joint_id]))
281  { // key and joint exists in humanoid
282  const JointModel & joint_model_humanoid =
283  humanoid.joints[humanoid.getJointId(model1.names[joint_id])];
284  BOOST_CHECK(
285  joint_model_humanoid.jointConfigSelector(humanoid_config->second)
286  == joint_model1.jointConfigSelector(config_vector));
287  BOOST_CHECK(
288  joint_model_humanoid.JointMappedConfigSelector(humanoid_config->second)
289  == joint_model1.JointMappedConfigSelector(config_vector));
290  // std::cerr<<"humanoid "<<config_name<<" "<<model1.names[joint_id]<<std::endl;
291  }
292  else if (
293  manipulator_config != manipulator.referenceConfigurations.end()
294  && manipulator.existJointName(model1.names[joint_id]))
295  { // key and joint exists in manipulator.
296  const JointModel & joint_model_manipulator =
297  manipulator.joints[manipulator.getJointId(model1.names[joint_id])];
298  BOOST_CHECK(
299  joint_model_manipulator.jointConfigSelector(manipulator_config->second)
300  == joint_model1.jointConfigSelector(config_vector));
301  BOOST_CHECK(
302  joint_model_manipulator.JointMappedConfigSelector(manipulator_config->second)
303  == joint_model1.JointMappedConfigSelector(config_vector));
304  // std::cerr<<"manipulator "<<config_name<<" "<<model1.names[joint_id]<<std::endl;
305  }
306  else
307  { // joint and key combo not found, should with neutral
308  BOOST_CHECK(
309  joint_model1.jointConfigSelector(neutral_config_vector)
310  == joint_model1.jointConfigSelector(config_vector));
311  BOOST_CHECK(
312  joint_model1.JointMappedConfigSelector(neutral_config_vector)
313  == joint_model1.JointMappedConfigSelector(config_vector));
314  // std::cerr<<"neutral "<<config_name<<" "<<model1.names[joint_id]<<std::endl;
315  }
316  }
317  }
318 
319  {
320  Model model2 = appendModel(humanoid, manipulator, fid, SE3::Identity());
321  Model model3;
322  appendModel(humanoid, manipulator, fid, SE3::Identity(), model3);
323  BOOST_CHECK(model1 == model2);
324  BOOST_CHECK(model1 == model3);
325  BOOST_CHECK(model2 == model3);
326  }
327 
328  Data data1(model1);
329  BOOST_CHECK(model1.check(data1));
330 
331  BOOST_TEST_MESSAGE(model1);
332 
333  // Second, append a model to a moving frame.
334  Model model;
335 
336  GeometryModel geomModel;
337  fid = humanoid.addFrame(Frame(
338  "humanoid/add_manipulator", humanoid.getJointId("humanoid/chest2_joint"),
339  humanoid.getFrameId("humanoid/chest2_joint"), aMb, OP_FRAME));
340 
341  // Append manipulator to chest2_joint of humanoid
342  appendModel(
343  humanoid, manipulator, geomHumanoid, geomManipulator, fid, SE3::Identity(), model, geomModel);
344 
345  neutral_config_vector.resize(model.nq);
346  neutral(model, neutral_config_vector);
347 
348  BOOST_CHECK(model.referenceConfigurations.size() == 3);
349  for (typename ConfigVectorMap::const_iterator config_it = model.referenceConfigurations.begin();
350  config_it != model.referenceConfigurations.end(); ++config_it)
351  {
352  const std::string & config_name = config_it->first;
353  const typename Model::ConfigVectorType & config_vector = config_it->second;
354 
355  typename ConfigVectorMap::const_iterator humanoid_config =
356  humanoid.referenceConfigurations.find(config_name);
357  typename ConfigVectorMap::const_iterator manipulator_config =
358  manipulator.referenceConfigurations.find(config_name);
359  for (JointIndex joint_id = 1; joint_id < model.joints.size(); ++joint_id)
360  {
361  const JointModel & joint_model = model.joints[joint_id];
362  if (
363  humanoid_config != humanoid.referenceConfigurations.end()
364  && humanoid.existJointName(model.names[joint_id]))
365  { // key and joint exists in humanoid
366  const JointModel & joint_model_humanoid =
367  humanoid.joints[humanoid.getJointId(model.names[joint_id])];
368  BOOST_CHECK(
369  joint_model_humanoid.jointConfigSelector(humanoid_config->second)
370  == joint_model.jointConfigSelector(config_vector));
371  BOOST_CHECK(
372  joint_model_humanoid.JointMappedConfigSelector(humanoid_config->second)
373  == joint_model.JointMappedConfigSelector(config_vector));
374  // std::cerr<<"humanoid "<<config_name<<" "<<model.names[joint_id]<<std::endl;
375  }
376  else if (
377  manipulator_config != manipulator.referenceConfigurations.end()
378  && manipulator.existJointName(model.names[joint_id]))
379  { // key and joint exists in manipulator.
380  const JointModel & joint_model_manipulator =
381  manipulator.joints[manipulator.getJointId(model.names[joint_id])];
382  BOOST_CHECK(
383  joint_model_manipulator.jointConfigSelector(manipulator_config->second)
384  == joint_model.jointConfigSelector(config_vector));
385  BOOST_CHECK(
386  joint_model_manipulator.JointMappedConfigSelector(manipulator_config->second)
387  == joint_model.JointMappedConfigSelector(config_vector));
388  // std::cerr<<"manipulator "<<config_name<<" "<<model.names[joint_id]<<std::endl;
389  }
390  else
391  { // joint and key combo not found, should with neutral
392  BOOST_CHECK(
393  joint_model.jointConfigSelector(neutral_config_vector)
394  == joint_model.jointConfigSelector(config_vector));
395  BOOST_CHECK(
396  joint_model.JointMappedConfigSelector(neutral_config_vector)
397  == joint_model.JointMappedConfigSelector(config_vector));
398  // std::cerr<<"neutral "<<config_name<<" "<<model.names[joint_id]<<std::endl;
399  }
400  }
401  }
402 
403  {
404  Model model2 = appendModel(humanoid, manipulator, fid, SE3::Identity());
405  Model model3;
406  appendModel(humanoid, manipulator, fid, SE3::Identity(), model3);
407  BOOST_CHECK(model == model2);
408  BOOST_CHECK(model == model3);
409  BOOST_CHECK(model2 == model3);
410  }
411 
412  BOOST_TEST_MESSAGE(model);
413 
414  Data data(model);
415  BOOST_CHECK(model.check(data));
416 
417  // Check the model
418  BOOST_CHECK_EQUAL(
419  model.getJointId("humanoid/chest2_joint"),
420  model.parents[model.getJointId("manipulator/shoulder1_joint")]);
421 
422  // check the joint order and the inertias
423  // All the joints of the manipulator should be at the end of the merged model
424  JointIndex hnj = (JointIndex)humanoid.njoints;
425  JointIndex chest2 = model.getJointId("humanoid/chest2_joint");
426  for (JointIndex jid = 1; jid < hnj; ++jid)
427  {
428  BOOST_TEST_MESSAGE("Checking joint " << jid << " " << model.names[jid]);
429  BOOST_CHECK_EQUAL(model.names[jid], humanoid.names[jid]);
430  if (jid != chest2)
431  BOOST_CHECK_EQUAL(model.inertias[jid], humanoid.inertias[jid]);
432  else
433  BOOST_CHECK_MESSAGE(
434  model.inertias[jid].isApprox(
435  manipulator.inertias[0].se3Action(aMb) + humanoid.inertias[jid]),
436  model.inertias[jid] << " != "
437  << manipulator.inertias[0].se3Action(aMb) + humanoid.inertias[jid]);
438 
439  BOOST_CHECK_EQUAL(model.jointPlacements[jid], humanoid.jointPlacements[jid]);
440  }
441  for (JointIndex jid = 1; jid < manipulator.joints.size() - 1; ++jid)
442  {
443  BOOST_TEST_MESSAGE("Checking joint " << hnj - 1 + jid << " " << model.names[hnj + jid]);
444  BOOST_CHECK_EQUAL(model.names[hnj - 1 + jid], manipulator.names[jid]);
445  BOOST_CHECK_EQUAL(model.inertias[hnj - 1 + jid], manipulator.inertias[jid]);
446  if (jid == 1)
447  BOOST_CHECK_EQUAL(
448  model.jointPlacements[hnj - 1 + jid], aMb * manipulator.jointPlacements[jid]);
449  else
450  BOOST_CHECK_EQUAL(model.jointPlacements[hnj - 1 + jid], manipulator.jointPlacements[jid]);
451  }
452  // Check the frames
453  for (FrameIndex fid = 1; fid < humanoid.frames.size(); ++fid)
454  {
455  const Frame &frame = humanoid.frames[fid], parent = humanoid.frames[frame.parentFrame];
456  BOOST_CHECK(model.existFrame(frame.name, frame.type));
457  const Frame &nframe = model.frames[model.getFrameId(frame.name, frame.type)],
458  nparent = model.frames[nframe.parentFrame];
459  BOOST_CHECK_EQUAL(parent.name, nparent.name);
460  BOOST_CHECK_EQUAL(frame.placement, nframe.placement);
461  }
462  for (FrameIndex fid = 1; fid < manipulator.frames.size(); ++fid)
463  {
464  const Frame &frame = manipulator.frames[fid], parent = manipulator.frames[frame.parentFrame];
465  BOOST_CHECK(model.existFrame(frame.name, frame.type));
466  const Frame &nframe = model.frames[model.getFrameId(frame.name, frame.type)],
467  nparent = model.frames[nframe.parentFrame];
468  if (frame.parentFrame > 0)
469  {
470  BOOST_CHECK_EQUAL(parent.name, nparent.name);
471  BOOST_CHECK_EQUAL(frame.placement, nframe.placement);
472  }
473  }
474 
475  {
476  Inertia inertia(2., Eigen::Vector3d(0.1, 0.1, 0.1), Eigen::Matrix3d::Identity());
477  Frame additional_frame("inertial_frame", 2, SE3::Identity(), FrameType::JOINT, inertia);
478  humanoid.addFrame(additional_frame);
479  double mass_humanoid = computeTotalMass(humanoid);
480  double mass_manipulator = computeTotalMass(manipulator);
481  double total_mass = mass_manipulator + mass_humanoid;
482 
483  Model model4;
484  GeometryModel geomModel4;
485  appendModel(
486  humanoid, manipulator, geomHumanoid, geomManipulator, 0, SE3::Identity(), model4, geomModel4);
487  BOOST_CHECK_CLOSE(computeTotalMass(model4), total_mass, 1e-6);
488  }
489  {
490  Model ff_model;
491  auto ff_id = ff_model.addJoint(0, JointModelFreeFlyer(), SE3::Identity(), "floating_base");
492  ff_model.addJointFrame(ff_id);
493  GeometryModel ff_geom_model = GeometryModel();
494  FrameIndex frame_id = ff_model.getFrameId("floating_base");
495  Model model4;
496  GeometryModel geomModel4;
497  appendModel(
498  ff_model, manipulator, ff_geom_model, geomManipulator, frame_id, SE3::Identity(), model4,
499  geomModel4);
500  BOOST_CHECK(model4.inertias[1] == model4.inertias[1]);
501  }
502 
503  {
504  Model model5, gripperModel;
505 
506  Inertia inertia(1., SE3::Vector3(0.5, 0., 0.0), SE3::Matrix3::Identity());
507  SE3 pos(1);
508 
509  pos.translation() = SE3::LinearType(0.1, 0., 0.);
510  JointIndex idx = gripperModel.addJoint(0, JointModelPX(), pos, "left_finger");
511  gripperModel.addJointFrame(idx);
512 
513  pos.translation() = SE3::LinearType(-0.1, 0., 0.);
514  idx = gripperModel.addJoint(0, JointModelPX(), pos, "right_finger");
515  gripperModel.addJointFrame(idx);
516 
517  SE3 transformManGr = SE3::Random();
518  FrameIndex fid = (FrameIndex)(manipulator.frames.size() - 1);
519  appendModel(manipulator, gripperModel, fid, transformManGr, model5);
520 
521  JointIndex jid5 = model5.getJointId("left_finger");
522  JointIndex jidG = gripperModel.getJointId("left_finger");
523  BOOST_CHECK(
524  model5.jointPlacements[jid5].isApprox(transformManGr * gripperModel.jointPlacements[jidG]));
525 
526  jid5 = model5.getJointId("right_finger");
527  jidG = gripperModel.getJointId("right_finger");
528  BOOST_CHECK(
529  model5.jointPlacements[jid5].isApprox(transformManGr * gripperModel.jointPlacements[jidG]));
530  }
531 }
532 #endif
533 
534 BOOST_AUTO_TEST_CASE(test_buildReducedModel_empty)
535 {
536  Model humanoid_model;
537  buildModels::humanoid(humanoid_model);
538 
539  static const std::vector<JointIndex> empty_index_vector;
540 
541  humanoid_model.lowerPositionLimit.head<3>().fill(-1.);
542  humanoid_model.upperPositionLimit.head<3>().fill(1.);
543 
544  humanoid_model.referenceConfigurations.insert(
545  std::pair<std::string, Eigen::VectorXd>("neutral", neutral(humanoid_model)));
546  Eigen::VectorXd reference_config_humanoid = randomConfiguration(humanoid_model);
547  Model humanoid_copy_model =
548  buildReducedModel(humanoid_model, empty_index_vector, reference_config_humanoid);
549 
550  BOOST_CHECK(humanoid_copy_model.names == humanoid_model.names);
551  BOOST_CHECK(humanoid_copy_model.joints == humanoid_model.joints);
552  BOOST_CHECK(humanoid_copy_model == humanoid_model);
553  BOOST_CHECK(
554  humanoid_copy_model.referenceConfigurations["neutral"].isApprox(neutral(humanoid_copy_model)));
555 
556  const std::vector<JointIndex> empty_joints_to_lock;
557 
558  const Model reduced_humanoid_model =
559  buildReducedModel(humanoid_model, empty_joints_to_lock, reference_config_humanoid);
560  BOOST_CHECK(reduced_humanoid_model.njoints == humanoid_model.njoints);
561  BOOST_CHECK(reduced_humanoid_model.frames == humanoid_model.frames);
562  BOOST_CHECK(reduced_humanoid_model.jointPlacements == humanoid_model.jointPlacements);
563  BOOST_CHECK(reduced_humanoid_model.joints == humanoid_model.joints);
564 
565  for (JointIndex joint_id = 1; joint_id < (JointIndex)reduced_humanoid_model.njoints; ++joint_id)
566  {
567  BOOST_CHECK(
568  reduced_humanoid_model.inertias[joint_id].isApprox(humanoid_model.inertias[joint_id]));
569  }
570 }
571 
572 BOOST_AUTO_TEST_CASE(test_buildReducedModel)
573 {
574  Model humanoid_model;
575  buildModels::humanoid(humanoid_model);
576 
577  static const std::vector<JointIndex> empty_index_vector;
578 
579  humanoid_model.lowerPositionLimit.head<3>().fill(-1.);
580  humanoid_model.upperPositionLimit.head<3>().fill(1.);
581 
582  humanoid_model.referenceConfigurations.insert(
583  std::pair<std::string, Eigen::VectorXd>("neutral", neutral(humanoid_model)));
584  Eigen::VectorXd reference_config_humanoid = randomConfiguration(humanoid_model);
585  Model humanoid_copy_model =
586  buildReducedModel(humanoid_model, empty_index_vector, reference_config_humanoid);
587 
588  BOOST_CHECK(humanoid_copy_model.names == humanoid_model.names);
589  BOOST_CHECK(humanoid_copy_model.joints == humanoid_model.joints);
590  BOOST_CHECK(humanoid_copy_model == humanoid_model);
591  BOOST_CHECK(
592  humanoid_copy_model.referenceConfigurations["neutral"].isApprox(neutral(humanoid_copy_model)));
593 
594  std::vector<JointIndex> joints_to_lock;
595  const std::string joint1_to_lock = "rarm_shoulder2_joint";
596  joints_to_lock.push_back(humanoid_model.getJointId(joint1_to_lock));
597  const std::string joint2_to_lock = "larm_shoulder2_joint";
598  joints_to_lock.push_back(humanoid_model.getJointId(joint2_to_lock));
599 
600  Model reduced_humanoid_model =
601  buildReducedModel(humanoid_model, joints_to_lock, reference_config_humanoid);
602  BOOST_CHECK(
603  reduced_humanoid_model.njoints == humanoid_model.njoints - (int)joints_to_lock.size());
604 
605  BOOST_CHECK(reduced_humanoid_model != humanoid_model);
606 
607  Model reduced_humanoid_model_other_signature;
609  humanoid_model, joints_to_lock, reference_config_humanoid,
610  reduced_humanoid_model_other_signature);
611  BOOST_CHECK(reduced_humanoid_model == reduced_humanoid_model_other_signature);
612 
613  // Check that forward kinematics give same results
614  Data data(humanoid_model), reduced_data(reduced_humanoid_model);
615  Eigen::VectorXd q = reference_config_humanoid;
616  Eigen::VectorXd reduced_q(reduced_humanoid_model.nq);
617 
618  for (JointIndex joint_id = 1; joint_id < (JointIndex)reduced_humanoid_model.njoints; ++joint_id)
619  {
620  const JointIndex reference_joint_id =
621  humanoid_model.getJointId(reduced_humanoid_model.names[joint_id]);
622 
623  reduced_humanoid_model.joints[joint_id].jointConfigSelector(reduced_q) =
624  humanoid_model.joints[reference_joint_id].jointConfigSelector(q);
625  reduced_humanoid_model.joints[joint_id].JointMappedConfigSelector(reduced_q) =
626  humanoid_model.joints[reference_joint_id].JointMappedConfigSelector(q);
627  }
628 
629  BOOST_CHECK(reduced_humanoid_model.referenceConfigurations["neutral"].isApprox(
630  neutral(reduced_humanoid_model)));
631 
632  framesForwardKinematics(humanoid_model, data, q);
633  framesForwardKinematics(reduced_humanoid_model, reduced_data, reduced_q);
634 
635  for (size_t frame_id = 0; frame_id < reduced_humanoid_model.frames.size(); ++frame_id)
636  {
637  const Frame & reduced_frame = reduced_humanoid_model.frames[frame_id];
638  switch (reduced_frame.type)
639  {
640  case JOINT:
641  case FIXED_JOINT: {
642  // May not be present in the original model
643  if (humanoid_model.existJointName(reduced_frame.name))
644  {
645  const JointIndex joint_id = humanoid_model.getJointId(reduced_frame.name);
646  BOOST_CHECK(data.oMi[joint_id].isApprox(reduced_data.oMf[frame_id]));
647  }
648  else if (humanoid_model.existFrame(reduced_frame.name))
649  {
650  const FrameIndex humanoid_frame_id = humanoid_model.getFrameId(reduced_frame.name);
651  BOOST_CHECK(data.oMf[humanoid_frame_id].isApprox(reduced_data.oMf[frame_id]));
652  }
653  else
654  {
655  BOOST_CHECK_MESSAGE(
656  false, "The frame " << reduced_frame.name << " is not presend in the humanoid_model");
657  }
658  break;
659  }
660  default: {
661  BOOST_CHECK(humanoid_model.existFrame(reduced_frame.name));
662  const FrameIndex humanoid_frame_id = humanoid_model.getFrameId(reduced_frame.name);
663  BOOST_CHECK(data.oMf[humanoid_frame_id].isApprox(reduced_data.oMf[frame_id]));
664  break;
665  }
666  }
667  }
668 }
669 
670 BOOST_AUTO_TEST_CASE(test_aligned_vector_of_model)
671 {
672  typedef PINOCCHIO_ALIGNED_STD_VECTOR(Model) VectorOfModels;
673 
674  VectorOfModels models;
675  for (size_t k = 0; k < 100; ++k)
676  {
677  models.push_back(Model());
678  buildModels::humanoidRandom(models[k]);
679  }
680 }
681 
682 #ifdef PINOCCHIO_WITH_HPP_FCL
683 BOOST_AUTO_TEST_CASE(test_buildReducedModel_with_geom)
684 {
685  Model humanoid_model;
686  buildModels::humanoid(humanoid_model);
687 
688  humanoid_model.lowerPositionLimit.head<3>().fill(-1.);
689  humanoid_model.upperPositionLimit.head<3>().fill(1.);
690  const Eigen::VectorXd reference_config_humanoid = randomConfiguration(humanoid_model);
691 
692  GeometryModel humanoid_geometry;
693  buildModels::humanoidGeometries(humanoid_model, humanoid_geometry);
694 
695  static const std::vector<JointIndex> empty_index_vector;
696 
697  Model humanoid_copy_model;
698  GeometryModel humanoid_copy_geometry;
700  humanoid_model, humanoid_geometry, empty_index_vector, reference_config_humanoid,
701  humanoid_copy_model, humanoid_copy_geometry);
702 
703  BOOST_CHECK(humanoid_copy_model == humanoid_model);
704  BOOST_CHECK(humanoid_copy_geometry == humanoid_geometry);
705 
706  std::vector<JointIndex> joints_to_lock;
707  const std::string joint1_to_lock = "rarm_shoulder2_joint";
708  joints_to_lock.push_back(humanoid_model.getJointId(joint1_to_lock));
709  const std::string joint2_to_lock = "larm_shoulder2_joint";
710  joints_to_lock.push_back(humanoid_model.getJointId(joint2_to_lock));
711 
712  Model reduced_humanoid_model;
713  GeometryModel reduced_humanoid_geometry;
715  humanoid_model, humanoid_geometry, joints_to_lock, reference_config_humanoid,
716  reduced_humanoid_model, reduced_humanoid_geometry);
717 
718  BOOST_CHECK(reduced_humanoid_geometry.ngeoms == humanoid_geometry.ngeoms);
719  BOOST_CHECK(reduced_humanoid_geometry.collisionPairs == humanoid_geometry.collisionPairs);
720  BOOST_CHECK(
721  reduced_humanoid_geometry.geometryObjects.size() == humanoid_geometry.geometryObjects.size());
722 
723  for (Index i = 0; i < humanoid_geometry.geometryObjects.size(); ++i)
724  {
725  const GeometryObject & go1 = humanoid_geometry.geometryObjects[i];
726  const GeometryObject & go2 = reduced_humanoid_geometry.geometryObjects[i];
727  BOOST_CHECK_EQUAL(go1.name, go2.name);
728  BOOST_CHECK_EQUAL(go1.geometry, go2.geometry);
729  BOOST_CHECK_EQUAL(go1.meshPath, go2.meshPath);
730  BOOST_CHECK_EQUAL(go1.meshScale, go2.meshScale);
731  BOOST_CHECK_EQUAL(go1.overrideMaterial, go2.overrideMaterial);
732  BOOST_CHECK_EQUAL(go1.meshColor, go2.meshColor);
733  BOOST_CHECK_EQUAL(go1.meshTexturePath, go2.meshTexturePath);
734  BOOST_CHECK_EQUAL(go1.parentFrame, go2.parentFrame);
735  BOOST_CHECK_EQUAL(
736  humanoid_model.frames[go1.parentFrame].name,
737  reduced_humanoid_model.frames[go2.parentFrame].name);
738  }
739 
740  Data data(humanoid_model), reduced_data(reduced_humanoid_model);
741  const Eigen::VectorXd q = reference_config_humanoid;
742  Eigen::VectorXd reduced_q(reduced_humanoid_model.nq);
743 
744  for (JointIndex joint_id = 1; joint_id < (JointIndex)reduced_humanoid_model.njoints; ++joint_id)
745  {
746  const JointIndex reference_joint_id =
747  humanoid_model.getJointId(reduced_humanoid_model.names[joint_id]);
748 
749  reduced_humanoid_model.joints[joint_id].jointConfigSelector(reduced_q) =
750  humanoid_model.joints[reference_joint_id].jointConfigSelector(q);
751  reduced_humanoid_model.joints[joint_id].JointMappedConfigSelector(reduced_q) =
752  humanoid_model.joints[reference_joint_id].JointMappedConfigSelector(q);
753  }
754 
755  framesForwardKinematics(humanoid_model, data, q);
756  framesForwardKinematics(reduced_humanoid_model, reduced_data, reduced_q);
757 
758  for (size_t frame_id = 0; frame_id < reduced_humanoid_model.frames.size(); ++frame_id)
759  {
760  const Frame & reduced_frame = reduced_humanoid_model.frames[frame_id];
761  switch (reduced_frame.type)
762  {
763  case JOINT:
764  case FIXED_JOINT: {
765  // May not be present in the original model
766  if (humanoid_model.existJointName(reduced_frame.name))
767  {
768  const JointIndex joint_id = humanoid_model.getJointId(reduced_frame.name);
769  BOOST_CHECK(data.oMi[joint_id].isApprox(reduced_data.oMf[frame_id]));
770  }
771  else if (humanoid_model.existFrame(reduced_frame.name))
772  {
773  const FrameIndex humanoid_frame_id = humanoid_model.getFrameId(reduced_frame.name);
774  BOOST_CHECK(data.oMf[humanoid_frame_id].isApprox(reduced_data.oMf[frame_id]));
775  }
776  else
777  {
778  BOOST_CHECK_MESSAGE(
779  false, "The frame " << reduced_frame.name << " is not presend in the humanoid_model");
780  }
781  break;
782  }
783  default: {
784  BOOST_CHECK(humanoid_model.existFrame(reduced_frame.name));
785  const FrameIndex humanoid_frame_id = humanoid_model.getFrameId(reduced_frame.name);
786  BOOST_CHECK(data.oMf[humanoid_frame_id].isApprox(reduced_data.oMf[frame_id]));
787  break;
788  }
789  }
790  }
791 
792  // Test GeometryObject placements
793  GeometryData geom_data(humanoid_geometry), reduded_geom_data(reduced_humanoid_geometry);
794  updateGeometryPlacements(humanoid_model, data, humanoid_geometry, geom_data);
796  reduced_humanoid_model, reduced_data, reduced_humanoid_geometry, reduded_geom_data);
797 
798  BOOST_CHECK(geom_data.oMg.size() == reduded_geom_data.oMg.size());
799  for (FrameIndex i = 0; i < geom_data.oMg.size(); ++i)
800  {
801  BOOST_CHECK(geom_data.oMg[i].isApprox(reduded_geom_data.oMg[i]));
802  }
803 
804  // Test other signature
805  std::vector<GeometryModel> full_geometry_models;
806  full_geometry_models.push_back(humanoid_geometry);
807  full_geometry_models.push_back(humanoid_geometry);
808  full_geometry_models.push_back(humanoid_geometry);
809 
810  std::vector<GeometryModel> reduced_geometry_models;
811 
812  Model reduced_humanoid_model_other_sig;
814  humanoid_model, full_geometry_models, joints_to_lock, reference_config_humanoid,
815  reduced_humanoid_model_other_sig, reduced_geometry_models);
816 
817  BOOST_CHECK(reduced_geometry_models[0] == reduced_humanoid_geometry);
818  BOOST_CHECK(reduced_geometry_models[1] == reduced_humanoid_geometry);
819  BOOST_CHECK(reduced_geometry_models[2] == reduced_humanoid_geometry);
820 }
821 #endif // PINOCCHIO_WITH_HPP_FCL
822 
823 BOOST_AUTO_TEST_CASE(test_findCommonAncestor)
824 {
825  Model model;
827 
828  {
829  size_t id_ancestor1, id_ancestor2;
830  JointIndex ancestor = findCommonAncestor(model, 0, 0, id_ancestor1, id_ancestor2);
831  BOOST_CHECK(ancestor == 0);
832  BOOST_CHECK(id_ancestor1 == 0);
833  BOOST_CHECK(id_ancestor2 == 0);
834  }
835 
836  {
837  size_t id_ancestor1, id_ancestor2;
838  JointIndex ancestor =
839  findCommonAncestor(model, 0, (JointIndex)(model.njoints - 1), id_ancestor1, id_ancestor2);
840  BOOST_CHECK(ancestor == 0);
841  BOOST_CHECK(id_ancestor1 == 0);
842  BOOST_CHECK(id_ancestor2 == 0);
843  }
844 
845  {
846  size_t id_ancestor1, id_ancestor2;
847  JointIndex ancestor =
848  findCommonAncestor(model, (JointIndex)(model.njoints - 1), 0, id_ancestor1, id_ancestor2);
849  BOOST_CHECK(ancestor == 0);
850  BOOST_CHECK(id_ancestor1 == 0);
851  BOOST_CHECK(id_ancestor2 == 0);
852  }
853 
854  {
855  size_t id_ancestor1, id_ancestor2;
856  JointIndex ancestor =
857  findCommonAncestor(model, (JointIndex)(model.njoints - 1), 1, id_ancestor1, id_ancestor2);
858  BOOST_CHECK(ancestor == 1);
859  BOOST_CHECK(model.supports[(JointIndex)(model.njoints - 1)][id_ancestor1] == ancestor);
860  BOOST_CHECK(model.supports[1][id_ancestor2] == ancestor);
861  }
862 }
863 
864 BOOST_AUTO_TEST_CASE(test_has_configuration_limit)
865 {
866  using namespace Eigen;
867 
868  // Test joint specific function hasConfigurationLimit
869  JointModelFreeFlyer test_joint_ff = JointModelFreeFlyer();
870  std::vector<bool> cf_limits_ff = test_joint_ff.hasConfigurationLimit();
871  std::vector<bool> cf_limits_tangent_ff = test_joint_ff.hasConfigurationLimitInTangent();
872  std::vector<bool> expected_cf_limits_ff({true, true, true, false, false, false, false});
873  std::vector<bool> expected_cf_limits_tangent_ff({true, true, true, false, false, false});
874  BOOST_CHECK(cf_limits_ff == expected_cf_limits_ff);
875  BOOST_CHECK(cf_limits_tangent_ff == expected_cf_limits_tangent_ff);
876 
877  JointModelPlanar test_joint_planar = JointModelPlanar();
878  std::vector<bool> cf_limits_planar = test_joint_planar.hasConfigurationLimit();
879  std::vector<bool> cf_limits_tangent_planar = test_joint_planar.hasConfigurationLimitInTangent();
880  std::vector<bool> expected_cf_limits_planar({true, true, false, false});
881  std::vector<bool> expected_cf_limits_tangent_planar({true, true, false});
882  BOOST_CHECK(cf_limits_planar == expected_cf_limits_planar);
883  BOOST_CHECK(cf_limits_tangent_planar == expected_cf_limits_tangent_planar);
884 
885  JointModelPX test_joint_p = JointModelPX();
886  std::vector<bool> cf_limits_prismatic = test_joint_p.hasConfigurationLimit();
887  std::vector<bool> cf_limits_tangent_prismatic = test_joint_p.hasConfigurationLimitInTangent();
888  std::vector<bool> expected_cf_limits_prismatic({true});
889  std::vector<bool> expected_cf_limits_tangent_prismatic({true});
890  BOOST_CHECK(cf_limits_prismatic == expected_cf_limits_prismatic);
891  BOOST_CHECK(cf_limits_tangent_prismatic == expected_cf_limits_tangent_prismatic);
892 
893  // Test model.hasConfigurationLimit() function
894  Model model;
895  JointIndex jointId = 0;
896 
897  jointId = model.addJoint(jointId, JointModelFreeFlyer(), SE3::Identity(), "Joint0");
898  jointId = model.addJoint(jointId, JointModelRZ(), SE3::Identity(), "Joint1");
899  jointId = model.addJoint(
900  jointId, JointModelRUBZ(), SE3(Matrix3d::Identity(), Vector3d(1.0, 0.0, 0.0)), "Joint2");
901 
902  std::vector<bool> expected_cf_limits_model(
903  {true, true, true, // translation of FF
904  false, false, false, false, // rotation of FF
905  true, // roational joint
906  false, false}); // unbounded rotational joint
907  std::vector<bool> model_cf_limits = model.hasConfigurationLimit();
908  BOOST_CHECK((model_cf_limits == expected_cf_limits_model));
909 
910  // Test model.hasConfigurationLimitInTangent() function
911  std::vector<bool> expected_cf_limits_tangent_model(
912  {true, true, true, // translation of FF
913  false, false, false, // rotation of FF
914  true, // roational joint
915  false}); // unbounded rotational joint
916  std::vector<bool> model_cf_limits_tangent = model.hasConfigurationLimitInTangent();
917  BOOST_CHECK((model_cf_limits_tangent == expected_cf_limits_tangent_model));
918 }
919 
922 BOOST_AUTO_TEST_CASE(test_transform_to_mimic)
923 {
924  Model model_exp;
925  // j1 -- j2 -- j3
926  // -- j4 -- j5
927  // with j5 mimicking j3
928  // with j2 mimicking j1
929  model_exp.addJoint(0, JointModelRX(), SE3::Identity(), "j1");
930  model_exp.addJoint(
931  1, JointModelMimic(JointModelRX(), model_exp.joints[1], 1., 0.), SE3::Identity(), "j2");
932  model_exp.addJoint(2, JointModelRX(), SE3::Identity(), "j3");
933  model_exp.addJoint(1, JointModelRX(), SE3::Identity(), "j4");
934  model_exp.addJoint(
935  4, JointModelMimic(JointModelRX(), model_exp.joints[3], 1., 0.), SE3::Identity(), "j5");
936 
937  // Model transformed with transformJointIntoMimic
938  Model model;
940  Model model_mimic2;
941  model.addJoint(0, JointModelRX(), SE3::Identity(), "j1");
942  model.addJoint(1, JointModelRX(), SE3::Identity(), "j2");
943  model.addJoint(2, JointModelRX(), SE3::Identity(), "j3");
944  model.addJoint(1, JointModelRX(), SE3::Identity(), "j4");
945  model.addJoint(4, JointModelRX(), SE3::Identity(), "j5");
946  transformJointIntoMimic(model, 3, 5, 1., 0., model_mimic);
947  transformJointIntoMimic(model_mimic, 1, 2, 1., 0., model_mimic2);
948 
949  BOOST_CHECK(model_mimic2 == model_exp);
950 }
951 
954 BOOST_AUTO_TEST_CASE(test_build_mimic_model)
955 {
956  Model model_exp;
957  // j1 -- j2 -- j3
958  // -- j4 -- j5
959  // with j5 mimicking j3
960  // with j2 mimicking j1
961  model_exp.addJoint(0, JointModelRX(), SE3::Identity(), "j1");
962  model_exp.addJoint(
963  1, JointModelMimic(JointModelRX(), model_exp.joints[1], 1., 0.), SE3::Identity(), "j2");
964  model_exp.addJoint(2, JointModelRX(), SE3::Identity(), "j3");
965  model_exp.addJoint(1, JointModelRX(), SE3::Identity(), "j4");
966  model_exp.addJoint(
967  4, JointModelMimic(JointModelRX(), model_exp.joints[3], 2., 1.), SE3::Identity(), "j5");
968 
969  // Model transformed with buildMimicModel
970  Model model;
972  model.addJoint(0, JointModelRX(), SE3::Identity(), "j1");
973  model.addJoint(1, JointModelRX(), SE3::Identity(), "j2");
974  model.addJoint(2, JointModelRX(), SE3::Identity(), "j3");
975  model.addJoint(1, JointModelRX(), SE3::Identity(), "j4");
976  model.addJoint(4, JointModelRX(), SE3::Identity(), "j5");
977  buildMimicModel(model, {1, 3}, {2, 5}, {1., 2.}, {0., 1.}, model_mimic);
978 
979  BOOST_CHECK(model_mimic == model_exp);
980 }
981 
982 BOOST_AUTO_TEST_CASE(test_cast_mimic)
983 {
984  Model humanoid_model, humanoid_mimic;
985  buildModels::humanoid(humanoid_model);
986  Data data(humanoid_model);
987  BOOST_CHECK(humanoid_model.check(data));
988 
989  JointIndex index_p = humanoid_model.getJointId("rleg_shoulder3_joint");
990  JointIndex index_s = humanoid_model.getJointId("lleg_shoulder3_joint");
991 
992  transformJointIntoMimic(humanoid_model, index_p, index_s, 2.0, 0.4, humanoid_mimic);
993  data = Data(humanoid_mimic);
994  BOOST_CHECK(humanoid_mimic.check(data));
995 
996  humanoid_mimic = humanoid_mimic.cast<double>();
997  data = Data(humanoid_mimic);
998  BOOST_CHECK(humanoid_mimic.check(data));
999  index_s = humanoid_mimic.getJointId("lleg_shoulder3_joint");
1000  BOOST_CHECK_EQUAL(humanoid_mimic.idx_qs[index_s], humanoid_mimic.joints[index_s].idx_q());
1001  BOOST_CHECK_EQUAL(humanoid_mimic.idx_vs[index_s], humanoid_mimic.joints[index_s].idx_v());
1002 
1003  ModelTpl<float> humanoid_mimic_f = humanoid_mimic.cast<float>();
1004  DataTpl<float> data_f(humanoid_mimic_f);
1005  BOOST_CHECK(humanoid_mimic_f.check(data_f));
1006  index_s = humanoid_mimic_f.getJointId("lleg_shoulder3_joint");
1007  BOOST_CHECK_EQUAL(humanoid_mimic_f.idx_qs[index_s], humanoid_mimic_f.joints[index_s].idx_q());
1008  BOOST_CHECK_EQUAL(humanoid_mimic_f.idx_vs[index_s], humanoid_mimic_f.joints[index_s].idx_v());
1009 }
1010 
1011 BOOST_AUTO_TEST_CASE(test_build_reduced_model_mimic)
1012 {
1013  typedef Model::JointCollection JC;
1014  typedef Model::JointIndex JointIndex;
1017  Model model;
1018 
1019  JointIndex ffidx = model.addJoint(0, JC::JointModelFreeFlyer(), SE3::Identity(), "root_joint");
1020  model.addJointFrame(ffidx);
1021  JointIndex idx = model.addJoint(ffidx, JointModelRX(), SE3::Identity(), "joint1");
1022  model.addJointFrame(idx);
1023  idx = model.addJoint(idx, JointModelRUBX(), SE3::Identity(), "joint2");
1024  model.addJointFrame(idx);
1025  idx = model.addJoint(idx, JointModelRX(), SE3::Identity(), "joint3");
1026  model.addJointFrame(idx);
1027 
1028  auto const mimic = JointModelMimic(
1029  boost::get<JointModelRX>(model.joints[model.getJointId("joint3")].toVariant()), 1., 0.);
1030  idx = model.addJoint(idx, mimic, SE3::Identity(), "joint4");
1031  model.addJointFrame(idx);
1032 
1033  // We lock a 2 DoF joint that is before the mimicked joint.
1034  const std::vector<JointIndex> joints_to_lock = {model.getJointId("joint2")};
1035  Model reduced_model;
1036  BOOST_CHECK_NO_THROW(buildReducedModel(model, joints_to_lock, neutral(model), reduced_model));
1037 }
1038 
1039 BOOST_AUTO_TEST_CASE(test_mimicked_mimicking_vectors)
1040 {
1041  Model man_model, man_mimic;
1042  buildModels::manipulator(man_model);
1043 
1044  // No mimic in model, so vectors should be empty
1045  BOOST_CHECK(man_model.mimicked_joints.size() == 0);
1046  BOOST_CHECK(man_model.mimicking_joints.size() == 0);
1047 
1048  buildMimicModel(man_model, {1, 3}, {2, 5}, {1., 2.}, {0., 1.}, man_mimic);
1049 
1050  // Vectors should be of size 2
1051  BOOST_CHECK(man_mimic.mimicked_joints.size() == 2);
1052  BOOST_CHECK(man_mimic.mimicking_joints.size() == 2);
1053  for (size_t i = 1; i < man_mimic.mimicked_joints.size(); i++)
1054  {
1055  BOOST_CHECK(man_mimic.mimicked_joints[i - 1] < man_mimic.mimicked_joints[i]);
1056  BOOST_CHECK(man_mimic.mimicking_joints[i - 1] < man_mimic.mimicking_joints[i]);
1057  }
1058 }
1059 
1060 BOOST_AUTO_TEST_SUITE_END()
pinocchio::buildModels::manipulator
void manipulator(ModelTpl< Scalar, Options, JointCollectionTpl > &model, const bool mimic=false)
Create a 6-DOF kinematic chain shoulder-elbow-wrist.
pinocchio::FrameIndex
Index FrameIndex
Definition: multibody/fwd.hpp:28
frames.hpp
pinocchio::GeometryObject::overrideMaterial
bool overrideMaterial
Decide whether to override the Material.
Definition: multibody/geometry-object.hpp:119
fwd.hpp
pinocchio::ModelTpl::names
std::vector< std::string > names
Name of the joints.
Definition: multibody/model.hpp:159
Eigen
pinocchio::ModelTpl::frames
FrameVector frames
Vector of operational frames registered on the model.
Definition: multibody/model.hpp:193
pinocchio::JointModelFreeFlyer
JointModelFreeFlyerTpl< context::Scalar > JointModelFreeFlyer
Definition: multibody/joint/fwd.hpp:110
pinocchio::DataTpl
Definition: context/generic.hpp:25
append-urdf-model-with-another-model.model2
model2
Definition: append-urdf-model-with-another-model.py:24
pinocchio::ModelTpl::existJointName
bool existJointName(const std::string &name) const
Check if a joint given by its name exists.
pinocchio::updateGeometryPlacements
void updateGeometryPlacements(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const GeometryModel &geom_model, GeometryData &geom_data, const Eigen::MatrixBase< ConfigVectorType > &q)
Apply a forward kinematics and update the placement of the geometry objects.
pinocchio::GeometryObject::meshScale
Eigen::Vector3d meshScale
Scaling vector applied to the GeometryObject::geometry object.
Definition: multibody/geometry-object.hpp:116
pinocchio::JointModelPlanarTpl::hasConfigurationLimit
const std::vector< bool > hasConfigurationLimit() const
Definition: joint-planar.hpp:580
model.hpp
pinocchio::JointModelPlanarTpl
Definition: multibody/joint/fwd.hpp:118
JointModelFreeFlyer
pinocchio::JointModelFreeFlyerTpl< pinocchio::context::Scalar > JointModelFreeFlyer
Definition: timings.cpp:31
pinocchio::ModelTpl::jointPlacements
SE3Vector jointPlacements
Vector of joint placements: placement of a joint i wrt its parent joint frame.
Definition: multibody/model.hpp:119
pinocchio::transformJointIntoMimic
void transformJointIntoMimic(const ModelTpl< Scalar, Options, JointCollectionTpl > &input_model, const JointIndex &index_mimicked, const JointIndex &index_mimicking, const Scalar &scaling, const Scalar &offset, ModelTpl< Scalar, Options, JointCollectionTpl > &output_model)
Transform of a joint of the model into a mimic joint. Keep the type of the joint as it was previously...
lambdas.parent
parent
Definition: lambdas.py:21
pinocchio::ModelTpl::lowerPositionLimit
ConfigVectorType lowerPositionLimit
Lower joint configuration limit.
Definition: multibody/model.hpp:187
pinocchio::ModelTpl::cast
CastType< NewScalar, ModelTpl >::type cast() const
pinocchio::SE3
context::SE3 SE3
Definition: spatial/fwd.hpp:59
inverse-kinematics.i
int i
Definition: inverse-kinematics.py:17
kinematics.hpp
pinocchio::GeometryObject::meshColor
Eigen::Vector4d meshColor
RGBA color value of the GeometryObject::geometry object.
Definition: multibody/geometry-object.hpp:122
pinocchio::buildMimicModel
void buildMimicModel(const ModelTpl< Scalar, Options, JointCollectionTpl > &input_model, const std::vector< JointIndex > &index_mimicked, const std::vector< JointIndex > &index_mimicking, const std::vector< Scalar > &scaling, const std::vector< Scalar > &offset, ModelTpl< Scalar, Options, JointCollectionTpl > &output_model)
Transform joints of a model into mimic joints.
model.hpp
pinocchio::JointModelRX
JointModelRevoluteTpl< context::Scalar, context::Options, 0 > JointModelRX
Definition: joint-revolute.hpp:877
setup.data
data
Definition: cmake/cython/setup.in.py:48
pinocchio::JointModelPlanarTpl::hasConfigurationLimitInTangent
const std::vector< bool > hasConfigurationLimitInTangent() const
Definition: joint-planar.hpp:585
pinocchio::ModelTpl::joints
JointModelVector joints
Vector of joint models.
Definition: multibody/model.hpp:122
pinocchio::FrameTpl::placement
SE3 placement
Position of kinematic element in parent joint frame.
Definition: model-item.hpp:39
pinocchio::buildReducedModel
void buildReducedModel(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, std::vector< JointIndex > list_of_joints_to_lock, const Eigen::MatrixBase< ConfigVectorType > &reference_configuration, ModelTpl< Scalar, Options, JointCollectionTpl > &reduced_model)
Build a reduced model from a given input model and a list of joint to lock.
pinocchio::ModelTpl::existFrame
bool existFrame(const std::string &name, const FrameType &type=(FrameType)(JOINT|FIXED_JOINT|BODY|OP_FRAME|SENSOR)) const
Checks if a frame given by its name exists.
pinocchio::JointModelPrismaticTpl::hasConfigurationLimit
const std::vector< bool > hasConfigurationLimit() const
Definition: joint-prismatic.hpp:689
pinocchio::JointModelFreeFlyerTpl
Definition: multibody/joint/fwd.hpp:110
pinocchio::GeometryObject
Definition: multibody/geometry-object.hpp:87
pinocchio::GeometryData
Definition: multibody/geometry.hpp:241
pinocchio::buildModels::humanoidRandom
void humanoidRandom(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true, bool mimic=false)
Create a humanoid kinematic tree with 6-DOF limbs and random joint placements.
BOOST_AUTO_TEST_CASE
BOOST_AUTO_TEST_CASE(test_model_subtree)
Definition: unittest/model.cpp:26
autodiff-rnea.f
f
Definition: autodiff-rnea.py:24
pinocchio::GeometryModel::addAllCollisionPairs
void addAllCollisionPairs()
Add all possible collision pairs.
pinocchio::cast
NewScalar cast(const Scalar &value)
Definition: utils/cast.hpp:13
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:315
pinocchio::computeTotalMass
Scalar computeTotalMass(const ModelTpl< Scalar, Options, JointCollectionTpl > &model)
Compute the total mass of the model and return it.
pinocchio::ModelTpl::addJoint
JointIndex addJoint(const JointIndex parent, const JointModel &joint_model, const SE3 &joint_placement, const std::string &joint_name)
Add a joint to the kinematic tree with infinite bounds.
pinocchio::appendModel
void appendModel(const ModelTpl< Scalar, Options, JointCollectionTpl > &modelA, const ModelTpl< Scalar, Options, JointCollectionTpl > &modelB, const FrameIndex frameInModelA, const SE3Tpl< Scalar, Options > &aMb, ModelTpl< Scalar, Options, JointCollectionTpl > &model)
Append a child model into a parent model, after a specific frame given by its index.
pinocchio::ModelTpl::check
bool check(const AlgorithmCheckerBase< D > &checker=AlgorithmCheckerBase< D >()) const
Check the validity of the attributes of Model with respect to the specification of some algorithms.
Definition: multibody/model.hpp:488
pinocchio::FrameTpl
A Plucker coordinate frame attached to a parent joint inside a kinematic tree.
Definition: multibody/frame.hpp:55
pinocchio::Index
PINOCCHIO_COMPILER_DIAGNOSTIC_POP typedef std::size_t Index
Definition: multibody/fwd.hpp:22
pinocchio::findCommonAncestor
JointIndex findCommonAncestor(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, JointIndex joint1_id, JointIndex joint2_id, size_t &index_ancestor_in_support1, size_t &index_ancestor_in_support2)
Computes the common ancestor between two joints belonging to the same kinematic tree.
pinocchio::JointModelPlanar
JointModelPlanarTpl< context::Scalar > JointModelPlanar
Definition: multibody/joint/fwd.hpp:118
pinocchio::ModelTpl::inertias
InertiaVector inertias
Vector of spatial inertias supported by each joint.
Definition: multibody/model.hpp:116
pinocchio::ModelTpl< Scalar >::JointIndex
pinocchio::JointIndex JointIndex
Definition: multibody/model.hpp:67
pinocchio::ModelTpl::referenceConfigurations
ConfigVectorMap referenceConfigurations
Map of reference configurations, indexed by user given names.
Definition: multibody/model.hpp:162
pinocchio::GeometryModel::geometryObjects
GeometryObjectVector geometryObjects
Vector of GeometryObjects used for collision computations.
Definition: multibody/geometry.hpp:221
dpendulum.p
p
Definition: dpendulum.py:7
append-urdf-model-with-another-model.parent_id
int parent_id
Definition: append-urdf-model-with-another-model.py:28
center-of-mass.hpp
joint-configuration.hpp
pinocchio::JointModelPrismaticTpl
Definition: multibody/joint/fwd.hpp:89
geometry.hpp
pinocchio::ModelItem::parentFrame
FrameIndex parentFrame
Index of the parent frame.
Definition: model-item.hpp:36
pinocchio::ModelTpl::mimicking_joints
std::vector< JointIndex > mimicking_joints
Vector of mimicking joints in the tree (with type MimicTpl)
Definition: multibody/model.hpp:151
pinocchio::Data
DataTpl< context::Scalar, context::Options > Data
Definition: multibody/fwd.hpp:34
pinocchio::ModelTpl::idx_qs
std::vector< int > idx_qs
Vector of starting index of the *i*th joint in the configuration space.
Definition: multibody/model.hpp:125
pinocchio::ModelTpl< Scalar >::IndexVector
std::vector< Index > IndexVector
Definition: multibody/model.hpp:70
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.
pinocchio::JointModelTpl< Scalar, Options, JointCollectionTpl >
data.hpp
pinocchio::JointModelRZ
JointModelRevoluteTpl< context::Scalar, context::Options, 2 > JointModelRZ
Definition: joint-revolute.hpp:885
pinocchio::JointModelFreeFlyerTpl::hasConfigurationLimitInTangent
const std::vector< bool > hasConfigurationLimitInTangent() const
Definition: joint-free-flyer.hpp:308
pinocchio::ModelTpl::mimicked_joints
std::vector< JointIndex > mimicked_joints
Vector of mimicked joints in the tree (can be any joint type) The i-th element of this vector corresp...
Definition: multibody/model.hpp:156
pinocchio::GeometryObject::meshPath
std::string meshPath
Absolute path to the mesh file (if the geometry pointee is also a Mesh)
Definition: multibody/geometry-object.hpp:113
pinocchio::Inertia
context::Inertia Inertia
Definition: spatial/fwd.hpp:64
pinocchio::GeometryObject::meshTexturePath
std::string meshTexturePath
Absolute path to the mesh texture file.
Definition: multibody/geometry-object.hpp:130
mimic_dynamics.model_mimic
model_mimic
Definition: mimic_dynamics.py:19
pinocchio::GeometryObject::geometry
CollisionGeometryPtr geometry
The FCL CollisionGeometry (might be a Mesh, a Geometry Primitive, etc.)
Definition: multibody/geometry-object.hpp:110
pos
pos
pinocchio::GeometryModel::ngeoms
Index ngeoms
The number of GeometryObjects.
Definition: multibody/geometry.hpp:218
pinocchio::ModelTpl::addJointFrame
FrameIndex addJointFrame(const JointIndex &joint_index, int previous_frame_index=-1)
Add a joint to the frame tree.
pinocchio::framesForwardKinematics
void framesForwardKinematics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
First calls the forwardKinematics on the model, then computes the placement of each frame....
contact-cholesky.frame
frame
Definition: contact-cholesky.py:24
pinocchio::JointModelPrismaticTpl::hasConfigurationLimitInTangent
const std::vector< bool > hasConfigurationLimitInTangent() const
Definition: joint-prismatic.hpp:694
pinocchio::Frame
FrameTpl< context::Scalar, context::Options > Frame
Definition: multibody/fwd.hpp:31
pinocchio::JOINT
@ JOINT
joint frame: attached to the child body of a joint (a.k.a. child frame)
Definition: multibody/frame.hpp:34
pinocchio::q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
Definition: joint-configuration.hpp:1083
pinocchio::JointModelRUBX
JointModelRevoluteUnboundedTpl< context::Scalar, context::Options, 0 > JointModelRUBX
Definition: joint-revolute-unbounded.hpp:239
contact-cholesky.frame_id
frame_id
Definition: contact-cholesky.py:19
pinocchio::ModelTpl::getJointId
JointIndex getJointId(const std::string &name) const
Return the index of a joint given by its name.
collisions.geom_data
geom_data
Definition: collisions.py:42
pinocchio::buildModels::humanoid
void humanoid(ModelTpl< Scalar, Options, JointCollectionTpl > &model, bool usingFF=true)
Create a 28-DOF kinematic chain of a floating humanoid robot.
check.hpp
pinocchio::OP_FRAME
@ OP_FRAME
operational frame: user-defined frames that are defined at runtime
Definition: multibody/frame.hpp:33
fill
fill
append-urdf-model-with-another-model.model1
model1
Definition: append-urdf-model-with-another-model.py:19
pinocchio::FrameTpl::parentFrame
FrameIndex parentFrame
Index of the parent frame.
Definition: model-item.hpp:36
pinocchio::GeometryModel::collisionPairs
CollisionPairVector collisionPairs
Vector of collision pairs.
Definition: multibody/geometry.hpp:224
pinocchio::ModelTpl< Scalar >::ConfigVectorMap
std::map< std::string, ConfigVectorType > ConfigVectorMap
Map between a string (key) and a configuration vector.
Definition: multibody/model.hpp:90
pinocchio::FIXED_JOINT
@ FIXED_JOINT
fixed joint frame: joint frame but for a fixed joint
Definition: multibody/frame.hpp:35
boost::fusion::append
result_of::push_front< V const, T >::type append(T const &t, V const &v)
Append the element T at the front of boost fusion vector V.
Definition: fusion.hpp:32
pinocchio::JointIndex
Index JointIndex
Definition: multibody/fwd.hpp:26
pinocchio::FrameTpl::type
FrameType type
Type of the frame.
Definition: multibody/frame.hpp:243
pinocchio::JointModelFreeFlyerTpl::hasConfigurationLimit
const std::vector< bool > hasConfigurationLimit() const
Definition: joint-free-flyer.hpp:303
pinocchio::JointModelPX
JointModelPrismaticTpl< context::Scalar, context::Options, 0 > JointModelPX
Definition: joint-prismatic.hpp:781
pinocchio::FrameTpl::name
std::string name
Name of the kinematic element.
Definition: model-item.hpp:25
pinocchio::GeometryModel
Definition: multibody/geometry.hpp:54
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:363
pinocchio::python::context::JointModelRUBX
JointModelRevoluteUnboundedTpl< Scalar, Options, 0 > JointModelRUBX
Definition: bindings/python/context/generic.hpp:85
pinocchio::ModelTpl< Scalar >::ConfigVectorType
VectorXs ConfigVectorType
Dense vectorized version of a joint configuration vector.
Definition: multibody/model.hpp:87
pinocchio::JointModelMimic
JointModelMimicTpl< context::Scalar > JointModelMimic
Definition: multibody/joint/fwd.hpp:155
sample-models.hpp
pinocchio::JointModelRUBZ
JointModelRevoluteUnboundedTpl< context::Scalar, context::Options, 2 > JointModelRUBZ
Definition: joint-revolute-unbounded.hpp:247
pinocchio::ModelTpl::idx_vs
std::vector< int > idx_vs
Starting index of the *i*th joint in the tangent configuration space.
Definition: multibody/model.hpp:131
append-urdf-model-with-another-model.joint_id
joint_id
Definition: append-urdf-model-with-another-model.py:34
pinocchio::ModelTpl< Scalar >
pinocchio::ModelTpl::upperPositionLimit
ConfigVectorType upperPositionLimit
Upper joint configuration limit.
Definition: multibody/model.hpp:190
pinocchio::PINOCCHIO_ALIGNED_STD_VECTOR
typedef PINOCCHIO_ALIGNED_STD_VECTOR(JointData) JointDataVector
pinocchio::ModelItem::name
std::string name
Name of the kinematic element.
Definition: model-item.hpp:25
pinocchio::ModelTpl::nq
int nq
Dimension of the configuration vector representation.
Definition: multibody/model.hpp:98
pinocchio::Model
ModelTpl< context::Scalar, context::Options > Model
Definition: multibody/fwd.hpp:33
n
Vec3f n
pinocchio::ModelTpl< Scalar >::JointCollection
JointCollectionTpl< Scalar, Options > JointCollection
Definition: multibody/model.hpp:57
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1082
pinocchio::python::context::JointModelRX
JointModelRevoluteTpl< Scalar, Options, 0 > JointModelRX
Definition: bindings/python/context/generic.hpp:73
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:33
pinocchio::ModelTpl::njoints
int njoints
Number of joints.
Definition: multibody/model.hpp:107


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