pinocchio  2.4.4
A fast and flexible implementation of Rigid Body Dynamics algorithms and their analytical derivatives
joint-free-flyer.hpp
1 //
2 // Copyright (c) 2015-2019 CNRS INRIA
3 // Copyright (c) 2015-2016 Wandercraft, 86 rue de Paris 91400 Orsay, France.
4 //
5 
6 #ifndef __pinocchio_joint_free_flyer_hpp__
7 #define __pinocchio_joint_free_flyer_hpp__
8 
9 #include "pinocchio/macros.hpp"
10 #include "pinocchio/spatial/inertia.hpp"
11 #include "pinocchio/spatial/explog.hpp"
12 #include "pinocchio/multibody/joint/joint-base.hpp"
13 #include "pinocchio/multibody/constraint.hpp"
14 #include "pinocchio/math/fwd.hpp"
15 #include "pinocchio/math/quaternion.hpp"
16 
17 namespace pinocchio
18 {
19 
20  template<typename Scalar, int Options> struct ConstraintIdentityTpl;
21 
22  template<typename _Scalar, int _Options>
23  struct traits< ConstraintIdentityTpl<_Scalar,_Options> >
24  {
25  typedef _Scalar Scalar;
26  enum { Options = _Options };
27  typedef Eigen::Matrix<Scalar,6,6,Options> Matrix6;
28  enum {
29  LINEAR = 0,
30  ANGULAR = 3
31  };
33  typedef Eigen::Matrix<Scalar,6,1,Options> JointForce;
34  typedef Eigen::Matrix<Scalar,6,6,Options> DenseBase;
35  typedef typename Matrix6::IdentityReturnType ConstMatrixReturnType;
36  typedef typename Matrix6::IdentityReturnType MatrixReturnType;
37  }; // traits ConstraintRevolute
38 
39 
40  template<typename _Scalar, int _Options>
42  : ConstraintBase< ConstraintIdentityTpl<_Scalar,_Options> >
43  {
44  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
45  PINOCCHIO_CONSTRAINT_TYPEDEF_TPL(ConstraintIdentityTpl)
46 
47  enum { NV = 6 };
48 
49  template<typename Vector6Like>
50  JointMotion __mult__(const Eigen::MatrixBase<Vector6Like> & vj) const
51  {
52  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector6Like,6);
53  return JointMotion(vj);
54  }
55 
56  template<typename S1, int O1>
58  se3Action(const SE3Tpl<S1,O1> & m) const
59  {
60  return m.toActionMatrix();
61 
62  }
63 
64  template<typename S1, int O1>
66  se3ActionInverse(const SE3Tpl<S1,O1> & m) const
67  {
68  return m.toActionMatrixInverse();
69  }
70 
71  int nv_impl() const { return NV; }
72 
74  {
75  template<typename Derived>
77  operator*(const ForceDense<Derived> & phi)
78  { return phi.toVector(); }
79 
80  /* [CRBA] MatrixBase operator* (Constraint::Transpose S, ForceSet::Block) */
81  template<typename MatrixDerived>
82  typename PINOCCHIO_EIGEN_REF_CONST_TYPE(MatrixDerived)
83  operator*(const Eigen::MatrixBase<MatrixDerived> & F)
84  {
85  return F.derived();
86  }
87  };
88 
89  TransposeConst transpose() const { return TransposeConst(); }
90  MatrixReturnType matrix_impl() const { return DenseBase::Identity(); }
91 
92  template<typename MotionDerived>
93  typename MotionDerived::ActionMatrixType
94  motionAction(const MotionBase<MotionDerived> & v) const
95  { return v.toActionMatrix(); }
96 
97  bool isEqual(const ConstraintIdentityTpl &) const { return true; }
98 
99  }; // struct ConstraintIdentityTpl
100 
101  template<typename Scalar, int Options, typename Vector6Like>
104  const Eigen::MatrixBase<Vector6Like> & v)
105  {
106  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector6Like,6);
107 // typedef typename ConstraintIdentityTpl<Scalar,Options>::Motion Motion;
109  return Motion(v.derived());
110  }
111 
112  /* [CRBA] ForceSet operator* (Inertia Y,Constraint S) */
113  template<typename S1, int O1, typename S2, int O2>
114  inline typename InertiaTpl<S1,O1>::Matrix6
116  {
117  return Y.matrix();
118  }
119 
120  /* [ABA] Y*S operator*/
121  template<typename Matrix6Like, typename S2, int O2>
122  inline typename PINOCCHIO_EIGEN_REF_CONST_TYPE(Matrix6Like)
123  operator*(const Eigen::MatrixBase<Matrix6Like> & Y, const ConstraintIdentityTpl<S2,O2> &)
124  {
125  return Y.derived();
126  }
127 
128  template<typename S1, int O1>
130  { typedef typename SE3Tpl<S1,O1>::ActionMatrixType ReturnType; };
131 
132  template<typename S1, int O1, typename MotionDerived>
133  struct MotionAlgebraAction< ConstraintIdentityTpl<S1,O1>,MotionDerived >
134  { typedef typename SE3Tpl<S1,O1>::ActionMatrixType ReturnType; };
135 
136  template<typename Scalar, int Options> struct JointFreeFlyerTpl;
137 
138  template<typename _Scalar, int _Options>
139  struct traits< JointFreeFlyerTpl<_Scalar,_Options> >
140  {
141  enum {
142  NQ = 7,
143  NV = 6
144  };
145  typedef _Scalar Scalar;
146  enum { Options = _Options };
153 
154  // [ABA]
155  typedef Eigen::Matrix<Scalar,6,NV,Options> U_t;
156  typedef Eigen::Matrix<Scalar,NV,NV,Options> D_t;
157  typedef Eigen::Matrix<Scalar,6,NV,Options> UD_t;
158 
159  PINOCCHIO_JOINT_DATA_BASE_ACCESSOR_DEFAULT_RETURN_TYPE
160 
161  typedef Eigen::Matrix<Scalar,NQ,1,Options> ConfigVector_t;
162  typedef Eigen::Matrix<Scalar,NV,1,Options> TangentVector_t;
163  };
164 
165  template<typename Scalar, int Options>
166  struct traits< JointDataFreeFlyerTpl<Scalar,Options> >
168 
169  template<typename Scalar, int Options>
170  struct traits< JointModelFreeFlyerTpl<Scalar,Options> >
172 
173  template<typename _Scalar, int _Options>
174  struct JointDataFreeFlyerTpl : public JointDataBase< JointDataFreeFlyerTpl<_Scalar,_Options> >
175  {
176  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
178  PINOCCHIO_JOINT_DATA_TYPEDEF_TEMPLATE(JointDerived);
179  PINOCCHIO_JOINT_DATA_BASE_DEFAULT_ACCESSOR
180 
181  Constraint_t S;
182  Transformation_t M;
183  Motion_t v;
184  Bias_t c;
185 
186  // [ABA] specific data
187  U_t U;
188  D_t Dinv;
189  UD_t UDinv;
190 
192  : M(Transformation_t::Identity())
193  , v(Motion_t::Zero())
194  , U(U_t::Zero())
195  , Dinv(D_t::Zero())
196  , UDinv(UD_t::Identity())
197  {}
198 
199  static std::string classname() { return std::string("JointDataFreeFlyer"); }
200  std::string shortname() const { return classname(); }
201 
202  }; // struct JointDataFreeFlyerTpl
203 
204  PINOCCHIO_JOINT_CAST_TYPE_SPECIALIZATION(JointModelFreeFlyerTpl);
205  template<typename _Scalar, int _Options>
207  : public JointModelBase< JointModelFreeFlyerTpl<_Scalar,_Options> >
208  {
209  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
211  PINOCCHIO_JOINT_TYPEDEF_TEMPLATE(JointDerived);
212 
214  using Base::id;
215  using Base::idx_q;
216  using Base::idx_v;
217  using Base::setIndexes;
218 
219  JointDataDerived createData() const { return JointDataDerived(); }
220 
221  template<typename ConfigVectorLike>
222  inline void forwardKinematics(Transformation_t & M, const Eigen::MatrixBase<ConfigVectorLike> & q_joint) const
223  {
224  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(ConfigVector_t,ConfigVectorLike);
225  typedef typename Eigen::Quaternion<typename ConfigVectorLike::Scalar,PINOCCHIO_EIGEN_PLAIN_TYPE(ConfigVectorLike)::Options> Quaternion;
226  typedef Eigen::Map<const Quaternion> ConstQuaternionMap;
227 
228  ConstQuaternionMap quat(q_joint.template tail<4>().data());
229  //assert(math::fabs(quat.coeffs().squaredNorm()-1.) <= sqrt(Eigen::NumTraits<typename V::Scalar>::epsilon())); TODO: check validity of the rhs precision
230  assert(math::fabs(quat.coeffs().squaredNorm()-1.) <= 1e-4);
231 
232  M.rotation(quat.matrix());
233  M.translation(q_joint.template head<3>());
234  }
235 
236  template<typename Vector3Derived, typename QuaternionDerived>
237  EIGEN_DONT_INLINE
238  void calc(JointDataDerived & data,
239  const typename Eigen::MatrixBase<Vector3Derived> & trans,
240  const typename Eigen::QuaternionBase<QuaternionDerived> & quat) const
241  {
242  data.M.translation(trans);
243  data.M.rotation(quat.matrix());
244  }
245 
246  template<typename ConfigVector>
247  EIGEN_DONT_INLINE
248  void calc(JointDataDerived & data,
249  const typename Eigen::PlainObjectBase<ConfigVector> & qs) const
250  {
251  typedef typename Eigen::Quaternion<typename ConfigVector::Scalar,ConfigVector::Options> Quaternion;
252  typedef Eigen::Map<const Quaternion> ConstQuaternionMap;
253 
254  typename ConfigVector::template ConstFixedSegmentReturnType<NQ>::Type q = qs.template segment<NQ>(idx_q());
255  ConstQuaternionMap quat(q.template tail<4>().data());
256 
257  calc(data,q.template head<3>(),quat);
258  }
259 
260  template<typename ConfigVector>
261  EIGEN_DONT_INLINE
262  void calc(JointDataDerived & data,
263  const typename Eigen::MatrixBase<ConfigVector> & qs) const
264  {
265  typedef typename Eigen::Quaternion<Scalar,Options> Quaternion;
266 
267  typename ConfigVector::template ConstFixedSegmentReturnType<NQ>::Type q = qs.template segment<NQ>(idx_q());
268  const Quaternion quat(q.template tail<4>());
269 
270  calc(data,q.template head<3>(),quat);
271  }
272 
273  template<typename ConfigVector, typename TangentVector>
274  EIGEN_DONT_INLINE
275  void calc(JointDataDerived & data,
276  const typename Eigen::MatrixBase<ConfigVector> & qs,
277  const typename Eigen::MatrixBase<TangentVector> & vs) const
278  {
279  calc(data,qs.derived());
280 
281  data.v = vs.template segment<NV>(idx_v());
282  }
283 
284  template<typename Matrix6Like>
285  void calc_aba(JointDataDerived & data,
286  const Eigen::MatrixBase<Matrix6Like> & I,
287  const bool update_I) const
288  {
289  data.U = I;
290 
291  // compute inverse
292 // data.Dinv.setIdentity();
293 // I.llt().solveInPlace(data.Dinv);
294  internal::PerformStYSInversion<Scalar>::run(I,data.Dinv);
295 
296  if (update_I)
297  PINOCCHIO_EIGEN_CONST_CAST(Matrix6Like,I).setZero();
298  }
299 
300  static std::string classname() { return std::string("JointModelFreeFlyer"); }
301  std::string shortname() const { return classname(); }
302 
304  template<typename NewScalar>
306  {
308  ReturnType res;
309  res.setIndexes(id(),idx_q(),idx_v());
310  return res;
311  }
312 
313  }; // struct JointModelFreeFlyerTpl
314 
315 } // namespace pinocchio
316 
317 #include <boost/type_traits.hpp>
318 
319 namespace boost
320 {
321  template<typename Scalar, int Options>
322  struct has_nothrow_constructor< ::pinocchio::JointModelFreeFlyerTpl<Scalar,Options> >
323  : public integral_constant<bool,true> {};
324 
325  template<typename Scalar, int Options>
326  struct has_nothrow_copy< ::pinocchio::JointModelFreeFlyerTpl<Scalar,Options> >
327  : public integral_constant<bool,true> {};
328 
329  template<typename Scalar, int Options>
330  struct has_nothrow_constructor< ::pinocchio::JointDataFreeFlyerTpl<Scalar,Options> >
331  : public integral_constant<bool,true> {};
332 
333  template<typename Scalar, int Options>
334  struct has_nothrow_copy< ::pinocchio::JointDataFreeFlyerTpl<Scalar,Options> >
335  : public integral_constant<bool,true> {};
336 }
337 
338 #endif // ifndef __pinocchio_joint_free_flyer_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...
Definition: casadi.hpp:13
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...
Return type of the ation of a Motion onto an object of type D.
Definition: motion.hpp:44
JointDataTpl< Scalar, Options, JointCollectionTpl > createData(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through CreateData visitor to create a JointDataTpl.
JointModelFreeFlyerTpl< NewScalar, Options > cast() const
ToVectorConstReturnType toVector() const
Return the force as an Eigen vector.
Definition: force-base.hpp:81
ActionMatrixType toActionMatrixInverse() const
The action matrix of .
Definition: se3-base.hpp:70
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.
std::string shortname(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointShortnameVisitor to get the shortname of the derived joint model...
ActionMatrixType toActionMatrix() const
The action matrix of .
Definition: se3-base.hpp:60
Main pinocchio namespace.
Definition: treeview.dox:24
Common traits structure to fully define base classes for CRTP.
Definition: fwd.hpp:35
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...
MultiplicationOp< InertiaTpl< Scalar, Options >, ConstraintDerived >::ReturnType operator*(const InertiaTpl< Scalar, Options > &Y, const ConstraintBase< ConstraintDerived > &constraint)
 .