5 #include "pinocchio/multibody/liegroup/liegroup.hpp" 6 #include "pinocchio/multibody/liegroup/liegroup-collection.hpp" 7 #include "pinocchio/multibody/liegroup/liegroup-generic.hpp" 8 #include "pinocchio/multibody/liegroup/cartesian-product-variant.hpp" 10 #include "pinocchio/multibody/joint/joint-generic.hpp" 12 #include <boost/test/unit_test.hpp> 13 #include <boost/utility/binary.hpp> 14 #include <boost/algorithm/string.hpp> 16 #define EIGEN_VECTOR_IS_APPROX(Va, Vb, precision) \ 17 BOOST_CHECK_MESSAGE((Va).isApprox(Vb, precision), \ 18 "check " #Va ".isApprox(" #Vb ") failed " \ 19 "[\n" << (Va).transpose() << "\n!=\n" << (Vb).transpose() << "\n]") 20 #define EIGEN_MATRIX_IS_APPROX(Va, Vb, precision) \ 21 BOOST_CHECK_MESSAGE((Va).isApprox(Vb, precision), \ 22 "check " #Va ".isApprox(" #Vb ") failed " \ 23 "[\n" << (Va) << "\n!=\n" << (Vb) << "\n]") 28 #define IFVERBOSE if(VERBOSE) 31 template<
typename Derived>
32 std::ostream& operator<< (std::ostream& os, const LieGroupBase<Derived>& lg)
34 return os << lg.name();
36 template<
typename LieGroupCollection>
37 std::ostream& operator<< (std::ostream& os, const LieGroupGenericTpl<LieGroupCollection>& lg)
39 return os << lg.name();
48 const Scalar prec = Eigen::NumTraits<Scalar>::dummy_precision();
49 BOOST_TEST_MESSAGE (
"Testing Joint over " << jmodel.shortname());
50 typedef typename T::ConfigVector_t ConfigVector_t;
51 typedef typename T::TangentVector_t TangentVector_t;
53 ConfigVector_t
q1(ConfigVector_t::Random (jmodel.nq()));
54 TangentVector_t q1_dot(TangentVector_t::Random (jmodel.nv()));
55 ConfigVector_t
q2(ConfigVector_t::Random (jmodel.nq()));
58 static ConfigVector_t Ones(ConfigVector_t::Ones(jmodel.nq()));
63 q1 = LieGroupType().randomConfiguration(-Ones, Ones);
65 typename T::JointDataDerived jdata = jmodel.createData();
68 jmodel.calc(jdata, q1, q1_dot);
72 q2 = LieGroupType().integrate(q1,q1_dot);
73 jmodel.calc(jdata,q2);
78 if(jmodel.shortname() !=
"JointModelSphericalZYX")
80 BOOST_CHECK_MESSAGE(M2.
isApprox(M2_exp), std::string(
"Error when integrating1 " + jmodel.shortname()));
84 ConfigVector_t q3 = LieGroupType().integrate(q2,-q1_dot);
85 jmodel.calc(jdata,q3);
88 BOOST_CHECK_MESSAGE(M3.
isApprox(M1), std::string(
"Error when integrating back " + jmodel.shortname()));
91 ConfigVector_t q_interpolate = LieGroupType().interpolate(q1,q2,0.);
92 BOOST_CHECK_MESSAGE(q_interpolate.isApprox(q1), std::string(
"Error when interpolating " + jmodel.shortname()));
94 q_interpolate = LieGroupType().interpolate(q1,q2,1.);
95 BOOST_CHECK_MESSAGE(q_interpolate.isApprox(q2), std::string(
"Error when interpolating " + jmodel.shortname()));
97 if(jmodel.shortname() !=
"JointModelSphericalZYX")
99 q_interpolate = LieGroupType().interpolate(q1,q2,u);
100 jmodel.calc(jdata,q_interpolate);
101 SE3 M_interpolate = jdata.M;
103 SE3 M_interpolate_expected = M1*
exp6(u*v1);
104 BOOST_CHECK_MESSAGE(M_interpolate_expected.
isApprox(M_interpolate,1e2*prec), std::string(
"Error when interpolating " + jmodel.shortname()));
108 TangentVector_t
zero = LieGroupType().difference(q1,q1);
109 BOOST_CHECK_MESSAGE (zero.isZero(), std::string (
"Error: difference between two equal configurations is not 0."));
110 zero = LieGroupType().difference(q2,q2);
111 BOOST_CHECK_MESSAGE (zero.isZero(), std::string (
"Error: difference between two equal configurations is not 0."));
114 TangentVector_t vdiff = LieGroupType().difference(q1,q2);
115 BOOST_CHECK_MESSAGE(vdiff.isApprox(q1_dot,1e2*prec), std::string(
"Error when differentiating " + jmodel.shortname()));
118 Scalar dist = LieGroupType().distance(q1,q2);
119 BOOST_CHECK_MESSAGE(dist > 0.,
"distance - wrong results");
120 BOOST_CHECK_SMALL(math::fabs(dist-q1_dot.norm()), 10*prec);
122 std::string error_prefix(
"LieGroup");
123 error_prefix +=
" on joint " + jmodel.shortname();
125 BOOST_CHECK_MESSAGE(jmodel.nq() ==
LieGroupType::NQ, std::string(error_prefix +
" - nq "));
126 BOOST_CHECK_MESSAGE(jmodel.nv() ==
LieGroupType::NV, std::string(error_prefix +
" - nv "));
130 LieGroupType().randomConfiguration(-1 * Ones, Ones).size(),
131 std::string(error_prefix +
" - RandomConfiguration dimensions "));
133 ConfigVector_t q_normalize(ConfigVector_t::Random());
135 if(jmodel.shortname() ==
"JointModelSpherical")
137 BOOST_CHECK_MESSAGE(!LieGroupType().
isNormalized(q_normalize_ref), std::string(error_prefix +
" - !isNormalized "));
138 q_normalize_ref /= q_normalize_ref.norm();
140 else if(jmodel.shortname() ==
"JointModelFreeFlyer")
142 BOOST_CHECK_MESSAGE(!LieGroupType().
isNormalized(q_normalize_ref), std::string(error_prefix +
" - !isNormalized "));
143 q_normalize_ref.template tail<4>() /= q_normalize_ref.template tail<4>().norm();
145 else if(boost::algorithm::istarts_with(jmodel.shortname(),
"JointModelRUB"))
147 BOOST_CHECK_MESSAGE(!LieGroupType().
isNormalized(q_normalize_ref), std::string(error_prefix +
" - !isNormalized "));
148 q_normalize_ref /= q_normalize_ref.norm();
150 else if(jmodel.shortname() ==
"JointModelPlanar")
152 BOOST_CHECK_MESSAGE(!LieGroupType().
isNormalized(q_normalize_ref), std::string(error_prefix +
" - !isNormalized "));
153 q_normalize_ref.template tail<2>() /= q_normalize_ref.template tail<2>().norm();
155 BOOST_CHECK_MESSAGE(LieGroupType().
isNormalized(q_normalize_ref), std::string(error_prefix +
" - isNormalized "));
156 LieGroupType().normalize(q_normalize);
157 BOOST_CHECK_MESSAGE(q_normalize.isApprox(q_normalize_ref), std::string(error_prefix +
" - normalize "));
162 template <
typename T>
166 jmodel.setIndexes(0,0,0);
167 typename T::JointDataDerived jdata = jmodel.createData();
176 pinocchio::JointModelRevoluteUnaligned::JointDataDerived jdata = jmodel.
createData();
185 pinocchio::JointModelPrismaticUnaligned::JointDataDerived jdata = jmodel.
createData();
193 template <
typename T>
196 typedef typename T::ConfigVector_t ConfigVector_t;
197 typedef typename T::TangentVector_t TangentVector_t;
198 typedef typename T::JacobianMatrix_t JacobianMatrix_t;
202 BOOST_TEST_MESSAGE (lg.name());
203 ConfigVector_t
q[2], q_dv[2];
206 TangentVector_t va, vb, dv;
207 JacobianMatrix_t
J[2];
210 lg.difference (q[0], q[1], va);
211 lg.template dDifference<ARG0> (q[0], q[1], J[0]);
212 lg.template dDifference<ARG1> (q[0], q[1], J[1]);
214 const Scalar
eps = 1e-6;
215 for (
int k = 0; k < 2; ++k) {
216 BOOST_TEST_MESSAGE (
"Checking J" << k <<
'\n' << J[k]);
220 for (
int i = 0;
i < dv.size(); ++
i)
223 lg.integrate (q[k], dv, q_dv[k]);
224 lg.difference (q_dv[0], q_dv[1], vb);
227 TangentVector_t J_dv = J[k].col(
i);
228 TangentVector_t vb_va = (vb - va) / eps;
237 template <
typename T>
241 template <
typename Scalar,
int Options>
246 typedef typename LG_t::ConfigVector_t ConfigVector_t;
247 typedef typename LG_t::JacobianMatrix_t JacobianMatrix_t;
248 typedef typename LG_t::ConstQuaternionMap_t ConstQuaternionMap_t;
256 ConstQuaternionMap_t quat0(q[0].
template tail<4>().
data()), quat1(q[1].
template tail<4>().
data());
257 JacobianMatrix_t
J[2];
259 lg.template dDifference<ARG0> (q[0], q[1], J[0]);
260 lg.template dDifference<ARG1> (q[0], q[1], J[1]);
264 _1m2 (om1.actInv (om0)) ;
269 ConfigVector_t q_interp = lg.interpolate(q[0],q[1],u);
270 ConstQuaternionMap_t quat_interp(q_interp.template tail<4>().data());
272 SE3 M0(quat0,q[0].
template head<3>());
273 SE3 M1(quat1,q[1].
template head<3>());
276 SE3 M_interp(quat_interp,q_interp.template head<3>());
278 BOOST_CHECK(M_u.
isApprox(M_interp));
281 template <
typename Scalar,
int Options>
292 typedef typename LG_t::ConfigVector_t ConfigVector_t;
293 typedef typename LG_t::JacobianMatrix_t JacobianMatrix_t;
300 JacobianMatrix_t
J[2];
302 lg.template dDifference<ARG0> (q[0], q[1], J[0]);
303 lg.template dDifference<ARG1> (q[0], q[1], J[1]);
308 JacobianMatrix_t X (JacobianMatrix_t::Identity());
309 X.template bottomRightCorner<3,3>() = oR1.transpose() * oR0;
314 template<
bool around_
identity>
316 template <
typename T>
319 typedef typename T::ConfigVector_t ConfigVector_t;
320 typedef typename T::TangentVector_t TangentVector_t;
321 typedef typename T::JacobianMatrix_t JacobianMatrix_t;
325 ConfigVector_t
q = lg.random();
326 TangentVector_t
v,
dq, dv;
335 ConfigVector_t q_v = lg.integrate (q, v);
337 JacobianMatrix_t Jq, Jv;
338 lg.dIntegrate_dq (q, v, Jq);
339 lg.dIntegrate_dv (q, v, Jv);
341 const Scalar
eps = 1e-6;
342 for (
int i = 0;
i < v.size(); ++
i)
345 ConfigVector_t q_dq = lg.integrate (q, dq);
347 ConfigVector_t q_dq_v = lg.integrate (q_dq, v);
348 TangentVector_t Jq_dq = Jq.col(
i);
350 TangentVector_t dI_dq = lg.difference (q_v, q_dq_v) /
eps;
353 ConfigVector_t q_v_dv = lg.integrate (q, (v+dv).eval());
354 TangentVector_t Jv_dv = Jv.col(
i);
356 TangentVector_t dI_dv = lg.difference (q_v, q_v_dv) /
eps;
365 template <
typename T>
368 typedef typename T::ConfigVector_t ConfigVector_t;
369 typedef typename T::TangentVector_t TangentVector_t;
370 typedef typename T::JacobianMatrix_t JacobianMatrix_t;
373 BOOST_TEST_MESSAGE (lg.name());
374 ConfigVector_t qa, qb (lg.nq());
376 TangentVector_t
v (lg.nv());
378 lg.integrate(qa,
v, qb);
380 JacobianMatrix_t Jd_qb, Ji_v;
382 lg.template dDifference<ARG1> (qa, qb, Jd_qb);
383 lg.template dIntegrate <ARG1> (qa,
v , Ji_v );
385 BOOST_CHECK_MESSAGE ((Jd_qb * Ji_v).isIdentity(),
391 Jd_qb * Ji_v <<
'\n');
397 template <
typename T>
400 typedef typename T::ConfigVector_t ConfigVector_t;
401 typedef typename T::TangentVector_t TangentVector_t;
405 ConfigVector_t
q = lg.random();
406 TangentVector_t dv(TangentVector_t::Zero(lg.nv()));
408 BOOST_TEST_MESSAGE (lg.name());
409 typedef Eigen::Matrix<Scalar,T::NQ,T::NV> JacobianCoeffs;
410 JacobianCoeffs Jintegrate(JacobianCoeffs::Zero(lg.nq(),lg.nv()));
411 lg.integrateCoeffWiseJacobian(q,Jintegrate);
412 JacobianCoeffs Jintegrate_fd(JacobianCoeffs::Zero(lg.nq(),lg.nv()));
414 const Scalar
eps = 1e-8;
415 for (
int i = 0;
i < lg.nv(); ++
i)
418 ConfigVector_t q_next(ConfigVector_t::Zero(lg.nq()));
419 lg.integrate(q, dv,q_next);
420 Jintegrate_fd.col(
i) = (q_next -
q)/eps;
429 BOOST_AUTO_TEST_SUITE ( BOOST_TEST_MODULE )
442 for (
int i = 0;
i < 20; ++
i)
443 boost::mpl::for_each<Variant::types>(
TestJoint());
455 typedef boost::mpl::vector< VectorSpaceOperationTpl<1,Scalar,Options>
462 VectorSpaceOperationTpl<2,Scalar,Options>,
463 SpecialOrthogonalOperationTpl<2,Scalar,Options>
467 SpecialOrthogonalOperationTpl<3,Scalar,Options>
470 for (
int i = 0;
i < 20; ++
i)
479 typedef boost::mpl::vector< VectorSpaceOperationTpl<1,Scalar,Options>
486 VectorSpaceOperationTpl<2,Scalar,Options>,
487 SpecialOrthogonalOperationTpl<2,Scalar,Options>
491 SpecialOrthogonalOperationTpl<3,Scalar,Options>
494 for (
int i = 0;
i < 20; ++
i)
506 typedef boost::mpl::vector< VectorSpaceOperationTpl<1,Scalar,Options>
513 VectorSpaceOperationTpl<2,Scalar,Options>,
514 SpecialOrthogonalOperationTpl<2,Scalar,Options>
518 SpecialOrthogonalOperationTpl<3,Scalar,Options>
521 for (
int i = 0;
i < 20; ++
i)
530 typedef boost::mpl::vector< VectorSpaceOperationTpl<1,Scalar,Options>
537 VectorSpaceOperationTpl<2,Scalar,Options>,
538 SpecialOrthogonalOperationTpl<2,Scalar,Options>
542 SpecialOrthogonalOperationTpl<3,Scalar,Options>
545 for (
int i = 0;
i < 20; ++
i)
549 typedef SpecialEuclideanOperationTpl<3,Scalar,Options>
LieGroup;
550 typedef LieGroup::ConfigVector_t ConfigVector_t;
553 ConfigVector_t
q = lg.random();
556 typedef Eigen::Matrix<Scalar,LieGroup::NQ,LieGroup::NV> JacobianCoeffs;
557 JacobianCoeffs Jintegrate(JacobianCoeffs::Zero(lg.nq(),lg.nv()));
558 lg.integrateCoeffWiseJacobian(q,Jintegrate);
568 VSO_t::ConfigVector_t
q,
569 lo(VSO_t::ConfigVector_t::Constant(-std::numeric_limits<double>::infinity())),
572 up(VSO_t::ConfigVector_t::Constant( 0));
576 VSO_t ().randomConfiguration(lo, up, q);
577 }
catch (
const std::runtime_error&) {
580 BOOST_CHECK_MESSAGE(error,
"Random configuration between infinite bounds should return an error");
590 BOOST_CHECK (vs1.
nq () == 1);
591 BOOST_CHECK (vs1.
nv () == 1);
592 BOOST_CHECK (vs1.
name () ==
"R^1");
598 BOOST_CHECK (vs2.nq () == 2);
599 BOOST_CHECK (vs2.nv () == 2);
600 BOOST_CHECK (vs2.name () ==
"R^2");
601 BOOST_CHECK (vs2.neutral () ==
neutral);
606 BOOST_CHECK (vs3.nq () == 3);
607 BOOST_CHECK (vs3.nv () == 3);
608 BOOST_CHECK (vs3.name () ==
"R^3");
609 BOOST_CHECK (vs3.neutral () ==
neutral);
612 neutral.resize (2); neutral [0] = 1; neutral [1] = 0;
613 BOOST_CHECK (so2.nq () == 2);
614 BOOST_CHECK (so2.nv () == 1);
615 BOOST_CHECK (so2.name () ==
"SO(2)");
616 BOOST_CHECK (so2.neutral () ==
neutral);
619 neutral.resize (4); neutral.setZero ();
621 BOOST_CHECK (so3.nq () == 4);
622 BOOST_CHECK (so3.nv () == 3);
623 BOOST_CHECK (so3.name () ==
"SO(3)");
624 BOOST_CHECK (so3.neutral () ==
neutral);
627 neutral.resize (4); neutral.setZero ();
629 BOOST_CHECK (se2.nq () == 4);
630 BOOST_CHECK (se2.nv () == 3);
631 BOOST_CHECK (se2.name () ==
"SE(2)");
632 BOOST_CHECK (se2.neutral () ==
neutral);
635 neutral.resize (7); neutral.setZero ();
637 BOOST_CHECK (se3.nq () == 7);
638 BOOST_CHECK (se3.nv () == 6);
639 BOOST_CHECK (se3.name () ==
"SE(3)");
640 BOOST_CHECK (se3.neutral () ==
neutral);
644 neutral.resize (4); neutral.setZero ();
646 BOOST_CHECK (r2xso2.nq () == 4);
647 BOOST_CHECK (r2xso2.nv () == 3);
648 BOOST_CHECK (r2xso2.name () ==
"R^2*SO(2)");
649 BOOST_CHECK (r2xso2.neutral () ==
neutral);
653 neutral.resize (7); neutral.setZero ();
655 BOOST_CHECK (r3xso3.nq () == 7);
656 BOOST_CHECK (r3xso3.nv () == 6);
657 BOOST_CHECK (r3xso3.name () ==
"R^3*SO(3)");
658 BOOST_CHECK (r3xso3.neutral () ==
neutral);
664 BOOST_CHECK(dim == 2);
666 BOOST_CHECK(dim == Eigen::Dynamic);
668 BOOST_CHECK(dim == Eigen::Dynamic);
671 template<
typename LieGroupCollection>
679 template<
typename Derived>
686 template<
typename Derived>
690 BOOST_CHECK(lg.
nq() ==
nq(lg_generic));
691 BOOST_CHECK(lg.
nv() ==
nv(lg_generic));
693 BOOST_CHECK(lg.
name() ==
name(lg_generic));
701 TangentVector_t
v = TangentVector_t::Random(lg.
nv());
702 ConfigVector_t qout_ref(lg.
nq());
705 ConfigVectorGeneric qout(lg.
nq());
706 integrate(lg_generic, ConfigVectorGeneric(q0), TangentVectorGeneric(v), qout);
707 BOOST_CHECK(qout.isApprox(qout_ref));
709 ConfigVector_t
q1 (
nq(lg_generic));
714 ConfigVector_t
q2(
nq(lg_generic));
726 template<
typename Lg1,
typename Lg2>
730 BOOST_CHECK_EQUAL(LieGroupGeneric(lg1), LieGroupGeneric(lg2));
733 template<
typename Lg1,
typename Lg2>
737 BOOST_CHECK_PREDICATE( std::not_equal_to<LieGroupGeneric>(),
738 (LieGroupGeneric(lg1))(LieGroupGeneric(lg2)) );
751 BOOST_AUTO_TEST_SUITE_END ()
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > const Eigen::MatrixBase< TangentVectorType > & v
Scalar distance(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorIn1 > &q0, const Eigen::MatrixBase< ConfigVectorIn2 > &q1)
Distance between two configuration vectors, namely .
void test_liegroup_variant_not_equal(Lg1 lg1, Lg2 lg2)
JointModelPrismaticTpl< double, 0, 0 > JointModelPX
int nv(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointNvVisitor to get the dimension of the joint tangent space...
ConfigVector_t neutral() const
Get neutral element as a vector.
Eigen::Matrix< Scalar, NQ, 1, Options > ConfigVector_t
void operator()(const T) const
bool isApprox(const Derived &other, const Scalar &prec=Eigen::NumTraits< Scalar >::dummy_precision()) const
void setIndexes(JointIndex id, int q, int v)
boost::python::object matrix()
int nq(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointNqVisitor to get the dimension of the joint configuration space...
void integrate(const Eigen::MatrixBase< ConfigIn_t > &q, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< ConfigOut_t > &qout) const
Integrate a joint's configuration with a tangent vector during one unit time duration.
void integrate(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q, const Eigen::MatrixBase< TangentVectorType > &v, const Eigen::MatrixBase< ReturnType > &qout)
Integrate a configuration vector for the specified model for a tangent vector during one unit time...
ConfigVector_t neutral() const
void difference(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorIn1 > &q0, const Eigen::MatrixBase< ConfigVectorIn2 > &q1, const Eigen::MatrixBase< ReturnType > &dvout)
Compute the tangent vector that must be integrated during one unit time to go from q0 to q1...
void operator()(const T) const
JointModelRevoluteUnboundedTpl< double, 0, 1 > JointModelRUBY
LieGroupGeneric::TangentVector_t TangentVector_t
JointDataDerived createData() const
JointModelRevoluteUnalignedTpl< double > JointModelRevoluteUnaligned
JointModelPrismaticTpl< double, 0, 1 > JointModelPY
Eigen::Matrix< Scalar, NV, 1, Options > TangentVector_t
std::string name(const LieGroupGenericTpl< LieGroupCollection > &lg)
Visit a LieGroupVariant to get the name of it.
bool isNormalized(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &q, const Scalar &prec=Eigen::NumTraits< Scalar >::dummy_precision())
Check whether a configuration vector is normalized within the given precision provided by prec...
static void test(const LieGroupBase< Derived > &lg, const LieGroupGenericTpl< LieGroupCollection > &lg_generic)
LieGroupGeneric::ConfigVector_t ConfigVector_t
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
void test_liegroup_variant_equal(Lg1 lg1, Lg2 lg2)
void specificTests(const CartesianProductOperation< VectorSpaceOperationTpl< 3, Scalar, Options >, SpecialOrthogonalOperationTpl< 3, Scalar, Options > >) const
void random(const Eigen::MatrixBase< Config_t > &qout) const
Generate a random joint configuration, normalizing quaternions when necessary.
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorIn1 > const Eigen::MatrixBase< ConfigVectorIn2 > const Scalar & u
JointModelRevoluteTpl< double, 0, 1 > JointModelRY
JointModelSphericalZYXTpl< double > JointModelSphericalZYX
Eigen::Quaternion< Scalar, Options > Quaternion
JointModelRevoluteTpl< double, 0, 2 > JointModelRZ
BOOST_AUTO_TEST_CASE(test_all)
#define EIGEN_VECTOR_IS_APPROX(Va, Vb, precision)
void specificTests(const T) const
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorIn1 > const Eigen::MatrixBase< ConfigVectorIn2 > & q1
void operator()(const T) const
void operator()(const pinocchio::JointModelRevoluteUnaligned &) const
JointModelTranslationTpl< double > JointModelTranslation
void specificTests(const SpecialEuclideanOperationTpl< 3, Scalar, Options >) const
Main pinocchio namespace.
JointModelSphericalTpl< double > JointModelSpherical
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorIn1 > & q0
void neutral(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ReturnType > &qout)
Return the neutral configuration element related to the model configuration space.
void operator()(const pinocchio::JointModelPrismaticUnaligned &) const
JointModelPrismaticTpl< double, 0, 2 > JointModelPZ
void random(const LieGroupGenericTpl< LieGroupCollection > &lg, const Eigen::MatrixBase< Config_t > &qout)
#define EIGEN_MATRIX_IS_APPROX(Va, Vb, precision)
void operator()(const T) const
std::string name() const
Get name of instance.
JointModelFreeFlyerTpl< double > JointModelFreeFlyer
void normalize(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &qout)
Normalize a configuration vector.
void test_lie_group_methods(T &jmodel, typename T::JointDataDerived &)
SE3Tpl< typename MotionDerived::Scalar, PINOCCHIO_EIGEN_PLAIN_TYPE(typename MotionDerived::Vector3)::Options > exp6(const MotionDense< MotionDerived > &nu)
Exp: se3 -> SE3.
void operator()(const T) const
LieGroupGenericTpl< LieGroupCollection > LieGroupGeneric
traits< SE3Tpl >::Matrix3 Matrix3
Index nv() const
Get dimension of Lie Group tangent space.
JointDataDerived createData() const
JointModelRevoluteUnboundedTpl< double, 0, 0 > JointModelRUBX
JointModelPrismaticUnalignedTpl< double > JointModelPrismaticUnaligned
Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > VectorXd
JointModelRevoluteTpl< double, 0, 0 > JointModelRX
static SE3Tpl Interpolate(const SE3Tpl &A, const SE3Tpl &B, const OtherScalar &alpha)
Linear interpolation on the SE3 manifold.
JointModelPlanarTpl< double > JointModelPlanar
Scalar distance(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1) const
Distance between two configurations of the joint.