pinocchio  1.2.6-7-g6de3e-dirty
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  operator Symmetric3Tpl () const
153  {
154  const Scalar & x = v[0], & y = v[1], & z = v[2];
155  return Symmetric3Tpl(-m*(y*y+z*z),
156  m* x*y,-m*(x*x+z*z),
157  m* x*z,m* y*z,-m*(x*x+y*y));
158  }
159  };
160 
161  friend AlphaSkewSquare operator* (const Scalar & m, const SkewSquare & sk)
162  { return AlphaSkewSquare(m,sk); }
163 
164  Symmetric3Tpl operator- (const AlphaSkewSquare & v) const
165  {
166  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
167  return Symmetric3Tpl(data_[0]+v.m*(y*y+z*z),
168  data_[1]-v.m* x*y, data_[2]+v.m*(x*x+z*z),
169  data_[3]-v.m* x*z, data_[4]-v.m* y*z,
170  data_[5]+v.m*(x*x+y*y));
171  }
172 
173  Symmetric3Tpl& operator-= (const AlphaSkewSquare & v)
174  {
175  const Scalar & x = v.v[0], & y = v.v[1], & z = v.v[2];
176  data_[0]+=v.m*(y*y+z*z);
177  data_[1]-=v.m* x*y; data_[2]+=v.m*(x*x+z*z);
178  data_[3]-=v.m* x*z; data_[4]-=v.m* y*z; data_[5]+=v.m*(x*x+y*y);
179  return *this;
180  }
181 
182  const Vector6 & data () const {return data_;}
183  Vector6 & data () {return data_;}
184 
185  // static Symmetric3Tpl SkewSq( const Vector3 & v )
186  // {
187  // const Scalar & x = v[0], & y = v[1], & z = v[2];
188  // return Symmetric3Tpl(-y*y-z*z,
189  // x*y, -x*x-z*z,
190  // x*z, y*z, -x*x-y*y );
191  // }
192 
193  /* Shoot a positive definite matrix. */
194  static Symmetric3Tpl RandomPositive()
195  {
196  Scalar
197  a = Scalar(std::rand())/RAND_MAX*2.0-1.0,
198  b = Scalar(std::rand())/RAND_MAX*2.0-1.0,
199  c = Scalar(std::rand())/RAND_MAX*2.0-1.0,
200  d = Scalar(std::rand())/RAND_MAX*2.0-1.0,
201  e = Scalar(std::rand())/RAND_MAX*2.0-1.0,
202  f = Scalar(std::rand())/RAND_MAX*2.0-1.0;
203  return Symmetric3Tpl(a*a+b*b+d*d,
204  a*b+b*c+d*e, b*b+c*c+e*e,
205  a*d+b*e+d*f, b*d+c*e+e*f, d*d+e*e+f*f );
206  }
207 
208  Matrix3 matrix() const
209  {
210  Matrix3 res;
211  res(0,0) = data_(0); res(0,1) = data_(1); res(0,2) = data_(3);
212  res(1,0) = data_(1); res(1,1) = data_(2); res(1,2) = data_(4);
213  res(2,0) = data_(3); res(2,1) = data_(4); res(2,2) = data_(5);
214  return res;
215  }
216  operator Matrix3 () const { return matrix(); }
217 
218  Scalar vtiv (const Vector3 & v) const
219  {
220  const Scalar & x = v[0];
221  const Scalar & y = v[1];
222  const Scalar & z = v[2];
223 
224  const Scalar xx = x*x;
225  const Scalar xy = x*y;
226  const Scalar xz = x*z;
227  const Scalar yy = y*y;
228  const Scalar yz = y*z;
229  const Scalar zz = z*z;
230 
231  return data_(0)*xx + data_(2)*yy + data_(5)*zz + 2.*(data_(1)*xy + data_(3)*xz + data_(4)*yz);
232  }
233 
234  Symmetric3Tpl operator+(const Symmetric3Tpl & s2) const
235  {
236  return Symmetric3Tpl((data_+s2.data_).eval());
237  }
238 
239  Symmetric3Tpl & operator+=(const Symmetric3Tpl & s2)
240  {
241  data_ += s2.data_; return *this;
242  }
243 
244  Vector3 operator*(const Vector3 &v) const
245  {
246  return Vector3(
247  data_(0) * v(0) + data_(1) * v(1) + data_(3) * v(2),
248  data_(1) * v(0) + data_(2) * v(1) + data_(4) * v(2),
249  data_(3) * v(0) + data_(4) * v(1) + data_(5) * v(2)
250  );
251  }
252 
253  // Matrix3 operator*(const Matrix3 &a) const
254  // {
255  // Matrix3 r;
256  // for(unsigned int i=0; i<3; ++i)
257  // {
258  // r(0,i) = data_(0) * a(0,i) + data_(1) * a(1,i) + data_(3) * a(2,i);
259  // r(1,i) = data_(1) * a(0,i) + data_(2) * a(1,i) + data_(4) * a(2,i);
260  // r(2,i) = data_(3) * a(0,i) + data_(4) * a(1,i) + data_(5) * a(2,i);
261  // }
262  // return r;
263  // }
264 
265  const Scalar& operator()(const int &i,const int &j) const
266  {
267  return ((i!=2)&&(j!=2)) ? data_[i+j] : data_[i+j+1];
268  }
269 
270  Symmetric3Tpl operator-(const Matrix3 &S) const
271  {
272  assert( (S-S.transpose()).isMuchSmallerThan(S) );
273  return Symmetric3Tpl( data_(0)-S(0,0),
274  data_(1)-S(1,0), data_(2)-S(1,1),
275  data_(3)-S(2,0), data_(4)-S(2,1), data_(5)-S(2,2) );
276  }
277 
278  Symmetric3Tpl operator+(const Matrix3 &S) const
279  {
280  assert( (S-S.transpose()).isMuchSmallerThan(S) );
281  return Symmetric3Tpl( data_(0)+S(0,0),
282  data_(1)+S(1,0), data_(2)+S(1,1),
283  data_(3)+S(2,0), data_(4)+S(2,1), data_(5)+S(2,2) );
284  }
285 
286  /* --- Symmetric R*S*R' and R'*S*R products --- */
287  public: //private:
288 
291  Matrix32 decomposeltI() const
292  {
293  Matrix32 L;
294  L <<
295  data_(0) - data_(5), data_(1),
296  data_(1), data_(2) - data_(5),
297  2*data_(3), data_(4) + data_(4);
298  return L;
299  }
300 
301  /* R*S*R' */
302  template<typename D>
303  Symmetric3Tpl rotate(const Eigen::MatrixBase<D> & R) const
304  {
305  EIGEN_STATIC_ASSERT_MATRIX_SPECIFIC_SIZE(D,3,3);
306  assert( (R.transpose()*R).isApprox(Matrix3::Identity()) );
307 
308  Symmetric3Tpl Sres;
309 
310  // 4 a
311  const Matrix32 L( decomposeltI() );
312 
313  // Y = R' L ===> (12 m + 8 a)
314  const Matrix2 Y( R.template block<2,3>(1,0) * L );
315 
316  // Sres= Y R ===> (16 m + 8a)
317  Sres.data_(1) = Y(0,0)*R(0,0) + Y(0,1)*R(0,1);
318  Sres.data_(2) = Y(0,0)*R(1,0) + Y(0,1)*R(1,1);
319  Sres.data_(3) = Y(1,0)*R(0,0) + Y(1,1)*R(0,1);
320  Sres.data_(4) = Y(1,0)*R(1,0) + Y(1,1)*R(1,1);
321  Sres.data_(5) = Y(1,0)*R(2,0) + Y(1,1)*R(2,1);
322 
323  // r=R' v ( 6m + 3a)
324  const Vector3 r(-R(0,0)*data_(4) + R(0,1)*data_(3),
325  -R(1,0)*data_(4) + R(1,1)*data_(3),
326  -R(2,0)*data_(4) + R(2,1)*data_(3));
327 
328  // Sres_11 (3a)
329  Sres.data_(0) = L(0,0) + L(1,1) - Sres.data_(2) - Sres.data_(5);
330 
331  // Sres + D + (Ev)x ( 9a)
332  Sres.data_(0) += data_(5);
333  Sres.data_(1) += r(2); Sres.data_(2)+= data_(5);
334  Sres.data_(3) +=-r(1); Sres.data_(4)+= r(0); Sres.data_(5) += data_(5);
335 
336  return Sres;
337  }
338 
339  protected:
340  Vector6 data_;
341 
342  };
343 
344 } // namespace se3
345 
346 #endif // ifndef __se3_symmetric3_hpp__
347 
Matrix32 decomposeltI() const
Computes L for a symmetric matrix A.
Definition: symmetric3.hpp:291