special-euclidean.hpp
Go to the documentation of this file.
1 //
2 // Copyright (c) 2016-2021 CNRS INRIA
3 //
4 
5 #ifndef __pinocchio_multibody_liegroup_special_euclidean_operation_hpp__
6 #define __pinocchio_multibody_liegroup_special_euclidean_operation_hpp__
7 
8 #include <limits>
9 
10 #include "pinocchio/macros.hpp"
17 
21 
22 namespace pinocchio
23 {
24  template<int Dim, typename Scalar, int Options = 0>
26  {
27  };
28 
29  template<int Dim, typename Scalar, int Options>
31  {
32  };
33 
34  template<typename _Scalar, int _Options>
35  struct traits<SpecialEuclideanOperationTpl<2, _Scalar, _Options>>
36  {
37  typedef _Scalar Scalar;
38  enum
39  {
40  Options = _Options,
41  NQ = 4,
42  NV = 3
43  };
44  };
45 
46  // SE(2)
47  template<typename _Scalar, int _Options>
48  struct SpecialEuclideanOperationTpl<2, _Scalar, _Options>
49  : public LieGroupBase<SpecialEuclideanOperationTpl<2, _Scalar, _Options>>
50  {
52 
56 
57  typedef Eigen::Matrix<Scalar, 2, 2, Options> Matrix2;
58  typedef Eigen::Matrix<Scalar, 2, 1, Options> Vector2;
59 
60  template<typename TangentVector, typename Matrix2Like, typename Vector2Like>
61  static void exp(
62  const Eigen::MatrixBase<TangentVector> & v,
63  const Eigen::MatrixBase<Matrix2Like> & R,
64  const Eigen::MatrixBase<Vector2Like> & t)
65  {
66  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(TangentVector_t, TangentVector);
67  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix2Like, 2, 2);
68  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector2Like, 2);
69 
70  typedef typename Matrix2Like::Scalar Scalar;
71  const Scalar omega = v(2);
72  Scalar cv, sv;
73  SINCOS(omega, &sv, &cv);
74  PINOCCHIO_EIGEN_CONST_CAST(Matrix2Like, R) << cv, -sv, sv, cv;
76 
77  {
78  typename PINOCCHIO_EIGEN_PLAIN_TYPE(Vector2Like) vcross(-v(1), v(0));
79  vcross -= -v(1) * R.col(0) + v(0) * R.col(1);
80  vcross /= omega;
81  Scalar omega_abs = math::fabs(omega);
82  PINOCCHIO_EIGEN_CONST_CAST(Vector2Like, t).coeffRef(0) = if_then_else(
83  internal::GT, omega_abs, Scalar(1e-14), // TODO: change hard coded value
84  vcross.coeff(0), v.coeff(0));
85 
86  PINOCCHIO_EIGEN_CONST_CAST(Vector2Like, t).coeffRef(1) = if_then_else(
87  internal::GT, omega_abs, Scalar(1e-14), // TODO: change hard coded value
88  vcross.coeff(1), v.coeff(1));
89  }
90  }
91 
92  template<typename Matrix2Like, typename Vector2Like, typename Matrix3Like>
93  static void toInverseActionMatrix(
94  const Eigen::MatrixBase<Matrix2Like> & R,
95  const Eigen::MatrixBase<Vector2Like> & t,
96  const Eigen::MatrixBase<Matrix3Like> & M,
97  const AssignmentOperatorType op)
98  {
99  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix2Like, 2, 2);
100  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector2Like, 2);
101  PINOCCHIO_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix3Like, M, 3, 3);
102  Matrix3Like & Mout = PINOCCHIO_EIGEN_CONST_CAST(Matrix3Like, M);
103  typedef typename Matrix3Like::Scalar Scalar;
104 
105  typename PINOCCHIO_EIGEN_PLAIN_TYPE(Vector2Like) tinv((R.transpose() * t).reverse());
106  tinv[0] *= Scalar(-1.);
107  switch (op)
108  {
109  case SETTO:
110  Mout.template topLeftCorner<2, 2>() = R.transpose();
111  Mout.template topRightCorner<2, 1>() = tinv;
112  Mout.template bottomLeftCorner<1, 2>().setZero();
113  Mout(2, 2) = (Scalar)1;
114  break;
115  case ADDTO:
116  Mout.template topLeftCorner<2, 2>() += R.transpose();
117  Mout.template topRightCorner<2, 1>() += tinv;
118  Mout(2, 2) += (Scalar)1;
119  break;
120  case RMTO:
121  Mout.template topLeftCorner<2, 2>() -= R.transpose();
122  Mout.template topRightCorner<2, 1>() -= tinv;
123  Mout(2, 2) -= (Scalar)1;
124  break;
125  default:
126  assert(false && "Wrong Op requesed value");
127  break;
128  }
129  }
130 
131  template<typename Matrix2Like, typename Vector2Like, typename TangentVector>
132  static void log(
133  const Eigen::MatrixBase<Matrix2Like> & R,
134  const Eigen::MatrixBase<Vector2Like> & p,
135  const Eigen::MatrixBase<TangentVector> & v)
136  {
137  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix2Like, 2, 2);
138  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector2Like, 2);
139  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(TangentVector_t, TangentVector);
140 
141  TangentVector & vout = PINOCCHIO_EIGEN_CONST_CAST(TangentVector, v);
142 
143  typedef typename Matrix2Like::Scalar Scalar1;
144 
145  Scalar1 t = SO2_t::log(R);
146  const Scalar1 tabs = math::fabs(t);
147  const Scalar1 t2 = t * t;
148  Scalar1 st, ct;
149  SINCOS(tabs, &st, &ct);
150  Scalar1 alpha;
152  internal::LT, tabs, Scalar(1e-4), // TODO: change hard coded value
153  static_cast<Scalar>(1 - t2 / 12 - t2 * t2 / 720),
154  static_cast<Scalar>(tabs * st / (2 * (1 - ct))));
155 
156  vout.template head<2>().noalias() = alpha * p;
157  vout(0) += t / 2 * p(1);
158  vout(1) += -t / 2 * p(0);
159  vout(2) = t;
160  }
161 
162  template<typename Matrix2Like, typename Vector2Like, typename JacobianOutLike>
163  static void Jlog(
164  const Eigen::MatrixBase<Matrix2Like> & R,
165  const Eigen::MatrixBase<Vector2Like> & p,
166  const Eigen::MatrixBase<JacobianOutLike> & J)
167  {
168  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix2Like, 2, 2);
169  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector2Like, 2);
170  EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(JacobianOutLike, JacobianMatrix_t);
171 
172  JacobianOutLike & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOutLike, J);
173 
174  typedef typename Matrix2Like::Scalar Scalar1;
175 
176  Scalar1 t = SO2_t::log(R);
177  const Scalar1 tabs = math::fabs(t);
178  Scalar1 alpha, alpha_dot;
179  Scalar1 t2 = t * t;
180  Scalar1 st, ct;
181  SINCOS(t, &st, &ct);
182  Scalar1 inv_2_1_ct = 0.5 / (1 - ct);
183 
185  internal::LT, tabs, Scalar(1e-4), static_cast<Scalar>(1 - t2 / 12),
186  static_cast<Scalar>(t * st * inv_2_1_ct));
187  alpha_dot = internal::if_then_else(
188  internal::LT, tabs, Scalar(1e-4), static_cast<Scalar>(-t / 6 - t2 * t / 180),
189  static_cast<Scalar>((st - t) * inv_2_1_ct));
190 
191  typename PINOCCHIO_EIGEN_PLAIN_TYPE(Matrix2Like) V;
192  V(0, 0) = V(1, 1) = alpha;
193  V(1, 0) = -t / 2;
194  V(0, 1) = -V(1, 0);
195 
196  Jout.template topLeftCorner<2, 2>().noalias() = V * R;
197  Jout.template topRightCorner<2, 1>() << alpha_dot * p[0] + p[1] / 2,
198  -p(0) / 2 + alpha_dot * p(1);
199  Jout.template bottomLeftCorner<1, 2>().setZero();
200  Jout(2, 2) = 1;
201  }
202 
207  static Index nq()
208  {
209  return NQ;
210  }
212  static Index nv()
213  {
214  return NV;
215  }
216 
218  {
220  n << Scalar(0), Scalar(0), Scalar(1), Scalar(0);
221  return n;
222  }
223 
224  static std::string name()
225  {
226  return std::string("SE(2)");
227  }
228 
229  template<class ConfigL_t, class ConfigR_t, class Tangent_t>
230  static void difference_impl(
231  const Eigen::MatrixBase<ConfigL_t> & q0,
232  const Eigen::MatrixBase<ConfigR_t> & q1,
233  const Eigen::MatrixBase<Tangent_t> & d)
234  {
237  Matrix2 R0, R1;
238  Vector2 t0, t1;
239  forwardKinematics(R0, t0, q0);
240  forwardKinematics(R1, t1, q1);
242  Matrix2 R(R0.transpose() * R1);
243  Vector2 t(R0.transpose() * (t1 - t0));
244 
245  log(R, t, d);
246  }
247 
248  template<ArgumentPosition arg, class ConfigL_t, class ConfigR_t, class JacobianOut_t>
250  const Eigen::MatrixBase<ConfigL_t> & q0,
251  const Eigen::MatrixBase<ConfigR_t> & q1,
252  const Eigen::MatrixBase<JacobianOut_t> & J) const
253  {
256  Matrix2 R0, R1;
257  Vector2 t0, t1;
258  forwardKinematics(R0, t0, q0);
259  forwardKinematics(R1, t1, q1);
261  Matrix2 R(R0.transpose() * R1);
262  Vector2 t(R0.transpose() * (t1 - t0));
263 
264  if (arg == ARG0)
265  {
268  JacobianMatrix_t J1;
269  Jlog(R, t, J1);
271 
272  // pcross = [ y1-y0, - (x1 - x0) ]
273  Vector2 pcross(q1(1) - q0(1), q0(0) - q1(0));
274 
275  JacobianOut_t & J0 = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J);
276  J0.template topLeftCorner<2, 2>().noalias() = -R.transpose();
277  J0.template topRightCorner<2, 1>().noalias() = R1.transpose() * pcross;
278  J0.template bottomLeftCorner<1, 2>().setZero();
279  J0(2, 2) = -1;
280  J0.applyOnTheLeft(J1);
281  }
282  else if (arg == ARG1)
283  {
284  Jlog(R, t, J);
285  }
286  }
287 
288  template<class ConfigIn_t, class Velocity_t, class ConfigOut_t>
289  static void integrate_impl(
290  const Eigen::MatrixBase<ConfigIn_t> & q,
291  const Eigen::MatrixBase<Velocity_t> & v,
292  const Eigen::MatrixBase<ConfigOut_t> & qout)
293  {
294  ConfigOut_t & out = PINOCCHIO_EIGEN_CONST_CAST(ConfigOut_t, qout);
295 
298  Matrix2 R0, R;
299  Vector2 t0, t;
300  forwardKinematics(R0, t0, q);
301  exp(v, R, t);
303 
304  out.template head<2>().noalias() = R0 * t + t0;
305  out.template tail<2>().noalias() = R0 * R.col(0);
306  }
307 
308  template<class Config_t, class Jacobian_t>
310  const Eigen::MatrixBase<Config_t> & q, const Eigen::MatrixBase<Jacobian_t> & J)
311  {
312  assert(J.rows() == nq() && J.cols() == nv() && "J is not of the right dimension");
313 
314  Jacobian_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J);
315  Jout.setZero();
316 
317  const typename Config_t::Scalar &c_theta = q(2), &s_theta = q(3);
318 
319  Jout.template topLeftCorner<2, 2>() << c_theta, -s_theta, s_theta, c_theta;
320  Jout.template bottomRightCorner<2, 1>() << -s_theta, c_theta;
321  }
322 
323  template<class Config_t, class Tangent_t, class JacobianOut_t>
324  static void dIntegrate_dq_impl(
325  const Eigen::MatrixBase<Config_t> & /*q*/,
326  const Eigen::MatrixBase<Tangent_t> & v,
327  const Eigen::MatrixBase<JacobianOut_t> & J,
328  const AssignmentOperatorType op = SETTO)
329  {
330  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J);
331 
334  Matrix2 R;
335  Vector2 t;
336  exp(v, R, t);
338 
339  toInverseActionMatrix(R, t, Jout, op);
340  }
341 
342  template<class Config_t, class Tangent_t, class JacobianOut_t>
343  static void dIntegrate_dv_impl(
344  const Eigen::MatrixBase<Config_t> & /*q*/,
345  const Eigen::MatrixBase<Tangent_t> & v,
346  const Eigen::MatrixBase<JacobianOut_t> & J,
347  const AssignmentOperatorType op = SETTO)
348  {
349  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J);
350  // TODO sparse version
352  nu.toVector() << v.template head<2>(), 0, 0, 0, v[2];
355  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
356  Jexp6(nu, Jtmp6);
358 
359  switch (op)
360  {
361  case SETTO:
362  Jout << Jtmp6.template topLeftCorner<2, 2>(), Jtmp6.template topRightCorner<2, 1>(),
363  Jtmp6.template bottomLeftCorner<1, 2>(), Jtmp6.template bottomRightCorner<1, 1>();
364  break;
365  case ADDTO:
366  Jout.template topLeftCorner<2, 2>() += Jtmp6.template topLeftCorner<2, 2>();
367  Jout.template topRightCorner<2, 1>() += Jtmp6.template topRightCorner<2, 1>();
368  Jout.template bottomLeftCorner<1, 2>() += Jtmp6.template bottomLeftCorner<1, 2>();
369  Jout.template bottomRightCorner<1, 1>() += Jtmp6.template bottomRightCorner<1, 1>();
370  break;
371  case RMTO:
372  Jout.template topLeftCorner<2, 2>() -= Jtmp6.template topLeftCorner<2, 2>();
373  Jout.template topRightCorner<2, 1>() -= Jtmp6.template topRightCorner<2, 1>();
374  Jout.template bottomLeftCorner<1, 2>() -= Jtmp6.template bottomLeftCorner<1, 2>();
375  Jout.template bottomRightCorner<1, 1>() -= Jtmp6.template bottomRightCorner<1, 1>();
376  break;
377  default:
378  assert(false && "Wrong Op requesed value");
379  break;
380  }
381  }
382 
383  template<class Config_t, class Tangent_t, class JacobianIn_t, class JacobianOut_t>
385  const Eigen::MatrixBase<Config_t> & /*q*/,
386  const Eigen::MatrixBase<Tangent_t> & v,
387  const Eigen::MatrixBase<JacobianIn_t> & Jin,
388  const Eigen::MatrixBase<JacobianOut_t> & J_out)
389  {
390  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J_out);
393  Matrix2 R;
394  Vector2 t;
395  exp(v, R, t);
397 
398  Vector2 tinv = (R.transpose() * t).reverse();
399  tinv[0] *= Scalar(-1.);
400 
401  Jout.template topRows<2>().noalias() = R.transpose() * Jin.template topRows<2>();
402  Jout.template topRows<2>().noalias() += tinv * Jin.template bottomRows<1>();
403  Jout.template bottomRows<1>() = Jin.template bottomRows<1>();
404  }
405 
406  template<class Config_t, class Tangent_t, class JacobianIn_t, class JacobianOut_t>
408  const Eigen::MatrixBase<Config_t> & /*q*/,
409  const Eigen::MatrixBase<Tangent_t> & v,
410  const Eigen::MatrixBase<JacobianIn_t> & Jin,
411  const Eigen::MatrixBase<JacobianOut_t> & J_out)
412  {
413  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J_out);
415  nu.toVector() << v.template head<2>(), 0, 0, 0, v[2];
416 
419  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
420  Jexp6(nu, Jtmp6);
422 
423  Jout.template topRows<2>().noalias() =
424  Jtmp6.template topLeftCorner<2, 2>() * Jin.template topRows<2>();
425  Jout.template topRows<2>().noalias() +=
426  Jtmp6.template topRightCorner<2, 1>() * Jin.template bottomRows<1>();
427  Jout.template bottomRows<1>().noalias() =
428  Jtmp6.template bottomLeftCorner<1, 2>() * Jin.template topRows<2>();
429  Jout.template bottomRows<1>().noalias() +=
430  Jtmp6.template bottomRightCorner<1, 1>() * Jin.template bottomRows<1>();
431  }
432 
433  template<class Config_t, class Tangent_t, class Jacobian_t>
435  const Eigen::MatrixBase<Config_t> & /*q*/,
436  const Eigen::MatrixBase<Tangent_t> & v,
437  const Eigen::MatrixBase<Jacobian_t> & J)
438  {
439  Jacobian_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J);
442  Matrix2 R;
443  Vector2 t;
444  exp(v, R, t);
446 
447  Vector2 tinv = (R.transpose() * t).reverse();
448  tinv[0] *= Scalar(-1);
449  // TODO: Aliasing
450  Jout.template topRows<2>() = R.transpose() * Jout.template topRows<2>();
451  // No Aliasing
452  Jout.template topRows<2>().noalias() += tinv * Jout.template bottomRows<1>();
453  }
454 
455  template<class Config_t, class Tangent_t, class Jacobian_t>
457  const Eigen::MatrixBase<Config_t> & /*q*/,
458  const Eigen::MatrixBase<Tangent_t> & v,
459  const Eigen::MatrixBase<Jacobian_t> & J)
460  {
461  Jacobian_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J);
463  nu.toVector() << v.template head<2>(), 0, 0, 0, v[2];
464 
467  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
468  Jexp6(nu, Jtmp6);
470  // TODO: Remove aliasing
471  Jout.template topRows<2>() =
472  Jtmp6.template topLeftCorner<2, 2>() * Jout.template topRows<2>();
473  Jout.template topRows<2>().noalias() +=
474  Jtmp6.template topRightCorner<2, 1>() * Jout.template bottomRows<1>();
475  Jout.template bottomRows<1>() =
476  Jtmp6.template bottomRightCorner<1, 1>() * Jout.template bottomRows<1>();
477  Jout.template bottomRows<1>().noalias() +=
478  Jtmp6.template bottomLeftCorner<1, 2>() * Jout.template topRows<2>();
479  }
480 
481  template<class Config_t>
482  static void normalize_impl(const Eigen::MatrixBase<Config_t> & qout)
483  {
484  pinocchio::normalize(qout.const_cast_derived().template tail<2>());
485  }
486 
487  template<class Config_t>
488  static bool isNormalized_impl(const Eigen::MatrixBase<Config_t> & qin, const Scalar & prec)
489  {
490  const Scalar norm = Scalar(qin.template tail<2>().norm());
491  using std::abs;
492  return abs(norm - Scalar(1.0)) < prec;
493  }
494 
495  template<class Config_t>
496  static void random_impl(const Eigen::MatrixBase<Config_t> & qout)
497  {
498  R2crossSO2_t().random(qout);
499  }
500 
501  template<class ConfigL_t, class ConfigR_t, class ConfigOut_t>
503  const Eigen::MatrixBase<ConfigL_t> & lower,
504  const Eigen::MatrixBase<ConfigR_t> & upper,
505  const Eigen::MatrixBase<ConfigOut_t> & qout)
506  {
507  R2crossSO2_t().randomConfiguration(lower, upper, qout);
508  }
509 
510  template<class ConfigL_t, class ConfigR_t>
512  const Eigen::MatrixBase<ConfigL_t> & q0,
513  const Eigen::MatrixBase<ConfigR_t> & q1,
514  const Scalar & prec)
515  {
516  return R2crossSO2_t().isSameConfiguration(q0, q1, prec);
517  }
518 
519  protected:
520  template<typename Matrix2Like, typename Vector2Like, typename Vector4Like>
521  static void forwardKinematics(
522  const Eigen::MatrixBase<Matrix2Like> & R,
523  const Eigen::MatrixBase<Vector2Like> & t,
524  const Eigen::MatrixBase<Vector4Like> & q)
525  {
526  EIGEN_STATIC_ASSERT_SAME_MATRIX_SIZE(Matrix2Like, Matrix2);
527  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(Vector2Like, Vector2);
528  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(ConfigVector_t, Vector4Like);
529 
530  PINOCCHIO_EIGEN_CONST_CAST(Vector2Like, t) = q.template head<2>();
531  const typename Vector4Like::Scalar &c_theta = q(2), &s_theta = q(3);
532  PINOCCHIO_EIGEN_CONST_CAST(Matrix2Like, R) << c_theta, -s_theta, s_theta, c_theta;
533  }
534  }; // struct SpecialEuclideanOperationTpl<2>
535 
536  template<typename _Scalar, int _Options>
537  struct traits<SpecialEuclideanOperationTpl<3, _Scalar, _Options>>
538  {
539  typedef _Scalar Scalar;
540  enum
541  {
542  Options = _Options,
543  NQ = 7,
544  NV = 6
545  };
546  };
547 
549  template<typename _Scalar, int _Options>
550  struct SpecialEuclideanOperationTpl<3, _Scalar, _Options>
551  : public LieGroupBase<SpecialEuclideanOperationTpl<3, _Scalar, _Options>>
552  {
554 
559 
560  typedef Eigen::Quaternion<Scalar, Options> Quaternion_t;
561  typedef Eigen::Map<Quaternion_t> QuaternionMap_t;
562  typedef Eigen::Map<const Quaternion_t> ConstQuaternionMap_t;
565  typedef typename Eigen::NumTraits<Scalar>::Real RealScalar;
566 
571  static Index nq()
572  {
573  return NQ;
574  }
576  static Index nv()
577  {
578  return NV;
579  }
580 
582  {
584  n.template head<6>().setZero();
585  n[6] = 1;
586  return n;
587  }
588 
589  static std::string name()
590  {
591  return std::string("SE(3)");
592  }
593 
594  template<class ConfigL_t, class ConfigR_t, class Tangent_t>
595  static void difference_impl(
596  const Eigen::MatrixBase<ConfigL_t> & q0,
597  const Eigen::MatrixBase<ConfigR_t> & q1,
598  const Eigen::MatrixBase<Tangent_t> & d)
599  {
600  typedef typename Tangent_t::Scalar Scalar;
601  ConstQuaternionMap_t quat0(q0.derived().template tail<4>().data());
604  ConstQuaternionMap_t quat1(q1.derived().template tail<4>().data());
607 
608  typedef Eigen::Matrix<Scalar, 3, 1, PINOCCHIO_EIGEN_PLAIN_TYPE(Tangent_t)::Options> Vector3;
609  const Vector3 dv_pre = q1.derived().template head<3>() - q0.derived().template head<3>();
610  const Vector3 dv = quat0.conjugate() * dv_pre;
611  PINOCCHIO_EIGEN_CONST_CAST(Tangent_t, d).noalias() =
612  log6(quat0.conjugate() * quat1, dv).toVector();
613  }
614 
617  template<ArgumentPosition arg, class ConfigL_t, class ConfigR_t, class JacobianOut_t>
618  static void dDifference_impl(
619  const Eigen::MatrixBase<ConfigL_t> & q0,
620  const Eigen::MatrixBase<ConfigR_t> & q1,
621  const Eigen::MatrixBase<JacobianOut_t> & J)
622  {
623  typedef typename SE3::Vector3 Vector3;
624 
625  ConstQuaternionMap_t quat0(q0.derived().template tail<4>().data());
628  ConstQuaternionMap_t quat1(q1.derived().template tail<4>().data());
631 
632  const Vector3 dv_pre = q1.derived().template head<3>() - q0.derived().template head<3>();
633  const Vector3 trans = quat0.conjugate() * dv_pre;
634 
635  const Quaternion_t quat_diff = quat0.conjugate() * quat1;
636 
637  const SE3 M(quat_diff, trans);
638 
639  if (arg == ARG0)
640  {
643  JacobianMatrix_t J1;
644  Jlog6(M, J1);
646 
647  const Vector3 p1_p0 = quat1.conjugate() * dv_pre;
648 
649  JacobianOut_t & J0 = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J);
650  J0.template bottomRightCorner<3, 3>() = J0.template topLeftCorner<3, 3>() =
651  -M.rotation().transpose();
652  J0.template topRightCorner<3, 3>().noalias() =
653  skew(p1_p0) * M.rotation().transpose(); // = R1.T * skew(q1_t - q0_t) * R0;
654  J0.template bottomLeftCorner<3, 3>().setZero();
655  J0.applyOnTheLeft(J1);
656  }
657  else if (arg == ARG1)
658  {
659  Jlog6(M, J);
660  }
661  }
662 
663  template<class ConfigIn_t, class Velocity_t, class ConfigOut_t>
664  static void integrate_impl(
665  const Eigen::MatrixBase<ConfigIn_t> & q,
666  const Eigen::MatrixBase<Velocity_t> & v,
667  const Eigen::MatrixBase<ConfigOut_t> & qout)
668  {
669  ConfigOut_t & out = PINOCCHIO_EIGEN_CONST_CAST(ConfigOut_t, qout);
670  Quaternion_t const quat(q.derived().template tail<4>());
673  QuaternionMap_t res_quat(out.template tail<4>().data());
674 
676 
677  typedef typename ConfigOut_t::Scalar Scalar;
678  enum
679  {
681  };
682 
683  Eigen::Matrix<Scalar, 7, 1, Options> expv;
684  quaternion::exp6(v, expv);
685 
686  out.template head<3>() = (quat * expv.template head<3>()) + q.derived().template head<3>();
687 
688  ConstQuaternionMap_t quat1(expv.template tail<4>().data());
689  res_quat = quat * quat1;
690 
691  const Scalar dot_product = res_quat.dot(quat);
692  for (Eigen::DenseIndex k = 0; k < 4; ++k)
693  {
694  res_quat.coeffs().coeffRef(k) = if_then_else(
695  internal::LT, dot_product, Scalar(0), static_cast<Scalar>(-res_quat.coeffs().coeff(k)),
696  res_quat.coeffs().coeff(k));
697  }
698 
699  // Norm of qs might be epsilon-different to 1, so M1.rotation might be epsilon-different to a
700  // rotation matrix. It is then safer to re-normalized after converting M1.rotation to
701  // quaternion.
705  }
706 
707  template<class Config_t, class Jacobian_t>
709  const Eigen::MatrixBase<Config_t> & q, const Eigen::MatrixBase<Jacobian_t> & J)
710  {
711  assert(J.rows() == nq() && J.cols() == nv() && "J is not of the right dimension");
712 
713  typedef typename PINOCCHIO_EIGEN_PLAIN_TYPE(Config_t) ConfigPlainType;
714  typedef typename PINOCCHIO_EIGEN_PLAIN_TYPE(Jacobian_t) JacobianPlainType;
715  typedef typename ConfigPlainType::Scalar Scalar;
716 
717  Jacobian_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J);
718  Jout.setZero();
719 
720  ConstQuaternionMap_t quat_map(q.derived().template tail<4>().data());
723  Jout.template topLeftCorner<3, 3>() = quat_map.toRotationMatrix();
724  // Jexp3(quat,Jout.template bottomRightCorner<4,3>());
725 
726  typedef Eigen::Matrix<Scalar, 4, 3, JacobianPlainType::Options | Eigen::RowMajor> Jacobian43;
730  Jacobian43 Jexp3QuatCoeffWise;
731 
732  Scalar theta;
733  typename SE3::Vector3 v = quaternion::log3(quat_map, theta);
734  quaternion::Jexp3CoeffWise(v, Jexp3QuatCoeffWise);
735  typename SE3::Matrix3 Jlog;
736  Jlog3(theta, v, Jlog);
738 
739  // std::cout << "Jexp3QuatCoeffWise\n" << Jexp3QuatCoeffWise << std::endl;
740  // std::cout << "Jlog\n" << Jlog << std::endl;
741 
742  // if(quat_map.w() >= 0.) // comes from the log3 for quaternions which may change the
743  // sign.
744  if (quat_map.coeffs()[3] >= Scalar(0)) // comes from the log3 for quaternions which may change
745  // the sign.
746  PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J).template bottomRightCorner<4, 3>().noalias() =
747  Jexp3QuatCoeffWise * Jlog;
748  else
749  PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J).template bottomRightCorner<4, 3>().noalias() =
750  -Jexp3QuatCoeffWise * Jlog;
751  }
752 
753  template<class Config_t, class Tangent_t, class JacobianOut_t>
754  static void dIntegrate_dq_impl(
755  const Eigen::MatrixBase<Config_t> & /*q*/,
756  const Eigen::MatrixBase<Tangent_t> & v,
757  const Eigen::MatrixBase<JacobianOut_t> & J,
758  const AssignmentOperatorType op = SETTO)
759  {
760  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J);
761 
762  switch (op)
763  {
764  case SETTO:
765  Jout = exp6(MotionRef<const Tangent_t>(v.derived())).toDualActionMatrix().transpose();
766  break;
767  case ADDTO:
768  Jout += exp6(MotionRef<const Tangent_t>(v.derived())).toDualActionMatrix().transpose();
769  break;
770  case RMTO:
771  Jout -= exp6(MotionRef<const Tangent_t>(v.derived())).toDualActionMatrix().transpose();
772  break;
773  default:
774  assert(false && "Wrong Op requesed value");
775  break;
776  }
777  }
778 
779  template<class Config_t, class Tangent_t, class JacobianOut_t>
780  static void dIntegrate_dv_impl(
781  const Eigen::MatrixBase<Config_t> & /*q*/,
782  const Eigen::MatrixBase<Tangent_t> & v,
783  const Eigen::MatrixBase<JacobianOut_t> & J,
784  const AssignmentOperatorType op = SETTO)
785  {
786  switch (op)
787  {
788  case SETTO:
789  Jexp6<SETTO>(MotionRef<const Tangent_t>(v.derived()), J.derived());
790  break;
791  case ADDTO:
792  Jexp6<ADDTO>(MotionRef<const Tangent_t>(v.derived()), J.derived());
793  break;
794  case RMTO:
795  Jexp6<RMTO>(MotionRef<const Tangent_t>(v.derived()), J.derived());
796  break;
797  default:
798  assert(false && "Wrong Op requesed value");
799  break;
800  }
801  }
802 
803  template<class Config_t, class Tangent_t, class JacobianIn_t, class JacobianOut_t>
805  const Eigen::MatrixBase<Config_t> & /*q*/,
806  const Eigen::MatrixBase<Tangent_t> & v,
807  const Eigen::MatrixBase<JacobianIn_t> & Jin,
808  const Eigen::MatrixBase<JacobianOut_t> & J_out)
809  {
810  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J_out);
811  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
812  Jtmp6 = exp6(MotionRef<const Tangent_t>(v.derived())).toDualActionMatrix().transpose();
813 
814  Jout.template topRows<3>().noalias() =
815  Jtmp6.template topLeftCorner<3, 3>() * Jin.template topRows<3>();
816  Jout.template topRows<3>().noalias() +=
817  Jtmp6.template topRightCorner<3, 3>() * Jin.template bottomRows<3>();
818  Jout.template bottomRows<3>().noalias() =
819  Jtmp6.template bottomRightCorner<3, 3>() * Jin.template bottomRows<3>();
820  }
821 
822  template<class Config_t, class Tangent_t, class JacobianIn_t, class JacobianOut_t>
824  const Eigen::MatrixBase<Config_t> & /*q*/,
825  const Eigen::MatrixBase<Tangent_t> & v,
826  const Eigen::MatrixBase<JacobianIn_t> & Jin,
827  const Eigen::MatrixBase<JacobianOut_t> & J_out)
828  {
829  JacobianOut_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(JacobianOut_t, J_out);
832  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
833  Jexp6<SETTO>(MotionRef<const Tangent_t>(v.derived()), Jtmp6);
835 
836  Jout.template topRows<3>().noalias() =
837  Jtmp6.template topLeftCorner<3, 3>() * Jin.template topRows<3>();
838  Jout.template topRows<3>().noalias() +=
839  Jtmp6.template topRightCorner<3, 3>() * Jin.template bottomRows<3>();
840  Jout.template bottomRows<3>().noalias() =
841  Jtmp6.template bottomRightCorner<3, 3>() * Jin.template bottomRows<3>();
842  }
843 
844  template<class Config_t, class Tangent_t, class Jacobian_t>
846  const Eigen::MatrixBase<Config_t> & /*q*/,
847  const Eigen::MatrixBase<Tangent_t> & v,
848  const Eigen::MatrixBase<Jacobian_t> & J_out)
849  {
850  Jacobian_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J_out);
851  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
852  Jtmp6 = exp6(MotionRef<const Tangent_t>(v.derived())).toDualActionMatrix().transpose();
853 
854  // Aliasing
855  Jout.template topRows<3>() =
856  Jtmp6.template topLeftCorner<3, 3>() * Jout.template topRows<3>();
857  Jout.template topRows<3>().noalias() +=
858  Jtmp6.template topRightCorner<3, 3>() * Jout.template bottomRows<3>();
859  Jout.template bottomRows<3>() =
860  Jtmp6.template bottomRightCorner<3, 3>() * Jout.template bottomRows<3>();
861  }
862 
863  template<class Config_t, class Tangent_t, class Jacobian_t>
865  const Eigen::MatrixBase<Config_t> & /*q*/,
866  const Eigen::MatrixBase<Tangent_t> & v,
867  const Eigen::MatrixBase<Jacobian_t> & J_out)
868  {
869  Jacobian_t & Jout = PINOCCHIO_EIGEN_CONST_CAST(Jacobian_t, J_out);
872  Eigen::Matrix<Scalar, 6, 6> Jtmp6;
873  Jexp6<SETTO>(MotionRef<const Tangent_t>(v.derived()), Jtmp6);
875 
876  Jout.template topRows<3>() =
877  Jtmp6.template topLeftCorner<3, 3>() * Jout.template topRows<3>();
878  Jout.template topRows<3>().noalias() +=
879  Jtmp6.template topRightCorner<3, 3>() * Jout.template bottomRows<3>();
880  Jout.template bottomRows<3>() =
881  Jtmp6.template bottomRightCorner<3, 3>() * Jout.template bottomRows<3>();
882  }
883 
884  template<class ConfigL_t, class ConfigR_t>
886  const Eigen::MatrixBase<ConfigL_t> & q0, const Eigen::MatrixBase<ConfigR_t> & q1)
887  {
891  difference_impl(q0, q1, t);
893  return t.squaredNorm();
894  }
895 
896  template<class Config_t>
897  static void normalize_impl(const Eigen::MatrixBase<Config_t> & qout)
898  {
899  pinocchio::normalize(qout.const_cast_derived().template tail<4>());
900  }
901 
902  template<class Config_t>
903  static bool isNormalized_impl(const Eigen::MatrixBase<Config_t> & qin, const Scalar & prec)
904  {
905  Scalar norm = Scalar(qin.template tail<4>().norm());
906  using std::abs;
907  return abs(norm - Scalar(1.0)) < prec;
908  }
909 
910  template<class Config_t>
911  static void random_impl(const Eigen::MatrixBase<Config_t> & qout)
912  {
913  R3crossSO3_t().random(qout);
914  }
915 
916  template<class ConfigL_t, class ConfigR_t, class ConfigOut_t>
918  const Eigen::MatrixBase<ConfigL_t> & lower,
919  const Eigen::MatrixBase<ConfigR_t> & upper,
920  const Eigen::MatrixBase<ConfigOut_t> & qout)
921  {
922  R3crossSO3_t().randomConfiguration(lower, upper, qout);
923  }
924 
925  template<class ConfigL_t, class ConfigR_t>
927  const Eigen::MatrixBase<ConfigL_t> & q0,
928  const Eigen::MatrixBase<ConfigR_t> & q1,
929  const Scalar & prec)
930  {
931  return R3crossSO3_t().isSameConfiguration(q0, q1, prec);
932  }
933  }; // struct SpecialEuclideanOperationTpl<3>
934 
935 } // namespace pinocchio
936 
937 #endif // ifndef __pinocchio_multibody_liegroup_special_euclidean_operation_hpp__
pinocchio::CartesianProductOperation
Definition: cartesian-product.hpp:40
cassie-simulation.qout
def qout
Definition: cassie-simulation.py:251
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::toInverseActionMatrix
static void toInverseActionMatrix(const Eigen::MatrixBase< Matrix2Like > &R, const Eigen::MatrixBase< Vector2Like > &t, const Eigen::MatrixBase< Matrix3Like > &M, const AssignmentOperatorType op)
Definition: special-euclidean.hpp:93
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dIntegrateTransport_dq_impl
static void dIntegrateTransport_dq_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianIn_t > &Jin, const Eigen::MatrixBase< JacobianOut_t > &J_out)
Definition: special-euclidean.hpp:384
pinocchio::SpecialEuclideanOperationTpl
Definition: special-euclidean.hpp:25
quaternion.hpp
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::integrate_impl
static void integrate_impl(const Eigen::MatrixBase< ConfigIn_t > &q, const Eigen::MatrixBase< Velocity_t > &v, const Eigen::MatrixBase< ConfigOut_t > &qout)
Definition: special-euclidean.hpp:289
fwd.hpp
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::isSameConfiguration_impl
static bool isSameConfiguration_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1, const Scalar &prec)
Definition: special-euclidean.hpp:511
pinocchio::SE3Tpl< Scalar, Options >::Vector3
traits< SE3Tpl >::Vector3 Vector3
Definition: spatial/se3-tpl.hpp:55
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::integrateCoeffWiseJacobian_impl
static void integrateCoeffWiseJacobian_impl(const Eigen::MatrixBase< Config_t > &q, const Eigen::MatrixBase< Jacobian_t > &J)
Definition: special-euclidean.hpp:708
static-if.hpp
pinocchio::RMTO
@ RMTO
Definition: fwd.hpp:134
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::log
static void log(const Eigen::MatrixBase< Matrix2Like > &R, const Eigen::MatrixBase< Vector2Like > &p, const Eigen::MatrixBase< TangentVector > &v)
Definition: special-euclidean.hpp:132
V
V
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::normalize_impl
static void normalize_impl(const Eigen::MatrixBase< Config_t > &qout)
Definition: special-euclidean.hpp:897
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::integrate_impl
static void integrate_impl(const Eigen::MatrixBase< ConfigIn_t > &q, const Eigen::MatrixBase< Velocity_t > &v, const Eigen::MatrixBase< ConfigOut_t > &qout)
Definition: special-euclidean.hpp:664
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dIntegrateTransport_dq_impl
static void dIntegrateTransport_dq_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianIn_t > &Jin, const Eigen::MatrixBase< JacobianOut_t > &J_out)
Definition: special-euclidean.hpp:804
pinocchio::SETTO
@ SETTO
Definition: fwd.hpp:132
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.
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dIntegrate_dv_impl
static void dIntegrate_dv_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianOut_t > &J, const AssignmentOperatorType op=SETTO)
Definition: special-euclidean.hpp:780
macros.hpp
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::nq
static Index nq()
Definition: special-euclidean.hpp:207
autodiff-rnea.dv
dv
Definition: autodiff-rnea.py:27
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dIntegrateTransport_dv_impl
static void dIntegrateTransport_dv_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< Jacobian_t > &J_out)
Definition: special-euclidean.hpp:864
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::ConstQuaternionMap_t
Eigen::Map< const Quaternion_t > ConstQuaternionMap_t
Definition: special-euclidean.hpp:562
cartesian-product.hpp
cassie-simulation.qin
qin
Definition: cassie-simulation.py:249
quat
quat
pinocchio::LieGroupBase< SpecialEuclideanOperationTpl< 2, _Scalar, _Options > >::TangentVector_t
Eigen::Matrix< Scalar, NV, 1, Options > TangentVector_t
Definition: liegroup-base.hpp:57
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dIntegrateTransport_dq_impl
static void dIntegrateTransport_dq_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< Jacobian_t > &J_out)
Definition: special-euclidean.hpp:845
pinocchio::Options
Options
Definition: joint-configuration.hpp:1082
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dIntegrate_dq_impl
static void dIntegrate_dq_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianOut_t > &J, const AssignmentOperatorType op=SETTO)
Definition: special-euclidean.hpp:324
pinocchio::internal::GT
@ GT
Definition: utils/static-if.hpp:21
pinocchio.explog.exp
def exp(x)
Definition: bindings/python/pinocchio/explog.py:13
pinocchio::Jexp6
void Jexp6(const MotionDense< MotionDerived > &nu, const Eigen::MatrixBase< Matrix6Like > &Jexp)
Derivative of exp6 Computed as the inverse of Jlog6.
Definition: spatial/explog.hpp:496
pinocchio::nv
int nv(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointNvVisitor to get the dimension of the joint tangent space.
PINOCCHIO_EIGEN_CONST_CAST
#define PINOCCHIO_EIGEN_CONST_CAST(TYPE, OBJ)
Macro for an automatic const_cast.
Definition: eigen-macros.hpp:51
pinocchio::SE3Tpl< Scalar, Options >
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::normalize_impl
static void normalize_impl(const Eigen::MatrixBase< Config_t > &qout)
Definition: special-euclidean.hpp:482
pinocchio::nq
int nq(const JointModelTpl< Scalar, Options, JointCollectionTpl > &jmodel)
Visit a JointModelTpl through JointNqVisitor to get the dimension of the joint configuration space.
pinocchio::LieGroupBase< SpecialEuclideanOperationTpl< 2, _Scalar, _Options > >::Index
int Index
Definition: liegroup-base.hpp:47
pinocchio::quaternion::isNormalized
bool isNormalized(const Eigen::QuaternionBase< Quaternion > &quat, const typename Quaternion::Coefficients::RealScalar &prec=Eigen::NumTraits< typename Quaternion::Coefficients::RealScalar >::dummy_precision())
Check whether the input quaternion is Normalized within the given precision.
Definition: math/quaternion.hpp:229
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::integrateCoeffWiseJacobian_impl
static void integrateCoeffWiseJacobian_impl(const Eigen::MatrixBase< Config_t > &q, const Eigen::MatrixBase< Jacobian_t > &J)
Definition: special-euclidean.hpp:309
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::isSameConfiguration_impl
static bool isSameConfiguration_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1, const Scalar &prec)
Definition: special-euclidean.hpp:926
PINOCCHIO_COMPILER_DIAGNOSTIC_POP
#define PINOCCHIO_COMPILER_DIAGNOSTIC_POP
Definition: include/pinocchio/macros.hpp:130
pinocchio::python::Scalar
context::Scalar Scalar
Definition: admm-solver.cpp:29
inverse-kinematics-3d.J
J
Definition: inverse-kinematics-3d.py:28
PINOCCHIO_LIE_GROUP_TPL_PUBLIC_INTERFACE
#define PINOCCHIO_LIE_GROUP_TPL_PUBLIC_INTERFACE(Derived)
Definition: liegroup-base.hpp:40
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::random_impl
static void random_impl(const Eigen::MatrixBase< Config_t > &qout)
Definition: special-euclidean.hpp:911
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::random_impl
static void random_impl(const Eigen::MatrixBase< Config_t > &qout)
Definition: special-euclidean.hpp:496
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::R3crossSO3_t
CartesianProductOperation< VectorSpaceOperationTpl< 3, Scalar, Options >, SpecialOrthogonalOperationTpl< 3, Scalar, Options > > R3crossSO3_t
Definition: special-euclidean.hpp:558
pinocchio::Jlog6
void Jlog6(const SE3Tpl< Scalar, Options > &M, const Eigen::MatrixBase< Matrix6Like > &Jlog)
Derivative of log6.
Definition: spatial/explog.hpp:668
R
R
pinocchio::internal::if_then_else
if_then_else_impl< LhsType, RhsType, ThenType, ElseType >::ReturnType if_then_else(const ComparisonOperators op, const LhsType &lhs_value, const RhsType &rhs_value, const ThenType &then_value, const ElseType &else_value)
Definition: utils/static-if.hpp:89
pinocchio::SINCOS
void SINCOS(const S1 &a, S2 *sa, S3 *ca)
Computes sin/cos values of a given input scalar.
Definition: sincos.hpp:27
pinocchio::traits< SpecialEuclideanOperationTpl< 2, _Scalar, _Options > >::Scalar
_Scalar Scalar
Definition: special-euclidean.hpp:37
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dIntegrate_dv_impl
static void dIntegrate_dv_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianOut_t > &J, const AssignmentOperatorType op=SETTO)
Definition: special-euclidean.hpp:343
pinocchio::ARG0
@ ARG0
Definition: fwd.hpp:123
pinocchio::traits< SpecialEuclideanOperationTpl< 3, _Scalar, _Options > >::Scalar
_Scalar Scalar
Definition: special-euclidean.hpp:539
PINOCCHIO_COMPILER_DIAGNOSTIC_IGNORED_MAYBE_UNINITIALIZED
#define PINOCCHIO_COMPILER_DIAGNOSTIC_IGNORED_MAYBE_UNINITIALIZED
Definition: include/pinocchio/macros.hpp:134
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::nv
static Index nv()
Get dimension of Lie Group tangent space.
Definition: special-euclidean.hpp:576
pinocchio::Jlog3
void Jlog3(const Scalar &theta, const Eigen::MatrixBase< Vector3Like > &log, const Eigen::MatrixBase< Matrix3Like > &Jlog)
Derivative of log3.
Definition: spatial/explog.hpp:240
pinocchio::SE3Tpl< Scalar, Options >::Matrix3
traits< SE3Tpl >::Matrix3 Matrix3
Definition: spatial/se3-tpl.hpp:56
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::name
static std::string name()
Definition: special-euclidean.hpp:224
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::randomConfiguration_impl
static void randomConfiguration_impl(const Eigen::MatrixBase< ConfigL_t > &lower, const Eigen::MatrixBase< ConfigR_t > &upper, const Eigen::MatrixBase< ConfigOut_t > &qout)
Definition: special-euclidean.hpp:917
dpendulum.p
p
Definition: dpendulum.py:7
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dDifference_impl
static void dDifference_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1, const Eigen::MatrixBase< JacobianOut_t > &J)
Definition: special-euclidean.hpp:618
matrix.hpp
reachable-workspace-with-collisions.alpha
float alpha
Definition: reachable-workspace-with-collisions.py:162
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::neutral
static ConfigVector_t neutral()
Definition: special-euclidean.hpp:217
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dIntegrateTransport_dv_impl
static void dIntegrateTransport_dv_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< Jacobian_t > &J)
Definition: special-euclidean.hpp:456
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::RealScalar
Eigen::NumTraits< Scalar >::Real RealScalar
Definition: special-euclidean.hpp:565
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::nq
static Index nq()
Definition: special-euclidean.hpp:571
pinocchio::LieGroupBase< SpecialEuclideanOperationTpl< 2, _Scalar, _Options > >::Scalar
traits< LieGroupDerived >::Scalar Scalar
Definition: liegroup-base.hpp:48
t2
dictionary t2
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::isNormalized_impl
static bool isNormalized_impl(const Eigen::MatrixBase< Config_t > &qin, const Scalar &prec)
Definition: special-euclidean.hpp:488
pinocchio::skew
void skew(const Eigen::MatrixBase< Vector3 > &v, const Eigen::MatrixBase< Matrix3 > &M)
Computes the skew representation of a given 3d vector, i.e. the antisymmetric matrix representation o...
Definition: skew.hpp:22
se3.hpp
pinocchio::exp6
SE3Tpl< typename MotionDerived::Scalar, PINOCCHIO_EIGEN_PLAIN_TYPE(typename MotionDerived::Vector3)::Options > exp6(const MotionDense< MotionDerived > &nu)
Exp: se3 -> SE3.
Definition: spatial/explog.hpp:347
pinocchio::SpecialOrthogonalOperationTpl
Definition: special-orthogonal.hpp:18
pinocchio::quaternion::log3
Eigen::Matrix< typename QuaternionLike::Scalar, 3, 1, PINOCCHIO_EIGEN_PLAIN_TYPE(typename QuaternionLike::Vector3)::Options > log3(const Eigen::QuaternionBase< QuaternionLike > &quat, typename QuaternionLike::Scalar &theta)
Same as log3 but with a unit quaternion as input.
Definition: explog-quaternion.hpp:211
special-orthogonal.hpp
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::neutral
static ConfigVector_t neutral()
Definition: special-euclidean.hpp:581
pinocchio::ARG1
@ ARG1
Definition: fwd.hpp:124
vector-space.hpp
pinocchio::q0
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorIn1 > & q0
Definition: joint-configuration.hpp:1137
pinocchio::AssignmentOperatorType
AssignmentOperatorType
Definition: fwd.hpp:130
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::SO2_t
SpecialOrthogonalOperationTpl< 2, Scalar, Options > SO2_t
Definition: special-euclidean.hpp:54
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::difference_impl
static void difference_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1, const Eigen::MatrixBase< Tangent_t > &d)
Definition: special-euclidean.hpp:595
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::R2crossSO2_t
CartesianProductOperation< R2_t, SO2_t > R2crossSO2_t
Definition: special-euclidean.hpp:55
PINOCCHIO_COMPILER_DIAGNOSTIC_PUSH
#define PINOCCHIO_COMPILER_DIAGNOSTIC_PUSH
macros for pragma push/pop/ignore deprecated warnings
Definition: include/pinocchio/macros.hpp:129
pinocchio::ADDTO
@ ADDTO
Definition: fwd.hpp:133
M
M
pinocchio::context::Vector3
Eigen::Matrix< Scalar, 3, 1, Options > Vector3
Definition: context/generic.hpp:53
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dIntegrateTransport_dq_impl
static void dIntegrateTransport_dq_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< Jacobian_t > &J)
Definition: special-euclidean.hpp:434
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::nv
static Index nv()
Get dimension of Lie Group tangent space.
Definition: special-euclidean.hpp:212
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::Matrix2
Eigen::Matrix< Scalar, 2, 2, Options > Matrix2
Definition: special-euclidean.hpp:57
pinocchio::LieGroupBase
Definition: liegroup-base.hpp:44
pinocchio::q
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > & q
Definition: joint-configuration.hpp:1083
pinocchio::LieGroupBase< SpecialEuclideanOperationTpl< 2, _Scalar, _Options > >::JacobianMatrix_t
Eigen::Matrix< Scalar, NV, NV, Options > JacobianMatrix_t
Definition: liegroup-base.hpp:58
pinocchio::v
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorType > const Eigen::MatrixBase< TangentVectorType > & v
Definition: joint-configuration.hpp:1084
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::Transformation_t
SE3Tpl< Scalar, Options > Transformation_t
Definition: special-euclidean.hpp:563
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::QuaternionMap_t
Eigen::Map< Quaternion_t > QuaternionMap_t
Definition: special-euclidean.hpp:561
PINOCCHIO_DEFAULT_QUATERNION_NORM_TOLERANCE_VALUE
#define PINOCCHIO_DEFAULT_QUATERNION_NORM_TOLERANCE_VALUE
Definition: math/quaternion.hpp:9
pinocchio::quaternion::firstOrderNormalize
void firstOrderNormalize(const Eigen::QuaternionBase< D > &q)
Definition: math/quaternion.hpp:89
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::squaredDistance_impl
static Scalar squaredDistance_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1)
Definition: special-euclidean.hpp:885
pinocchio::log6
MotionTpl< Scalar, Options > log6(const SE3Tpl< Scalar, Options > &M)
Log: SE3 -> se3.
Definition: spatial/explog.hpp:435
pinocchio::internal::LT
@ LT
Definition: utils/static-if.hpp:17
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::isNormalized_impl
static bool isNormalized_impl(const Eigen::MatrixBase< Config_t > &qin, const Scalar &prec)
Definition: special-euclidean.hpp:903
pinocchio::quaternion::exp6
void exp6(const MotionDense< MotionDerived > &motion, Eigen::MatrixBase< Config_t > &qout)
The se3 -> SE3 exponential map, using quaternions to represent the output rotation.
Definition: explog-quaternion.hpp:92
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::Quaternion_t
Eigen::Quaternion< Scalar, Options > Quaternion_t
Definition: special-euclidean.hpp:560
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dIntegrateTransport_dv_impl
static void dIntegrateTransport_dv_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianIn_t > &Jin, const Eigen::MatrixBase< JacobianOut_t > &J_out)
Definition: special-euclidean.hpp:823
pinocchio::q1
JointCollectionTpl const Eigen::MatrixBase< ConfigVectorIn1 > const Eigen::MatrixBase< ConfigVectorIn2 > & q1
Definition: joint-configuration.hpp:1138
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::R2_t
VectorSpaceOperationTpl< 2, Scalar, Options > R2_t
Definition: special-euclidean.hpp:53
pinocchio.explog.log
def log(x)
Definition: bindings/python/pinocchio/explog.py:29
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::difference_impl
static void difference_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1, const Eigen::MatrixBase< Tangent_t > &d)
Definition: special-euclidean.hpp:230
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::Jlog
static void Jlog(const Eigen::MatrixBase< Matrix2Like > &R, const Eigen::MatrixBase< Vector2Like > &p, const Eigen::MatrixBase< JacobianOutLike > &J)
Definition: special-euclidean.hpp:163
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::SE3
SE3Tpl< Scalar, Options > SE3
Definition: special-euclidean.hpp:564
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::dIntegrate_dq_impl
static void dIntegrate_dq_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianOut_t > &J, const AssignmentOperatorType op=SETTO)
Definition: special-euclidean.hpp:754
pinocchio::traits
Common traits structure to fully define base classes for CRTP.
Definition: fwd.hpp:71
t
Transform3f t
pinocchio::VectorSpaceOperationTpl
Definition: vector-space.hpp:16
pinocchio::PINOCCHIO_EIGEN_PLAIN_TYPE
PINOCCHIO_EIGEN_PLAIN_TYPE(ConfigVectorType) integrate(const ModelTpl< Scalar
Integrate a configuration vector for the specified model for a tangent vector during one unit time.
PINOCCHIO_ASSERT_MATRIX_SPECIFIC_SIZE
#define PINOCCHIO_ASSERT_MATRIX_SPECIFIC_SIZE(type, M, nrows, ncols)
Ensure that a matrix (or vector) is of correct size (compile-time and run-time assertion)
Definition: include/pinocchio/macros.hpp:75
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::exp
static void exp(const Eigen::MatrixBase< TangentVector > &v, const Eigen::MatrixBase< Matrix2Like > &R, const Eigen::MatrixBase< Vector2Like > &t)
Definition: special-euclidean.hpp:61
pinocchio::SpecialEuclideanOperationTpl< 3, _Scalar, _Options >::name
static std::string name()
Definition: special-euclidean.hpp:589
pinocchio::quaternion::Jexp3CoeffWise
void Jexp3CoeffWise(const Eigen::MatrixBase< Vector3Like > &v, const Eigen::MatrixBase< Matrix43Like > &Jexp)
Derivative of where is a small perturbation of at identity.
Definition: explog-quaternion.hpp:292
pinocchio::SE3
SE3Tpl< context::Scalar, context::Options > SE3
Definition: spatial/fwd.hpp:64
pinocchio::MotionTpl
Definition: context/casadi.hpp:27
dcrba.NV
NV
Definition: dcrba.py:536
dpendulum.NQ
int NQ
Definition: dpendulum.py:9
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::forwardKinematics
static void forwardKinematics(const Eigen::MatrixBase< Matrix2Like > &R, const Eigen::MatrixBase< Vector2Like > &t, const Eigen::MatrixBase< Vector4Like > &q)
Definition: special-euclidean.hpp:521
pinocchio::MotionRef
Definition: context/casadi.hpp:38
pinocchio::LieGroupBase< SpecialEuclideanOperationTpl< 2, _Scalar, _Options > >::ConfigVector_t
Eigen::Matrix< Scalar, NQ, 1, Options > ConfigVector_t
Definition: liegroup-base.hpp:56
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::Vector2
Eigen::Matrix< Scalar, 2, 1, Options > Vector2
Definition: special-euclidean.hpp:58
liegroup-base.hpp
d
FCL_REAL d
n
Vec3f n
pinocchio::normalize
void normalize(const ModelTpl< Scalar, Options, JointCollectionTpl > &model, const Eigen::MatrixBase< ConfigVectorType > &qout)
Normalize a configuration vector.
Definition: joint-configuration.hpp:887
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dIntegrateTransport_dv_impl
static void dIntegrateTransport_dv_impl(const Eigen::MatrixBase< Config_t > &, const Eigen::MatrixBase< Tangent_t > &v, const Eigen::MatrixBase< JacobianIn_t > &Jin, const Eigen::MatrixBase< JacobianOut_t > &J_out)
Definition: special-euclidean.hpp:407
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::dDifference_impl
void dDifference_impl(const Eigen::MatrixBase< ConfigL_t > &q0, const Eigen::MatrixBase< ConfigR_t > &q1, const Eigen::MatrixBase< JacobianOut_t > &J) const
Definition: special-euclidean.hpp:249
pinocchio
Main pinocchio namespace.
Definition: timings.cpp:27
pinocchio::SpecialEuclideanOperationTpl< 2, _Scalar, _Options >::randomConfiguration_impl
static void randomConfiguration_impl(const Eigen::MatrixBase< ConfigL_t > &lower, const Eigen::MatrixBase< ConfigR_t > &upper, const Eigen::MatrixBase< ConfigOut_t > &qout)
Definition: special-euclidean.hpp:502


pinocchio
Author(s):
autogenerated on Sun Dec 22 2024 03:41:12