joint-translation.cpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2015-2020 CNRS INRIA
3 // Copyright (c) 2015 Wandercraft, 86 rue de Paris 91400 Orsay, France.
4 //
5 
6 #include "pinocchio/math/fwd.hpp"
13 
14 #include <boost/test/unit_test.hpp>
15 #include <iostream>
16 
17 using namespace pinocchio;
18 
19 template<typename D>
21  Model & model,
22  const JointModelBase<D> & jmodel,
24  const SE3 & joint_placement,
25  const std::string & joint_name,
26  const Inertia & Y)
27 {
29 
30  idx = model.addJoint(parent_id, jmodel, joint_placement, joint_name);
31  model.appendBodyToJoint(idx, Y);
32 }
33 
34 BOOST_AUTO_TEST_SUITE(JointTranslation)
35 
37 {
38  typedef TransformTranslationTpl<double, 0> TransformTranslation;
39  typedef SE3::Vector3 Vector3;
40 
41  const Vector3 displacement(Vector3::Random());
42  SE3 Mplain, Mrand(SE3::Random());
43 
44  TransformTranslation Mtrans(displacement);
45  Mplain = Mtrans;
46  BOOST_CHECK(Mplain.translation().isApprox(displacement));
47  BOOST_CHECK(Mplain.rotation().isIdentity());
48  BOOST_CHECK((Mrand * Mplain).isApprox(Mrand * Mtrans));
49 
50  SE3 M(SE3::Random());
52 
54  Motion mp_dense(mp);
55 
56  BOOST_CHECK(M.act(mp).isApprox(M.act(mp_dense)));
57  BOOST_CHECK(M.actInv(mp).isApprox(M.actInv(mp_dense)));
58 
59  BOOST_CHECK(v.cross(mp).isApprox(v.cross(mp_dense)));
60 }
61 
63 {
64  using namespace pinocchio;
65  typedef SE3::Vector3 Vector3;
66  typedef Eigen::Matrix<double, 6, 1> Vector6;
67  typedef Eigen::Matrix<double, 7, 1> VectorFF;
68  typedef SE3::Matrix3 Matrix3;
69 
70  Model modelTranslation, modelFreeflyer;
71 
72  Inertia inertia(1., Vector3(0.5, 0., 0.0), Matrix3::Identity());
73  SE3 pos(1);
74  pos.translation() = SE3::LinearType(1., 0., 0.);
75 
77  modelTranslation, JointModelTranslation(), 0, SE3::Identity(), "translation", inertia);
78  addJointAndBody(modelFreeflyer, JointModelFreeFlyer(), 0, SE3::Identity(), "free-flyer", inertia);
79 
80  Data dataTranslation(modelTranslation);
81  Data dataFreeFlyer(modelFreeflyer);
82 
83  Eigen::VectorXd q = Eigen::VectorXd::Ones(modelTranslation.nq);
84  VectorFF qff;
85  qff << 1, 1, 1, 0, 0, 0, 1;
86  Eigen::VectorXd v = Eigen::VectorXd::Ones(modelTranslation.nv);
87  Vector6 vff;
88  vff << 1, 1, 1, 0, 0, 0;
89  Eigen::VectorXd tauTranslation = Eigen::VectorXd::Ones(modelTranslation.nv);
90  Eigen::VectorXd tauff(6);
91  tauff << 1, 1, 1, 0, 0, 0;
92  Eigen::VectorXd aTranslation = Eigen::VectorXd::Ones(modelTranslation.nv);
93  Eigen::VectorXd aff(vff);
94 
95  forwardKinematics(modelTranslation, dataTranslation, q, v);
96  forwardKinematics(modelFreeflyer, dataFreeFlyer, qff, vff);
97 
98  computeAllTerms(modelTranslation, dataTranslation, q, v);
99  computeAllTerms(modelFreeflyer, dataFreeFlyer, qff, vff);
100 
101  BOOST_CHECK(dataFreeFlyer.oMi[1].isApprox(dataTranslation.oMi[1]));
102  BOOST_CHECK(dataFreeFlyer.liMi[1].isApprox(dataTranslation.liMi[1]));
103  BOOST_CHECK(dataFreeFlyer.Ycrb[1].matrix().isApprox(dataTranslation.Ycrb[1].matrix()));
104  BOOST_CHECK(dataFreeFlyer.f[1].toVector().isApprox(dataTranslation.f[1].toVector()));
105 
106  Eigen::VectorXd nle_expected_ff(3);
107  nle_expected_ff << dataFreeFlyer.nle[0], dataFreeFlyer.nle[1], dataFreeFlyer.nle[2];
108  BOOST_CHECK(nle_expected_ff.isApprox(dataTranslation.nle));
109  BOOST_CHECK(dataFreeFlyer.com[0].isApprox(dataTranslation.com[0]));
110 
111  // InverseDynamics == rnea
112  tauTranslation = rnea(modelTranslation, dataTranslation, q, v, aTranslation);
113  tauff = rnea(modelFreeflyer, dataFreeFlyer, qff, vff, aff);
114 
115  Vector3 tau_expected;
116  tau_expected << tauff(0), tauff(1), tauff(2);
117  BOOST_CHECK(tauTranslation.isApprox(tau_expected));
118 
119  // ForwardDynamics == aba
120  Eigen::VectorXd aAbaTranslation =
121  aba(modelTranslation, dataTranslation, q, v, tauTranslation, Convention::WORLD);
122  Eigen::VectorXd aAbaFreeFlyer =
123  aba(modelFreeflyer, dataFreeFlyer, qff, vff, tauff, Convention::WORLD);
124  Vector3 a_expected;
125  a_expected << aAbaFreeFlyer[0], aAbaFreeFlyer[1], aAbaFreeFlyer[2];
126  BOOST_CHECK(aAbaTranslation.isApprox(a_expected));
127 
128  // crba
129  crba(modelTranslation, dataTranslation, q, Convention::WORLD);
130  crba(modelFreeflyer, dataFreeFlyer, qff, Convention::WORLD);
131 
132  Eigen::Matrix<double, 3, 3> M_expected(dataFreeFlyer.M.topLeftCorner<3, 3>());
133 
134  BOOST_CHECK(dataTranslation.M.isApprox(M_expected));
135 
136  // Jacobian
137  Eigen::Matrix<double, 6, Eigen::Dynamic> jacobian_translation;
138  jacobian_translation.resize(6, 3);
139  jacobian_translation.setZero();
140  Eigen::Matrix<double, 6, Eigen::Dynamic> jacobian_ff;
141  jacobian_ff.resize(6, 6);
142  jacobian_ff.setZero();
143  computeJointJacobians(modelTranslation, dataTranslation, q);
144  computeJointJacobians(modelFreeflyer, dataFreeFlyer, qff);
145  getJointJacobian(modelTranslation, dataTranslation, 1, LOCAL, jacobian_translation);
146  getJointJacobian(modelFreeflyer, dataFreeFlyer, 1, LOCAL, jacobian_ff);
147 
148  Eigen::Matrix<double, 6, 3> jacobian_expected;
149  jacobian_expected << jacobian_ff.col(0), jacobian_ff.col(1), jacobian_ff.col(2);
150  BOOST_CHECK(jacobian_translation.isApprox(jacobian_expected));
151 }
152 
153 BOOST_AUTO_TEST_SUITE_END()
pinocchio::InertiaTpl< context::Scalar, context::Options >
pinocchio::DataTpl::M
MatrixXs M
The joint space inertia matrix (a square matrix of dim model.nv).
Definition: multibody/data.hpp:199
pinocchio::addJointAndBody
void addJointAndBody(Model &model, const JointModelBase< D > &jmodel, const Model::JointIndex parent_id, const SE3 &joint_placement, const std::string &name, const Inertia &Y)
Definition: model-generator.hpp:11
pinocchio::SE3Tpl< context::Scalar, context::Options >::Vector3
traits< SE3Tpl >::Vector3 Vector3
Definition: spatial/se3-tpl.hpp:55
pinocchio::JointModelFreeFlyer
JointModelFreeFlyerTpl< context::Scalar > JointModelFreeFlyer
Definition: multibody/joint/fwd.hpp:110
pinocchio::DataTpl
Definition: context/generic.hpp:25
pinocchio::JointModelTranslation
JointModelTranslationTpl< context::Scalar > JointModelTranslation
Definition: multibody/joint/fwd.hpp:126
pinocchio::forwardKinematics
void forwardKinematics(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, DataTpl< Scalar, Options, JointCollectionTpl > &data, const Eigen::MatrixBase< ConfigVectorType > &q)
Update the joint placements according to the current joint configuration.
compute-all-terms.hpp
pinocchio::JointModelBase
Definition: joint-model-base.hpp:75
pinocchio::SE3Tpl< context::Scalar, context::Options >
pinocchio::Convention::WORLD
@ WORLD
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::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....
rnea.hpp
aba.hpp
pinocchio::SE3Tpl< context::Scalar, context::Options >::Matrix3
traits< SE3Tpl >::Matrix3 Matrix3
Definition: spatial/se3-tpl.hpp:56
pinocchio::SE3Tpl::translation
ConstLinearRef translation() const
Definition: se3-base.hpp:52
pinocchio::ModelTpl::JointIndex
pinocchio::JointIndex JointIndex
Definition: multibody/model.hpp:68
append-urdf-model-with-another-model.parent_id
int parent_id
Definition: append-urdf-model-with-another-model.py:28
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::SE3Tpl::rotation
ConstAngularRef rotation() const
Definition: se3-base.hpp:48
fwd.hpp
pos
pos
M
M
pinocchio::context::Vector3
Eigen::Matrix< Scalar, 3, 1, Options > Vector3
Definition: context/generic.hpp:53
pinocchio::MotionTpl::Random
static MotionTpl Random()
Definition: motion-tpl.hpp:140
pinocchio::q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
Definition: joint-configuration.hpp:1117
BOOST_AUTO_TEST_CASE
BOOST_AUTO_TEST_CASE(spatial)
Definition: joint-translation.cpp:36
pinocchio::v
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > const Eigen::MatrixBase< TangentVectorType > & v
Definition: joint-configuration.hpp:1118
append-urdf-model-with-another-model.joint_name
string joint_name
Definition: append-urdf-model-with-another-model.py:33
pinocchio::SE3Tpl< context::Scalar, context::Options >::Random
static SE3Tpl Random()
Definition: spatial/se3-tpl.hpp:154
append-urdf-model-with-another-model.joint_placement
joint_placement
Definition: append-urdf-model-with-another-model.py:29
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::ModelTpl::nv
int nv
Dimension of the velocity vector space.
Definition: multibody/model.hpp:102
pinocchio::SE3Tpl< context::Scalar, context::Options >::Identity
static SE3Tpl Identity()
Definition: spatial/se3-tpl.hpp:136
pinocchio::MotionTranslationTpl
Definition: joint-translation.hpp:19
Y
Y
joints.hpp
pinocchio::MotionTpl
Definition: context/casadi.hpp:27
jacobian.hpp
pinocchio::ModelTpl
Definition: context/generic.hpp:20
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...
pinocchio::ModelTpl::nq
int nq
Dimension of the configuration vector representation.
Definition: multibody/model.hpp:99
pinocchio::model
JointCollectionTpl & model
Definition: joint-configuration.hpp:1116
pinocchio::TransformTranslationTpl
Definition: joint-translation.hpp:196
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:27
pinocchio::DataTpl::nle
VectorXs nle
Vector of Non Linear Effects (dim model.nv). It corresponds to concatenation of the Coriolis,...
Definition: multibody/data.hpp:179
BOOST_CHECK
#define BOOST_CHECK(check)
Definition: overview-urdf.cpp:34
isApprox
bool isApprox(const Box &s1, const Box &s2, const FCL_REAL tol)
pinocchio::LOCAL
@ LOCAL
Definition: multibody/fwd.hpp:50


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