pinocchio  2.1.3
spatial/symmetric3.hpp
1 //
2 // Copyright (c) 2014-2017 CNRS
3 //
4 
5 #ifndef __pinocchio_symmetric3_hpp__
6 #define __pinocchio_symmetric3_hpp__
7 
8 #include "pinocchio/macros.hpp"
9 
10 namespace pinocchio
11 {
12 
13  template<typename _Scalar, int _Options>
14  class Symmetric3Tpl
15  {
16  public:
17  typedef _Scalar Scalar;
18  enum { Options = _Options };
19  typedef Eigen::Matrix<Scalar,3,1,Options> Vector3;
20  typedef Eigen::Matrix<Scalar,6,1,Options> Vector6;
21  typedef Eigen::Matrix<Scalar,3,3,Options> Matrix3;
22  typedef Eigen::Matrix<Scalar,2,2,Options> Matrix2;
23  typedef Eigen::Matrix<Scalar,3,2,Options> Matrix32;
24 
25  EIGEN_MAKE_ALIGNED_OPERATOR_NEW
26 
27  public:
28  Symmetric3Tpl(): m_data() {}
29 
30  template<typename Sc,int N,int Opt>
31  explicit Symmetric3Tpl(const Eigen::Matrix<Sc,N,N,Opt> & I)
32  {
33  EIGEN_STATIC_ASSERT(N==3,THIS_METHOD_IS_ONLY_FOR_MATRICES_OF_A_SPECIFIC_SIZE)
34  assert( (I-I.transpose()).isMuchSmallerThan(I) );
35  m_data(0) = I(0,0);
36  m_data(1) = I(1,0); m_data(2) = I(1,1);
37  m_data(3) = I(2,0); m_data(4) = I(2,1); m_data(5) = I(2,2);
38  }
39 
40  explicit Symmetric3Tpl(const Vector6 & I) : m_data(I) {}
41 
42  Symmetric3Tpl(const Scalar & a0, const Scalar & a1, const Scalar & a2,
43  const Scalar & a3, const Scalar & a4, const Scalar & a5)
44  { m_data << a0,a1,a2,a3,a4,a5; }
45 
46  static Symmetric3Tpl Zero() { return Symmetric3Tpl(Vector6::Zero()); }
47  void setZero() { m_data.setZero(); }
48 
49  static Symmetric3Tpl Random() { return RandomPositive(); }
50  void setRandom()
51  {
52  Scalar
53  a = Scalar(std::rand())/RAND_MAX*2.0-1.0,
54  b = Scalar(std::rand())/RAND_MAX*2.0-1.0,
55  c = Scalar(std::rand())/RAND_MAX*2.0-1.0,
56  d = Scalar(std::rand())/RAND_MAX*2.0-1.0,
57  e = Scalar(std::rand())/RAND_MAX*2.0-1.0,
58  f = Scalar(std::rand())/RAND_MAX*2.0-1.0;
59 
60  m_data << a, b, c, d, e, f;
61  }
62 
63  static Symmetric3Tpl Identity() { return Symmetric3Tpl(1, 0, 1, 0, 0, 1); }
64  void setIdentity() { m_data << 1, 0, 1, 0, 0, 1; }
65 
66  /* Requiered by Inertia::operator== */
67  bool operator==(const Symmetric3Tpl & other) const
68  { return m_data == other.m_data; }
69 
70  bool operator!=(const Symmetric3Tpl & other) const
71  { return !(*this == other); }
72 
73  bool isApprox(const Symmetric3Tpl & other,
74  const Scalar & prec = Eigen::NumTraits<Scalar>::dummy_precision()) const
75  { return m_data.isApprox(other.m_data,prec); }
76 
77  void fill(const Scalar value) { m_data.fill(value); }
78 
79  struct SkewSquare
80  {
81  const Vector3 & v;
82  SkewSquare( const Vector3 & v ) : v(v) {}
83  operator Symmetric3Tpl () const
84  {
85  const Scalar & x = v[0], & y = v[1], & z = v[2];
86  return Symmetric3Tpl( -y*y-z*z,
87  x*y , -x*x-z*z,
88  x*z , y*z , -x*x-y*y );
89  }
90  }; // struct SkewSquare
91 
92  Symmetric3Tpl operator- (const SkewSquare & v) const
93  {
94  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
95  return Symmetric3Tpl(m_data[0]+y*y+z*z,
96  m_data[1]-x*y,m_data[2]+x*x+z*z,
97  m_data[3]-x*z,m_data[4]-y*z,m_data[5]+x*x+y*y);
98  }
99 
100  Symmetric3Tpl& operator-= (const SkewSquare & v)
101  {
102  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
103  m_data[0]+=y*y+z*z;
104  m_data[1]-=x*y; m_data[2]+=x*x+z*z;
105  m_data[3]-=x*z; m_data[4]-=y*z; m_data[5]+=x*x+y*y;
106  return *this;
107  }
108 
109  template<typename D>
110  friend Matrix3 operator- (const Symmetric3Tpl & S, const Eigen::MatrixBase <D> & M)
111  {
112  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(D,3,3);
113  Matrix3 result (S.matrix());
114  result -= M;
115 
116  return result;
117  }
118 
120  {
121  const Scalar & m;
122  const Vector3 & v;
123 
124  AlphaSkewSquare(const Scalar & m, const SkewSquare & v) : m(m),v(v.v) {}
125  AlphaSkewSquare(const Scalar & m, const Vector3 & v) : m(m),v(v) {}
126 
127  operator Symmetric3Tpl () const
128  {
129  const Scalar & x = v[0], & y = v[1], & z = v[2];
130  return Symmetric3Tpl(-m*(y*y+z*z),
131  m* x*y,-m*(x*x+z*z),
132  m* x*z,m* y*z,-m*(x*x+y*y));
133  }
134  };
135 
136  friend AlphaSkewSquare operator* (const Scalar & m, const SkewSquare & sk)
137  { return AlphaSkewSquare(m,sk); }
138 
139  Symmetric3Tpl operator- (const AlphaSkewSquare & v) const
140  {
141  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
142  return Symmetric3Tpl(m_data[0]+v.m*(y*y+z*z),
143  m_data[1]-v.m* x*y, m_data[2]+v.m*(x*x+z*z),
144  m_data[3]-v.m* x*z, m_data[4]-v.m* y*z,
145  m_data[5]+v.m*(x*x+y*y));
146  }
147 
148  Symmetric3Tpl& operator-= (const AlphaSkewSquare & v)
149  {
150  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
151  m_data[0]+=v.m*(y*y+z*z);
152  m_data[1]-=v.m* x*y; m_data[2]+=v.m*(x*x+z*z);
153  m_data[3]-=v.m* x*z; m_data[4]-=v.m* y*z; m_data[5]+=v.m*(x*x+y*y);
154  return *this;
155  }
156 
157  const Vector6 & data () const {return m_data;}
158  Vector6 & data () {return m_data;}
159 
160  // static Symmetric3Tpl SkewSq( const Vector3 & v )
161  // {
162  // const Scalar & x = v[0], & y = v[1], & z = v[2];
163  // return Symmetric3Tpl(-y*y-z*z,
164  // x*y, -x*x-z*z,
165  // x*z, y*z, -x*x-y*y );
166  // }
167 
168  /* Shoot a positive definite matrix. */
169  static Symmetric3Tpl RandomPositive()
170  {
171  Scalar
172  a = Scalar(std::rand())/RAND_MAX*2.0-1.0,
173  b = Scalar(std::rand())/RAND_MAX*2.0-1.0,
174  c = Scalar(std::rand())/RAND_MAX*2.0-1.0,
175  d = Scalar(std::rand())/RAND_MAX*2.0-1.0,
176  e = Scalar(std::rand())/RAND_MAX*2.0-1.0,
177  f = Scalar(std::rand())/RAND_MAX*2.0-1.0;
178  return Symmetric3Tpl(a*a+b*b+d*d,
179  a*b+b*c+d*e, b*b+c*c+e*e,
180  a*d+b*e+d*f, b*d+c*e+e*f, d*d+e*e+f*f );
181  }
182 
183  Matrix3 matrix() const
184  {
185  Matrix3 res;
186  res(0,0) = m_data(0); res(0,1) = m_data(1); res(0,2) = m_data(3);
187  res(1,0) = m_data(1); res(1,1) = m_data(2); res(1,2) = m_data(4);
188  res(2,0) = m_data(3); res(2,1) = m_data(4); res(2,2) = m_data(5);
189  return res;
190  }
191  operator Matrix3 () const { return matrix(); }
192 
193  Scalar vtiv (const Vector3 & v) const
194  {
195  const Scalar & x = v[0];
196  const Scalar & y = v[1];
197  const Scalar & z = v[2];
198 
199  const Scalar xx = x*x;
200  const Scalar xy = x*y;
201  const Scalar xz = x*z;
202  const Scalar yy = y*y;
203  const Scalar yz = y*z;
204  const Scalar zz = z*z;
205 
206  return m_data(0)*xx + m_data(2)*yy + m_data(5)*zz + 2.*(m_data(1)*xy + m_data(3)*xz + m_data(4)*yz);
207  }
208 
219  template<typename Vector3, typename Matrix3>
220  static void vxs(const Eigen::MatrixBase<Vector3> & v,
221  const Symmetric3Tpl & S3,
222  const Eigen::MatrixBase<Matrix3> & M)
223  {
224  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3,3);
225  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix3,3,3);
226 
227  const Scalar & a = S3.data()[0];
228  const Scalar & b = S3.data()[1];
229  const Scalar & c = S3.data()[2];
230  const Scalar & d = S3.data()[3];
231  const Scalar & e = S3.data()[4];
232  const Scalar & f = S3.data()[5];
233 
234 
235  const typename Vector3::RealScalar & v0 = v[0];
236  const typename Vector3::RealScalar & v1 = v[1];
237  const typename Vector3::RealScalar & v2 = v[2];
238 
239  Matrix3 & M_ = const_cast<Eigen::MatrixBase<Matrix3> &>(M).derived();
240  M_(0,0) = d * v1 - b * v2;
241  M_(1,0) = a * v2 - d * v0;
242  M_(2,0) = b * v0 - a * v1;
243 
244  M_(0,1) = e * v1 - c * v2;
245  M_(1,1) = b * v2 - e * v0;
246  M_(2,1) = c * v0 - b * v1;
247 
248  M_(0,2) = f * v1 - e * v2;
249  M_(1,2) = d * v2 - f * v0;
250  M_(2,2) = e * v0 - d * v1;
251  }
252 
263  template<typename Vector3>
264  Matrix3 vxs(const Eigen::MatrixBase<Vector3> & v) const
265  {
266  Matrix3 M;
267  vxs(v,*this,M);
268  return M;
269  }
270 
280  template<typename Vector3, typename Matrix3>
281  static void svx(const Eigen::MatrixBase<Vector3> & v,
282  const Symmetric3Tpl & S3,
283  const Eigen::MatrixBase<Matrix3> & M)
284  {
285  EIGEN_STATIC_ASSERT_VECTOR_SPECIFIC_SIZE(Vector3,3);
286  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(Matrix3,3,3);
287 
288  const Scalar & a = S3.data()[0];
289  const Scalar & b = S3.data()[1];
290  const Scalar & c = S3.data()[2];
291  const Scalar & d = S3.data()[3];
292  const Scalar & e = S3.data()[4];
293  const Scalar & f = S3.data()[5];
294 
295  const typename Vector3::RealScalar & v0 = v[0];
296  const typename Vector3::RealScalar & v1 = v[1];
297  const typename Vector3::RealScalar & v2 = v[2];
298 
299  Matrix3 & M_ = const_cast<Eigen::MatrixBase<Matrix3> &>(M).derived();
300  M_(0,0) = b * v2 - d * v1;
301  M_(1,0) = c * v2 - e * v1;
302  M_(2,0) = e * v2 - f * v1;
303 
304  M_(0,1) = d * v0 - a * v2;
305  M_(1,1) = e * v0 - b * v2;
306  M_(2,1) = f * v0 - d * v2;
307 
308  M_(0,2) = a * v1 - b * v0;
309  M_(1,2) = b * v1 - c * v0;
310  M_(2,2) = d * v1 - e * v0;
311  }
312 
321  template<typename Vector3>
322  Matrix3 svx(const Eigen::MatrixBase<Vector3> & v) const
323  {
324  Matrix3 M;
325  svx(v,*this,M);
326  return M;
327  }
328 
329  Symmetric3Tpl operator+(const Symmetric3Tpl & s2) const
330  {
331  return Symmetric3Tpl((m_data+s2.m_data).eval());
332  }
333 
334  Symmetric3Tpl & operator+=(const Symmetric3Tpl & s2)
335  {
336  m_data += s2.m_data; return *this;
337  }
338 
339  template<typename V3in, typename V3out>
340  static void rhsMult(const Symmetric3Tpl & S3,
341  const Eigen::MatrixBase<V3in> & vin,
342  const Eigen::MatrixBase<V3out> & vout)
343  {
344  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(V3in,Vector3);
345  EIGEN_STATIC_ASSERT_SAME_VECTOR_SIZE(V3out,Vector3);
346 
347  Eigen::MatrixBase<V3out> & vout_ = const_cast<Eigen::MatrixBase<V3out>&>(vout);
348 
349  vout_[0] = S3.m_data(0) * vin[0] + S3.m_data(1) * vin[1] + S3.m_data(3) * vin[2];
350  vout_[1] = S3.m_data(1) * vin[0] + S3.m_data(2) * vin[1] + S3.m_data(4) * vin[2];
351  vout_[2] = S3.m_data(3) * vin[0] + S3.m_data(4) * vin[1] + S3.m_data(5) * vin[2];
352  }
353 
354  template<typename V3>
355  Vector3 operator*(const Eigen::MatrixBase<V3> & v) const
356  {
357  Vector3 res;
358  rhsMult(*this,v,res);
359  return res;
360  }
361 
362  // Matrix3 operator*(const Matrix3 &a) const
363  // {
364  // Matrix3 r;
365  // for(unsigned int i=0; i<3; ++i)
366  // {
367  // r(0,i) = m_data(0) * a(0,i) + m_data(1) * a(1,i) + m_data(3) * a(2,i);
368  // r(1,i) = m_data(1) * a(0,i) + m_data(2) * a(1,i) + m_data(4) * a(2,i);
369  // r(2,i) = m_data(3) * a(0,i) + m_data(4) * a(1,i) + m_data(5) * a(2,i);
370  // }
371  // return r;
372  // }
373 
374  const Scalar& operator()(const int &i,const int &j) const
375  {
376  return ((i!=2)&&(j!=2)) ? m_data[i+j] : m_data[i+j+1];
377  }
378 
379  Symmetric3Tpl operator-(const Matrix3 &S) const
380  {
381  assert( (S-S.transpose()).isMuchSmallerThan(S) );
382  return Symmetric3Tpl( m_data(0)-S(0,0),
383  m_data(1)-S(1,0), m_data(2)-S(1,1),
384  m_data(3)-S(2,0), m_data(4)-S(2,1), m_data(5)-S(2,2) );
385  }
386 
387  Symmetric3Tpl operator+(const Matrix3 &S) const
388  {
389  assert( (S-S.transpose()).isMuchSmallerThan(S) );
390  return Symmetric3Tpl( m_data(0)+S(0,0),
391  m_data(1)+S(1,0), m_data(2)+S(1,1),
392  m_data(3)+S(2,0), m_data(4)+S(2,1), m_data(5)+S(2,2) );
393  }
394 
395  /* --- Symmetric R*S*R' and R'*S*R products --- */
396  public: //private:
397 
400  Matrix32 decomposeltI() const
401  {
402  Matrix32 L;
403  L <<
404  m_data(0) - m_data(5), m_data(1),
405  m_data(1), m_data(2) - m_data(5),
406  2*m_data(3), m_data(4) + m_data(4);
407  return L;
408  }
409 
410  /* R*S*R' */
411  template<typename D>
412  Symmetric3Tpl rotate(const Eigen::MatrixBase<D> & R) const
413  {
414  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(D,3,3);
415  assert( (R.transpose()*R).isApprox(Matrix3::Identity()) );
416 
417  Symmetric3Tpl Sres;
418 
419  // 4 a
420  const Matrix32 L( decomposeltI() );
421 
422  // Y = R' L ===> (12 m + 8 a)
423  const Matrix2 Y( R.template block<2,3>(1,0) * L );
424 
425  // Sres= Y R ===> (16 m + 8a)
426  Sres.m_data(1) = Y(0,0)*R(0,0) + Y(0,1)*R(0,1);
427  Sres.m_data(2) = Y(0,0)*R(1,0) + Y(0,1)*R(1,1);
428  Sres.m_data(3) = Y(1,0)*R(0,0) + Y(1,1)*R(0,1);
429  Sres.m_data(4) = Y(1,0)*R(1,0) + Y(1,1)*R(1,1);
430  Sres.m_data(5) = Y(1,0)*R(2,0) + Y(1,1)*R(2,1);
431 
432  // r=R' v ( 6m + 3a)
433  const Vector3 r(-R(0,0)*m_data(4) + R(0,1)*m_data(3),
434  -R(1,0)*m_data(4) + R(1,1)*m_data(3),
435  -R(2,0)*m_data(4) + R(2,1)*m_data(3));
436 
437  // Sres_11 (3a)
438  Sres.m_data(0) = L(0,0) + L(1,1) - Sres.m_data(2) - Sres.m_data(5);
439 
440  // Sres + D + (Ev)x ( 9a)
441  Sres.m_data(0) += m_data(5);
442  Sres.m_data(1) += r(2); Sres.m_data(2)+= m_data(5);
443  Sres.m_data(3) +=-r(1); Sres.m_data(4)+= r(0); Sres.m_data(5) += m_data(5);
444 
445  return Sres;
446  }
447 
449  template<typename NewScalar>
451  {
452  return Symmetric3Tpl<NewScalar,Options>(m_data.template cast<NewScalar>());
453  }
454 
455  protected:
456  Vector6 m_data;
457 
458  };
459 
460 } // namespace pinocchio
461 
462 #endif // ifndef __pinocchio_symmetric3_hpp__
463 
Matrix3 svx(const Eigen::MatrixBase< Vector3 > &v) const
Performs the operation .
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...
Matrix32 decomposeltI() const
Computes L for a symmetric matrix A.
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...
Main pinocchio namespace.
Definition: treeview.dox:24
static void svx(const Eigen::MatrixBase< Vector3 > &v, const Symmetric3Tpl &S3, const Eigen::MatrixBase< Matrix3 > &M)
Performs the operation .
Symmetric3Tpl< NewScalar, Options > cast() const