pinocchio  UNKNOWN
symmetric3.hpp
1 //
2 // Copyright (c) 2014-2017 CNRS
3 //
4 // This file is part of Pinocchio
5 // Pinocchio is free software: you can redistribute it
6 // and/or modify it under the terms of the GNU Lesser General Public
7 // License as published by the Free Software Foundation, either version
8 // 3 of the License, or (at your option) any later version.
9 //
10 // Pinocchio is distributed in the hope that it will be
11 // useful, but WITHOUT ANY WARRANTY; without even the implied warranty
12 // of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // General Lesser Public License for more details. You should have
14 // received a copy of the GNU Lesser General Public License along with
15 // Pinocchio If not, see
16 // <http://www.gnu.org/licenses/>.
17 //
18 // This file is originally copied from metapod/tools/spatial/lti.hh.
19 // Authors: Olivier Stasse (LAAS, CNRS) and Sébastien Barthélémy (Aldebaran Robotics)
20 // The file was modified in pinocchio by Nicolas Mansard (LAAS, CNRS)
21 //
22 // metapod is free software, distributed under the terms of the GNU Lesser
23 // General Public License as published by
24 // the Free Software Foundation, either version 3 of the License, or
25 // (at your option) any later version.
26 
27 #ifndef __se3_symmetric3_hpp__
28 #define __se3_symmetric3_hpp__
29 
30 #include "pinocchio/macros.hpp"
31 
32 namespace se3
33 {
34 
35  template<typename _Scalar, int _Options>
36  class Symmetric3Tpl
37  {
38  public:
39  typedef _Scalar Scalar;
40  enum { Options = _Options };
41  typedef Eigen::Matrix<Scalar,3,1,Options> Vector3;
42  typedef Eigen::Matrix<Scalar,6,1,Options> Vector6;
43  typedef Eigen::Matrix<Scalar,3,3,Options> Matrix3;
44  typedef Eigen::Matrix<Scalar,2,2,Options> Matrix2;
45  typedef Eigen::Matrix<Scalar,3,2,Options> Matrix32;
46 
47  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
48 
49  public:
50  Symmetric3Tpl(): data_() {}
51 
52 // template<typename D>
53 // explicit Symmetric3Tpl(const Eigen::MatrixBase<D> & I)
54 // {
55 // EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(D,3,3);
56 // assert( (I-I.transpose()).isMuchSmallerThan(I) );
57 // data_(0) = I(0,0);
58 // data_(1) = I(1,0); data_(2) = I(1,1);
59 // data_(3) = I(2,0); data_(4) = I(2,1); data_(5) = I(2,2);
60 // }
61  template<typename Sc,int N,int Opt>
62  explicit Symmetric3Tpl(const Eigen::Matrix<Sc,N,N,Opt> & I)
63  {
64  EIGEN_STATIC_ASSERT(N==3,THIS_METHOD_IS_ONLY_FOR_MATRICES_OF_A_SPECIFIC_SIZE)
65  assert( (I-I.transpose()).isMuchSmallerThan(I) );
66  data_(0) = I(0,0);
67  data_(1) = I(1,0); data_(2) = I(1,1);
68  data_(3) = I(2,0); data_(4) = I(2,1); data_(5) = I(2,2);
69  }
70 
71  explicit Symmetric3Tpl(const Vector6 & I) : data_(I) {}
72 
73  Symmetric3Tpl(const Scalar & a0, const Scalar & a1, const Scalar & a2,
74  const Scalar & a3, const Scalar & a4, const Scalar & a5)
75  { data_ << a0,a1,a2,a3,a4,a5; }
76 
77  static Symmetric3Tpl Zero() { return Symmetric3Tpl(Vector6::Zero()); }
78  void setZero() { data_.setZero(); }
79 
80  static Symmetric3Tpl Random() { return RandomPositive(); }
81  void setRandom()
82  {
83  Scalar
84  a = Scalar(std::rand())/RAND_MAX*2.0-1.0,
85  b = Scalar(std::rand())/RAND_MAX*2.0-1.0,
86  c = Scalar(std::rand())/RAND_MAX*2.0-1.0,
87  d = Scalar(std::rand())/RAND_MAX*2.0-1.0,
88  e = Scalar(std::rand())/RAND_MAX*2.0-1.0,
89  f = Scalar(std::rand())/RAND_MAX*2.0-1.0;
90 
91  data_ << a, b, c, d, e, f;
92  }
93 
94  static Symmetric3Tpl Identity() { return Symmetric3Tpl(1, 0, 1, 0, 0, 1); }
95  void setIdentity() { data_ << 1, 0, 1, 0, 0, 1; }
96 
97  /* Requiered by Inertia::operator== */
98  bool operator== (const Symmetric3Tpl & S2) const { return data_ == S2.data_; }
99 
100  bool isApprox(const Symmetric3Tpl & other,
101  const Scalar & prec = Eigen::NumTraits<Scalar>::dummy_precision()) const
102  { return data_.isApprox(other.data_,prec); }
103 
104  void fill(const Scalar value) { data_.fill(value); }
105 
106  struct SkewSquare
107  {
108  const Vector3 & v;
109  SkewSquare( const Vector3 & v ) : v(v) {}
110  operator Symmetric3Tpl () const
111  {
112  const Scalar & x = v[0], & y = v[1], & z = v[2];
113  return Symmetric3Tpl( -y*y-z*z,
114  x*y , -x*x-z*z,
115  x*z , y*z , -x*x-y*y );
116  }
117  }; // struct SkewSquare
118 
119  Symmetric3Tpl operator- (const SkewSquare & v) const
120  {
121  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
122  return Symmetric3Tpl(data_[0]+y*y+z*z,
123  data_[1]-x*y,data_[2]+x*x+z*z,
124  data_[3]-x*z,data_[4]-y*z,data_[5]+x*x+y*y);
125  }
126 
127  Symmetric3Tpl& operator-= (const SkewSquare & v)
128  {
129  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
130  data_[0]+=y*y+z*z;
131  data_[1]-=x*y; data_[2]+=x*x+z*z;
132  data_[3]-=x*z; data_[4]-=y*z; data_[5]+=x*x+y*y;
133  return *this;
134  }
135 
136  template<typename D>
137  friend Matrix3 operator- (const Symmetric3Tpl & S, const Eigen::MatrixBase <D> & M)
138  {
139  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(D,3,3);
140  Matrix3 result (S.matrix());
141  result -= M;
142 
143  return result;
144  }
145 
147  {
148  const Scalar & m;
149  const Vector3 & v;
150 
151  AlphaSkewSquare(const Scalar & m, const SkewSquare & v) : m(m),v(v.v) {}
152  AlphaSkewSquare(const Scalar & m, const Vector3 & v) : m(m),v(v) {}
153 
154  operator Symmetric3Tpl () const
155  {
156  const Scalar & x = v[0], & y = v[1], & z = v[2];
157  return Symmetric3Tpl(-m*(y*y+z*z),
158  m* x*y,-m*(x*x+z*z),
159  m* x*z,m* y*z,-m*(x*x+y*y));
160  }
161  };
162 
163  friend AlphaSkewSquare operator* (const Scalar & m, const SkewSquare & sk)
164  { return AlphaSkewSquare(m,sk); }
165 
166  Symmetric3Tpl operator- (const AlphaSkewSquare & v) const
167  {
168  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
169  return Symmetric3Tpl(data_[0]+v.m*(y*y+z*z),
170  data_[1]-v.m* x*y, data_[2]+v.m*(x*x+z*z),
171  data_[3]-v.m* x*z, data_[4]-v.m* y*z,
172  data_[5]+v.m*(x*x+y*y));
173  }
174 
175  Symmetric3Tpl& operator-= (const AlphaSkewSquare & v)
176  {
177  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
178  data_[0]+=v.m*(y*y+z*z);
179  data_[1]-=v.m* x*y; data_[2]+=v.m*(x*x+z*z);
180  data_[3]-=v.m* x*z; data_[4]-=v.m* y*z; data_[5]+=v.m*(x*x+y*y);
181  return *this;
182  }
183 
184  const Vector6 & data () const {return data_;}
185  Vector6 & data () {return data_;}
186 
187  // static Symmetric3Tpl SkewSq( const Vector3 & v )
188  // {
189  // const Scalar & x = v[0], & y = v[1], & z = v[2];
190  // return Symmetric3Tpl(-y*y-z*z,
191  // x*y, -x*x-z*z,
192  // x*z, y*z, -x*x-y*y );
193  // }
194 
195  /* Shoot a positive definite matrix. */
196  static Symmetric3Tpl RandomPositive()
197  {
198  Scalar
199  a = Scalar(std::rand())/RAND_MAX*2.0-1.0,
200  b = Scalar(std::rand())/RAND_MAX*2.0-1.0,
201  c = Scalar(std::rand())/RAND_MAX*2.0-1.0,
202  d = Scalar(std::rand())/RAND_MAX*2.0-1.0,
203  e = Scalar(std::rand())/RAND_MAX*2.0-1.0,
204  f = Scalar(std::rand())/RAND_MAX*2.0-1.0;
205  return Symmetric3Tpl(a*a+b*b+d*d,
206  a*b+b*c+d*e, b*b+c*c+e*e,
207  a*d+b*e+d*f, b*d+c*e+e*f, d*d+e*e+f*f );
208  }
209 
210  Matrix3 matrix() const
211  {
212  Matrix3 res;
213  res(0,0) = data_(0); res(0,1) = data_(1); res(0,2) = data_(3);
214  res(1,0) = data_(1); res(1,1) = data_(2); res(1,2) = data_(4);
215  res(2,0) = data_(3); res(2,1) = data_(4); res(2,2) = data_(5);
216  return res;
217  }
218  operator Matrix3 () const { return matrix(); }
219 
220  Scalar vtiv (const Vector3 & v) const
221  {
222  const Scalar & x = v[0];
223  const Scalar & y = v[1];
224  const Scalar & z = v[2];
225 
226  const Scalar xx = x*x;
227  const Scalar xy = x*y;
228  const Scalar xz = x*z;
229  const Scalar yy = y*y;
230  const Scalar yz = y*z;
231  const Scalar zz = z*z;
232 
233  return data_(0)*xx + data_(2)*yy + data_(5)*zz + 2.*(data_(1)*xy + data_(3)*xz + data_(4)*yz);
234  }
235 
246  template<typename Vector3, typename Matrix3>
247  static void vxs(const Eigen::MatrixBase<Vector3> & v,
248  const Symmetric3Tpl & S3,
249  const Eigen::MatrixBase<Matrix3> & M)
250  {
251  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3,3);
252  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix3,3,3);
253 
254  const Scalar & a = S3.data()[0];
255  const Scalar & b = S3.data()[1];
256  const Scalar & c = S3.data()[2];
257  const Scalar & d = S3.data()[3];
258  const Scalar & e = S3.data()[4];
259  const Scalar & f = S3.data()[5];
260 
261 
262  const typename Vector3::RealScalar & v0 = v[0];
263  const typename Vector3::RealScalar & v1 = v[1];
264  const typename Vector3::RealScalar & v2 = v[2];
265 
266  Matrix3 & M_ = const_cast<Eigen::MatrixBase<Matrix3> &>(M).derived();
267  M_(0,0) = d * v1 - b * v2;
268  M_(1,0) = a * v2 - d * v0;
269  M_(2,0) = b * v0 - a * v1;
270 
271  M_(0,1) = e * v1 - c * v2;
272  M_(1,1) = b * v2 - e * v0;
273  M_(2,1) = c * v0 - b * v1;
274 
275  M_(0,2) = f * v1 - e * v2;
276  M_(1,2) = d * v2 - f * v0;
277  M_(2,2) = e * v0 - d * v1;
278  }
279 
290  template<typename Vector3>
291  Matrix3 vxs(const Eigen::MatrixBase<Vector3> & v) const
292  {
293  Matrix3 M;
294  vxs(v,*this,M);
295  return M;
296  }
297 
307  template<typename Vector3, typename Matrix3>
308  static void svx(const Eigen::MatrixBase<Vector3> & v,
309  const Symmetric3Tpl & S3,
310  const Eigen::MatrixBase<Matrix3> & M)
311  {
312  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3,3);
313  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix3,3,3);
314 
315  const Scalar & a = S3.data()[0];
316  const Scalar & b = S3.data()[1];
317  const Scalar & c = S3.data()[2];
318  const Scalar & d = S3.data()[3];
319  const Scalar & e = S3.data()[4];
320  const Scalar & f = S3.data()[5];
321 
322  const typename Vector3::RealScalar & v0 = v[0];
323  const typename Vector3::RealScalar & v1 = v[1];
324  const typename Vector3::RealScalar & v2 = v[2];
325 
326  Matrix3 & M_ = const_cast<Eigen::MatrixBase<Matrix3> &>(M).derived();
327  M_(0,0) = b * v2 - d * v1;
328  M_(1,0) = c * v2 - e * v1;
329  M_(2,0) = e * v2 - f * v1;
330 
331  M_(0,1) = d * v0 - a * v2;
332  M_(1,1) = e * v0 - b * v2;
333  M_(2,1) = f * v0 - d * v2;
334 
335  M_(0,2) = a * v1 - b * v0;
336  M_(1,2) = b * v1 - c * v0;
337  M_(2,2) = d * v1 - e * v0;
338  }
339 
348  template<typename Vector3>
349  Matrix3 svx(const Eigen::MatrixBase<Vector3> & v) const
350  {
351  Matrix3 M;
352  svx(v,*this,M);
353  return M;
354  }
355 
356  Symmetric3Tpl operator+(const Symmetric3Tpl & s2) const
357  {
358  return Symmetric3Tpl((data_+s2.data_).eval());
359  }
360 
361  Symmetric3Tpl & operator+=(const Symmetric3Tpl & s2)
362  {
363  data_ += s2.data_; return *this;
364  }
365 
366  template<typename V3in, typename V3out>
367  static void rhsMult(const Symmetric3Tpl & S3,
368  const Eigen::MatrixBase<V3in> & vin,
369  const Eigen::MatrixBase<V3out> & vout)
370  {
371  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(V3in,Vector3);
372  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(V3out,Vector3);
373 
374  Eigen::MatrixBase<V3out> & vout_ = const_cast<Eigen::MatrixBase<V3out>&>(vout);
375 
376  vout_[0] = S3.data_(0) * vin[0] + S3.data_(1) * vin[1] + S3.data_(3) * vin[2];
377  vout_[1] = S3.data_(1) * vin[0] + S3.data_(2) * vin[1] + S3.data_(4) * vin[2];
378  vout_[2] = S3.data_(3) * vin[0] + S3.data_(4) * vin[1] + S3.data_(5) * vin[2];
379  }
380 
381  template<typename V3>
382  Vector3 operator*(const Eigen::MatrixBase<V3> & v) const
383  {
384  Vector3 res;
385  rhsMult(*this,v,res);
386  return res;
387  }
388 
389  // Matrix3 operator*(const Matrix3 &a) const
390  // {
391  // Matrix3 r;
392  // for(unsigned int i=0; i<3; ++i)
393  // {
394  // r(0,i) = data_(0) * a(0,i) + data_(1) * a(1,i) + data_(3) * a(2,i);
395  // r(1,i) = data_(1) * a(0,i) + data_(2) * a(1,i) + data_(4) * a(2,i);
396  // r(2,i) = data_(3) * a(0,i) + data_(4) * a(1,i) + data_(5) * a(2,i);
397  // }
398  // return r;
399  // }
400 
401  const Scalar& operator()(const int &i,const int &j) const
402  {
403  return ((i!=2)&&(j!=2)) ? data_[i+j] : data_[i+j+1];
404  }
405 
406  Symmetric3Tpl operator-(const Matrix3 &S) const
407  {
408  assert( (S-S.transpose()).isMuchSmallerThan(S) );
409  return Symmetric3Tpl( data_(0)-S(0,0),
410  data_(1)-S(1,0), data_(2)-S(1,1),
411  data_(3)-S(2,0), data_(4)-S(2,1), data_(5)-S(2,2) );
412  }
413 
414  Symmetric3Tpl operator+(const Matrix3 &S) const
415  {
416  assert( (S-S.transpose()).isMuchSmallerThan(S) );
417  return Symmetric3Tpl( data_(0)+S(0,0),
418  data_(1)+S(1,0), data_(2)+S(1,1),
419  data_(3)+S(2,0), data_(4)+S(2,1), data_(5)+S(2,2) );
420  }
421 
422  /* --- Symmetric R*S*R' and R'*S*R products --- */
423  public: //private:
424 
427  Matrix32 decomposeltI() const
428  {
429  Matrix32 L;
430  L <<
431  data_(0) - data_(5), data_(1),
432  data_(1), data_(2) - data_(5),
433  2*data_(3), data_(4) + data_(4);
434  return L;
435  }
436 
437  /* R*S*R' */
438  template<typename D>
439  Symmetric3Tpl rotate(const Eigen::MatrixBase<D> & R) const
440  {
441  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(D,3,3);
442  assert( (R.transpose()*R).isApprox(Matrix3::Identity()) );
443 
444  Symmetric3Tpl Sres;
445 
446  // 4 a
447  const Matrix32 L( decomposeltI() );
448 
449  // Y = R' L ===> (12 m + 8 a)
450  const Matrix2 Y( R.template block<2,3>(1,0) * L );
451 
452  // Sres= Y R ===> (16 m + 8a)
453  Sres.data_(1) = Y(0,0)*R(0,0) + Y(0,1)*R(0,1);
454  Sres.data_(2) = Y(0,0)*R(1,0) + Y(0,1)*R(1,1);
455  Sres.data_(3) = Y(1,0)*R(0,0) + Y(1,1)*R(0,1);
456  Sres.data_(4) = Y(1,0)*R(1,0) + Y(1,1)*R(1,1);
457  Sres.data_(5) = Y(1,0)*R(2,0) + Y(1,1)*R(2,1);
458 
459  // r=R' v ( 6m + 3a)
460  const Vector3 r(-R(0,0)*data_(4) + R(0,1)*data_(3),
461  -R(1,0)*data_(4) + R(1,1)*data_(3),
462  -R(2,0)*data_(4) + R(2,1)*data_(3));
463 
464  // Sres_11 (3a)
465  Sres.data_(0) = L(0,0) + L(1,1) - Sres.data_(2) - Sres.data_(5);
466 
467  // Sres + D + (Ev)x ( 9a)
468  Sres.data_(0) += data_(5);
469  Sres.data_(1) += r(2); Sres.data_(2)+= data_(5);
470  Sres.data_(3) +=-r(1); Sres.data_(4)+= r(0); Sres.data_(5) += data_(5);
471 
472  return Sres;
473  }
474 
475  protected:
476  Vector6 data_;
477 
478  };
479 
480 } // namespace se3
481 
482 #endif // ifndef __se3_symmetric3_hpp__
483 
static void svx(const Eigen::MatrixBase< Vector3 > &v, const Symmetric3Tpl &S3, const Eigen::MatrixBase< Matrix3 > &M)
Performs the operation .
Definition: symmetric3.hpp:308
Matrix3 vxs(const Eigen::MatrixBase< Vector3 > &v) const
Performs the operation . This operation is equivalent to applying the cross product of v on each colu...
Definition: symmetric3.hpp:291
Matrix32 decomposeltI() const
Computes L for a symmetric matrix A.
Definition: symmetric3.hpp:427
static void vxs(const Eigen::MatrixBase< Vector3 > &v, const Symmetric3Tpl &S3, const Eigen::MatrixBase< Matrix3 > &M)
Performs the operation . This operation is equivalent to applying the cross product of v on each colu...
Definition: symmetric3.hpp:247
Matrix3 svx(const Eigen::MatrixBase< Vector3 > &v) const
Performs the operation .
Definition: symmetric3.hpp:349