pinocchio  2.1.3
joint-spherical-ZYX.hpp
1 //
2 // Copyright (c) 2015-2018 CNRS
3 // Copyright (c) 2015-2016 Wandercraft, 86 rue de Paris 91400 Orsay, France.
4 //
5 
6 #ifndef __pinocchio_joint_spherical_ZYX_hpp__
7 #define __pinocchio_joint_spherical_ZYX_hpp__
8 #include <iostream>
9 #include "pinocchio/macros.hpp"
10 #include "pinocchio/multibody/joint/joint-base.hpp"
11 #include "pinocchio/multibody/joint/joint-spherical.hpp"
12 #include "pinocchio/multibody/constraint.hpp"
13 #include "pinocchio/math/sincos.hpp"
14 #include "pinocchio/spatial/inertia.hpp"
15 #include "pinocchio/spatial/skew.hpp"
16 
17 namespace pinocchio
18 {
19  template<typename Scalar, int Options> struct ConstraintSphericalZYXTpl;
20 
21  template <typename _Scalar, int _Options>
22  struct traits< ConstraintSphericalZYXTpl<_Scalar,_Options> >
23  {
24  typedef _Scalar Scalar;
25  enum { Options = _Options };
26  typedef Eigen::Matrix<Scalar,3,1,Options> Vector3;
27  typedef Eigen::Matrix<Scalar,4,1,Options> Vector4;
28  typedef Eigen::Matrix<Scalar,6,1,Options> Vector6;
29  typedef Eigen::Matrix<Scalar,3,3,Options> Matrix3;
30  typedef Eigen::Matrix<Scalar,4,4,Options> Matrix4;
31  typedef Eigen::Matrix<Scalar,6,6,Options> Matrix6;
32  typedef Matrix3 Angular_t;
33  typedef Vector3 Linear_t;
34  typedef const Matrix3 ConstAngular_t;
35  typedef const Vector3 ConstLinear_t;
36  typedef Matrix6 ActionMatrix_t;
37  typedef Eigen::Quaternion<Scalar,Options> Quaternion_t;
42  enum {
43  LINEAR = 0,
44  ANGULAR = 3
45  };
47  typedef Eigen::Matrix<Scalar,3,1,Options> JointForce;
48  typedef Eigen::Matrix<Scalar,6,3,Options> DenseBase;
49 
50  typedef DenseBase MatrixReturnType;
51  typedef const DenseBase ConstMatrixReturnType;
52  }; // struct traits struct ConstraintRotationalSubspace
53 
54  template<typename _Scalar, int _Options>
55  struct ConstraintSphericalZYXTpl : public ConstraintBase< ConstraintSphericalZYXTpl<_Scalar,_Options> >
56  {
57  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
58 
59  SPATIAL_TYPEDEF_TEMPLATE(ConstraintSphericalZYXTpl);
60 
61  enum { NV = 3, Options = _Options };
63  typedef typename traits<ConstraintSphericalZYXTpl>::JointForce JointForce;
64  typedef typename traits<ConstraintSphericalZYXTpl>::DenseBase DenseBase;
65 
67 
68  template<typename Matrix3Like>
69  ConstraintSphericalZYXTpl(const Eigen::MatrixBase<Matrix3Like> & subspace)
70  : S_minimal(subspace)
71  { EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix3Like,3,3); }
72 
73  template<typename Vector3Like>
74  JointMotion __mult__(const Eigen::MatrixBase<Vector3Like> & v) const
75  {
76  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3Like,3);
77  return JointMotion(S_minimal * v);
78  }
79 
80  Matrix3 & operator() () { return S_minimal; }
81  const Matrix3 & operator() () const { return S_minimal; }
82 
83  int nv_impl() const { return NV; }
84 
86  {
87  const ConstraintSphericalZYXTpl & ref;
88  ConstraintTranspose(const ConstraintSphericalZYXTpl & ref) : ref(ref) {}
89 
90  template<typename Derived>
91 #if EIGEN_VERSION_AT_LEAST(3,2,90)
92  const typename Eigen::Product<
93  Eigen::Transpose<const Matrix3>,
94  Eigen::Block<const typename ForceDense<Derived>::Vector6,3,1>
95  >
96 #else
97  const typename Eigen::ProductReturnType<
98  Eigen::Transpose<const Matrix3>,
99  // typename Motion::ConstAngular_t::Base /* This feature leads currently to a bug */
100  Eigen::Block<const typename ForceDense<Derived>::Vector6,3,1>
101  >::Type
102 #endif
103  operator* (const ForceDense<Derived> & phi) const
104  {
105  return ref.S_minimal.transpose() * phi.angular();
106  }
107 
108  /* [CRBA] MatrixBase operator* (Constraint::Transpose S, ForceSet::Block) */
109  template<typename D>
110 #if EIGEN_VERSION_AT_LEAST(3,2,90)
111  const typename Eigen::Product<
112  typename Eigen::Transpose<const Matrix3>,
113  typename Eigen::MatrixBase<const D>::template NRowsBlockXpr<3>::Type
114  >
115 #else
116  const typename Eigen::ProductReturnType<
117  typename Eigen::Transpose<const Matrix3>,
118  typename Eigen::MatrixBase<const D>::template NRowsBlockXpr<3>::Type
119  >::Type
120 #endif
121  operator* (const Eigen::MatrixBase<D> & F) const
122  {
123  EIGEN_STATIC_ASSERT(D::RowsAtCompileTime==6,THIS_METHOD_IS_ONLY_FOR_MATRICES_OF_A_SPECIFIC_SIZE)
124  return ref.S_minimal.transpose () * F.template middleRows<3>(ANGULAR);
125  }
126  }; // struct ConstraintTranspose
127 
128  ConstraintTranspose transpose () const { return ConstraintTranspose(*this); }
129 
130  DenseBase matrix_impl() const
131  {
132  DenseBase S;
133  S.template middleRows<3>(LINEAR).setZero();
134  S.template middleRows<3>(ANGULAR) = S_minimal;
135  return S;
136  }
137 
138  // const typename Eigen::ProductReturnType<
139  // const ConstraintDense,
140  // const Matrix3
141  // >::Type
142  template<typename S1, int O1>
143  Eigen::Matrix<Scalar,6,3,Options>
144  se3Action(const SE3Tpl<S1,O1> & m) const
145  {
146  // Eigen::Matrix <Scalar,6,3,Options> X_subspace;
147  // X_subspace.template block <3,3> (Motion::LINEAR, 0) = skew (m.translation ()) * m.rotation ();
148  // X_subspace.template block <3,3> (Motion::ANGULAR, 0) = m.rotation ();
149  //
150  // return (X_subspace * S_minimal).eval();
151 
152  Eigen::Matrix<Scalar,6,3,Options> result;
153  result.template middleRows<3>(ANGULAR).noalias() = m.rotation () * S_minimal;
154  for(int k = 0; k < 3; ++k)
155  result.template middleRows<3>(LINEAR).col(k) =
156  m.translation().cross(result.template middleRows<3>(Motion::ANGULAR).col(k));
157 
158  return result;
159  }
160 
161  template<typename MotionDerived>
162  DenseBase motionAction(const MotionDense<MotionDerived> & m) const
163  {
164  const typename MotionDerived::ConstLinearType v = m.linear();
165  const typename MotionDerived::ConstAngularType w = m.angular();
166 
167  DenseBase res;
168  cross(v,S_minimal,res.template middleRows<3>(LINEAR));
169  cross(w,S_minimal,res.template middleRows<3>(ANGULAR));
170 
171  return res;
172  }
173 
174  // data
175  Matrix3 S_minimal;
176 
177  }; // struct ConstraintSphericalZYXTpl
178 
179  /* [CRBA] ForceSet operator* (Inertia Y,Constraint S) */
180  template <typename S1, int O1, typename S2, int O2>
181  Eigen::Matrix<S1,6,3,O1>
182  operator*(const InertiaTpl<S1,O1> & Y,
184  {
185  typedef typename InertiaTpl<S1,O1>::Symmetric3 Symmetric3;
186  typedef ConstraintSphericalZYXTpl<S2,O2> Constraint;
187  Eigen::Matrix<S1,6,3,O1> M;
188  alphaSkew (-Y.mass(),Y.lever(),M.template middleRows<3>(Constraint::LINEAR));
189  M.template middleRows<3>(Constraint::ANGULAR) = (Y.inertia () -
190  typename Symmetric3::AlphaSkewSquare(Y.mass (), Y.lever ())).matrix();
191 
192  return (M * S.S_minimal).eval();
193  }
194 
195  /* [ABA] Y*S operator (Inertia Y,Constraint S) */
196  // inline Eigen::Matrix<double,6,3>
197  template<typename Matrix6Like, typename S2, int O2>
198 #if EIGEN_VERSION_AT_LEAST(3,2,90)
199  const typename Eigen::Product<
200  typename Eigen::internal::remove_const<typename SizeDepType<3>::ColsReturn<Matrix6Like>::ConstType>::type,
202  >
203 #else
204  const typename Eigen::ProductReturnType<
205  typename Eigen::internal::remove_const<typename SizeDepType<3>::ColsReturn<Matrix6Like>::ConstType>::type,
207  >::Type
208 #endif
209  operator*(const Eigen::MatrixBase<Matrix6Like> & Y,
211  {
212  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix6Like,6,6);
213  return Y.derived().template middleCols<3>(Inertia::ANGULAR) * S.S_minimal;
214  }
215 
216  namespace internal
217  {
218  template<typename S1, int O1>
219  struct SE3GroupAction< ConstraintSphericalZYXTpl<S1,O1> >
220  {
221 // typedef const typename Eigen::ProductReturnType<
222 // Eigen::Matrix <double,6,3,0>,
223 // Eigen::Matrix <double,3,3,0>
224 // >::Type Type;
225  typedef Eigen::Matrix<S1,6,3,O1> ReturnType;
226  };
227 
228  template<typename S1, int O1, typename MotionDerived>
229  struct MotionAlgebraAction< ConstraintSphericalZYXTpl<S1,O1>, MotionDerived >
230  {
231  typedef Eigen::Matrix<S1,6,3,O1> ReturnType;
232  };
233  }
234 
235  template<typename Scalar, int Options> struct JointSphericalZYXTpl;
236 
237  template<typename _Scalar, int _Options>
238  struct traits< JointSphericalZYXTpl<_Scalar,_Options> >
239  {
240  enum {
241  NQ = 3,
242  NV = 3
243  };
244  typedef _Scalar Scalar;
245  enum { Options = _Options };
252  typedef Eigen::Matrix<Scalar,6,NV,Options> F_t;
253 
254  // [ABA]
255  typedef Eigen::Matrix<Scalar,6,NV,Options> U_t;
256  typedef Eigen::Matrix<Scalar,NV,NV,Options> D_t;
257  typedef Eigen::Matrix<Scalar,6,NV,Options> UD_t;
258 
259  PINOCCHIO_JOINT_DATA_BASE_ACCESSOR_DEFAULT_RETURN_TYPE
260 
261  typedef Eigen::Matrix<Scalar,NQ,1,Options> ConfigVector_t;
262  typedef Eigen::Matrix<Scalar,NV,1,Options> TangentVector_t;
263  };
264 
265  template<typename Scalar, int Options>
266  struct traits< JointDataSphericalZYXTpl<Scalar,Options> >
268 
269  template<typename Scalar, int Options>
270  struct traits< JointModelSphericalZYXTpl<Scalar,Options> >
272 
273  template<typename _Scalar, int _Options>
274  struct JointDataSphericalZYXTpl : public JointDataBase< JointDataSphericalZYXTpl<_Scalar,_Options> >
275  {
276  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
278  PINOCCHIO_JOINT_DATA_TYPEDEF_TEMPLATE;
279  PINOCCHIO_JOINT_DATA_BASE_DEFAULT_ACCESSOR
280 
281  Constraint_t S;
282  Transformation_t M;
283  Motion_t v;
284  Bias_t c;
285 
286  F_t F;
287 
288  // [ABA] specific data
289  U_t U;
290  D_t Dinv;
291  UD_t UDinv;
292  D_t StU;
293 
294  JointDataSphericalZYXTpl () : M(1), U(), Dinv(), UDinv() {}
295 
296  static std::string classname() { return std::string("JointDataSphericalZYX"); }
297  std::string shortname() const { return classname(); }
298 
299  }; // strcut JointDataSphericalZYXTpl
300 
301  PINOCCHIO_JOINT_CAST_TYPE_SPECIALIZATION(JointModelSphericalZYXTpl);
302  template<typename _Scalar, int _Options>
304  : public JointModelBase< JointModelSphericalZYXTpl<_Scalar,_Options> >
305  {
306  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
308  PINOCCHIO_JOINT_TYPEDEF_TEMPLATE;
309 
311  using Base::id;
312  using Base::idx_q;
313  using Base::idx_v;
314  using Base::setIndexes;
315 
316  JointDataDerived createData() const { return JointDataDerived(); }
317 
318  template<typename ConfigVector>
319  void calc(JointDataDerived & data,
320  const typename Eigen::MatrixBase<ConfigVector> & qs) const
321  {
322  typename ConfigVector::template ConstFixedSegmentReturnType<NQ>::Type & q = qs.template segment<NQ>(idx_q());
323 
324  typedef typename ConfigVector::Scalar S2;
325 
326  S2 c0,s0; SINCOS(q(0), &s0, &c0);
327  S2 c1,s1; SINCOS(q(1), &s1, &c1);
328  S2 c2,s2; SINCOS(q(2), &s2, &c2);
329 
330  data.M.rotation () << c0 * c1,
331  c0 * s1 * s2 - s0 * c2,
332  c0 * s1 * c2 + s0 * s2,
333  s0 * c1,
334  s0 * s1 * s2 + c0 * c2,
335  s0 * s1 * c2 - c0 * s2,
336  -s1,
337  c1 * s2,
338  c1 * c2;
339 
340  data.S.S_minimal
341  << -s1, Scalar(0), Scalar(1),
342  c1 * s2, c2, Scalar(0),
343  c1 * c2, -s2, Scalar(0);
344  }
345 
346  template<typename ConfigVector, typename TangentVector>
347  void calc(JointDataDerived & data,
348  const typename Eigen::MatrixBase<ConfigVector> & qs,
349  const typename Eigen::MatrixBase<TangentVector> & vs) const
350  {
351  typename ConfigVector::template ConstFixedSegmentReturnType<NQ>::Type & q = qs.template segment<NQ>(idx_q());
352 
353  typedef typename ConfigVector::Scalar S2;
354 
355  S2 c0,s0; SINCOS(q(0), &s0, &c0);
356  S2 c1,s1; SINCOS(q(1), &s1, &c1);
357  S2 c2,s2; SINCOS(q(2), &s2, &c2);
358 
359  data.M.rotation () << c0 * c1,
360  c0 * s1 * s2 - s0 * c2,
361  c0 * s1 * c2 + s0 * s2,
362  s0 * c1,
363  s0 * s1 * s2 + c0 * c2,
364  s0 * s1 * c2 - c0 * s2,
365  -s1,
366  c1 * s2,
367  c1 * c2;
368 
369  data.S.S_minimal
370  << -s1, Scalar(0), Scalar(1),
371  c1 * s2, c2, Scalar(0),
372  c1 * c2, -s2, Scalar(0);
373 
374  typename TangentVector::template ConstFixedSegmentReturnType<NV>::Type
375  & q_dot = vs.template segment<NV>(idx_v());
376 
377  data.v().noalias() = data.S.S_minimal * q_dot;
378 
379  data.c()(0) = -c1 * q_dot(0) * q_dot(1);
380  data.c()(1) = -s1 * s2 * q_dot(0) * q_dot(1) + c1 * c2 * q_dot(0) * q_dot(2) - s2 * q_dot(1) * q_dot(2);
381  data.c()(2) = -s1 * c2 * q_dot(0) * q_dot(1) - c1 * s2 * q_dot(0) * q_dot(2) - c2 * q_dot(1) * q_dot(2);
382  }
383 
384  template<typename Matrix6Like>
385  void calc_aba(JointDataDerived & data, const Eigen::MatrixBase<Matrix6Like> & I, const bool update_I) const
386  {
387  data.U.noalias() = I.template middleCols<3>(Motion::ANGULAR) * data.S.S_minimal;
388  data.StU.noalias() = data.S.S_minimal.transpose() * data.U.template middleRows<3>(Motion::ANGULAR);
389 
390  // compute inverse
391  data.Dinv.setIdentity();
392  data.StU.llt().solveInPlace(data.Dinv);
393 
394  data.UDinv.noalias() = data.U * data.Dinv;
395 
396  if (update_I)
397  PINOCCHIO_EIGEN_CONST_CAST(Matrix6Like,I) -= data.UDinv * data.U.transpose();
398  }
399 
400  Scalar finiteDifferenceIncrement() const
401  {
402  using math::sqrt;
403  return 2.*sqrt(sqrt(Eigen::NumTraits<Scalar>::epsilon()));
404  }
405 
406  static std::string classname() { return std::string("JointModelSphericalZYX"); }
407  std::string shortname() const { return classname(); }
408 
410  template<typename NewScalar>
412  {
414  ReturnType res;
415  res.setIndexes(id(),idx_q(),idx_v());
416  return res;
417  }
418 
419  }; // struct JointModelSphericalZYXTpl
420 
421 } // namespace pinocchio
422 
423 #include <boost/type_traits.hpp>
424 
425 namespace boost
426 {
427  template<typename Scalar, int Options>
428  struct has_nothrow_constructor< ::pinocchio::JointModelSphericalZYXTpl<Scalar,Options> >
429  : public integral_constant<bool,true> {};
430 
431  template<typename Scalar, int Options>
432  struct has_nothrow_copy< ::pinocchio::JointModelSphericalZYXTpl<Scalar,Options> >
433  : public integral_constant<bool,true> {};
434 
435  template<typename Scalar, int Options>
436  struct has_nothrow_constructor< ::pinocchio::JointDataSphericalZYXTpl<Scalar,Options> >
437  : public integral_constant<bool,true> {};
438 
439  template<typename Scalar, int Options>
440  struct has_nothrow_copy< ::pinocchio::JointDataSphericalZYXTpl<Scalar,Options> >
441  : public integral_constant<bool,true> {};
442 }
443 
444 #endif // ifndef __pinocchio_joint_spherical_ZYX_hpp__
int idx_q(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointIdxQVisitor to get the index in the full model configuration space...
int idx_v(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointIdxVVisitor to get the index in the full model tangent space corre...
JointDataTpl< Scalar, Options, JointCollectionTpl > createData(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through CreateData visitor to create a JointDataTpl.
ModelTpl< Scalar, Options, JointCollectionTpl >::TangentVectorType finiteDifferenceIncrement(const ModelTpl< Scalar, Options, JointCollectionTpl > &model)
Computes the finite difference increments for each degree of freedom according to the current joint c...
void SINCOS(const Scalar &a, Scalar *sa, Scalar *ca)
Computes sin/cos values of a given input scalar.
Definition: sincos.hpp:28
std::string shortname(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointShortnameVisitor to get the shortname of the derived joint model...
void cross(const Eigen::MatrixBase< Vector3 > &v, const Eigen::MatrixBase< Matrix3xIn > &Min, const Eigen::MatrixBase< Matrix3xOut > &Mout)
Applies the cross product onto the columns of M.
Definition: skew.hpp:211
ConstAngularType angular() const
Return the angular part of the force vector.
Definition: force-base.hpp:40
Main pinocchio namespace.
Definition: treeview.dox:24
JointModelSphericalZYXTpl< NewScalar, Options > cast() const
void alphaSkew(const Scalar alpha, const Eigen::MatrixBase< Vector3 > &v, const Eigen::MatrixBase< Matrix3 > &M)
Computes the skew representation of a given 3d vector multiplied by a given scalar. i.e. the antisymmetric matrix representation of the cross product operator ( )
Definition: skew.hpp:124
Common traits structure to fully define base classes for CRTP.
Definition: spatial/fwd.hpp:29
void calc_aba(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel, JointDataTpl< Scalar, Options, JointCollectionTpl > &jdata, const Eigen::MatrixBase< Matrix6Type > &I, const bool update_I)
Visit a JointModelTpl and the corresponding JointDataTpl through JointCalcAbaVisitor to...