Main MRPT website > C++ reference for MRPT 1.5.6
eigen_plugins.h
Go to the documentation of this file.
1 /* +---------------------------------------------------------------------------+
2  | Mobile Robot Programming Toolkit (MRPT) |
3  | http://www.mrpt.org/ |
4  | |
5  | Copyright (c) 2005-2017, Individual contributors, see AUTHORS file |
6  | See: http://www.mrpt.org/Authors - All rights reserved. |
7  | Released under BSD License. See details in http://www.mrpt.org/License |
8  +---------------------------------------------------------------------------+ */
9 
10 // -------------------------------------------------------------------------
11 // Note: This file will be included within the body of Eigen::MatrixBase
12 // -------------------------------------------------------------------------
13 public:
14  /** @name MRPT plugin: Types
15  * @{ */
16  // size is constant
17  enum { static_size = RowsAtCompileTime*ColsAtCompileTime };
18  /** @} */
19 
20 
21  /** @name MRPT plugin: Basic iterators. These iterators are intended for 1D matrices only, i.e. column or row vectors.
22  * @{ */
23  typedef Scalar* iterator;
24  typedef const Scalar* const_iterator;
25 
26  EIGEN_STRONG_INLINE iterator begin() { return derived().data(); }
27  EIGEN_STRONG_INLINE iterator end() { return (&(derived().data()[size()-1]))+1; }
28  EIGEN_STRONG_INLINE const_iterator begin() const { return derived().data(); }
29  EIGEN_STRONG_INLINE const_iterator end() const { return (&(derived().data()[size()-1]))+1; }
30 
31  /** @} */
32 
33 
34  /** @name MRPT plugin: Set/get/load/save and other miscelaneous methods
35  * @{ */
36 
37  /*! Fill all the elements with a given value */
38  EIGEN_STRONG_INLINE void fill(const Scalar v) { derived().setConstant(v); }
39 
40  /*! Fill all the elements with a given value */
41  EIGEN_STRONG_INLINE void assign(const Scalar v) { derived().setConstant(v); }
42  /*! Resize to N and set all the elements to a given value */
43  EIGEN_STRONG_INLINE void assign(size_t N, const Scalar v) { derived().resize(N); derived().setConstant(v); }
44 
45  /** Get number of rows */
46  EIGEN_STRONG_INLINE size_t getRowCount() const { return rows(); }
47  /** Get number of columns */
48  EIGEN_STRONG_INLINE size_t getColCount() const { return cols(); }
49 
50  /** Make the matrix an identity matrix (the diagonal values can be 1.0 or any other value) */
51  EIGEN_STRONG_INLINE void unit(const size_t nRows, const Scalar diag_vals) {
52  if (diag_vals==1)
53  derived().setIdentity(nRows,nRows);
54  else {
55  derived().setZero(nRows,nRows);
56  derived().diagonal().setConstant(diag_vals);
57  }
58  }
59 
60  /** Make the matrix an identity matrix */
61  EIGEN_STRONG_INLINE void unit() { derived().setIdentity(); }
62  /** Make the matrix an identity matrix */
63  EIGEN_STRONG_INLINE void eye() { derived().setIdentity(); }
64 
65  /** Set all elements to zero */
66  EIGEN_STRONG_INLINE void zeros() { derived().setZero(); }
67  /** Resize and set all elements to zero */
68  EIGEN_STRONG_INLINE void zeros(const size_t row,const size_t col) { derived().setZero(row,col); }
69 
70  /** Resize matrix and set all elements to one */
71  EIGEN_STRONG_INLINE void ones(const size_t row, const size_t col) { derived().setOnes(row,col); }
72  /** Set all elements to one */
73  EIGEN_STRONG_INLINE void ones() { derived().setOnes(); }
74 
75  /** Fast but unsafe method to obtain a pointer to a given row of the matrix (Use only in time critical applications)
76  * VERY IMPORTANT WARNING: You must be aware of the memory layout, either Column or Row-major ordering.
77  */
78  EIGEN_STRONG_INLINE Scalar * get_unsafe_row(size_t row) { return &derived().coeffRef(row,0); }
79  EIGEN_STRONG_INLINE const Scalar* get_unsafe_row(size_t row) const { return &derived().coeff(row,0); }
80 
81  /** Read-only access to one element (Use with caution, bounds are not checked!) */
82  EIGEN_STRONG_INLINE Scalar get_unsafe(const size_t row, const size_t col) const {
83 #ifdef _DEBUG
84  return derived()(row,col);
85 #else
86  return derived().coeff(row,col);
87 #endif
88  }
89  /** Reference access to one element (Use with caution, bounds are not checked!) */
90  EIGEN_STRONG_INLINE Scalar& get_unsafe(const size_t row, const size_t col) { //-V659
91 #ifdef _DEBUG
92  return derived()(row,col);
93 #else
94  return derived().coeffRef(row,col);
95 #endif
96  }
97  /** Sets an element (Use with caution, bounds are not checked!) */
98  EIGEN_STRONG_INLINE void set_unsafe(const size_t row, const size_t col, const Scalar val) {
99 #ifdef _DEBUG
100  derived()(row,col) = val;
101 #else
102  derived().coeffRef(row,col) = val;
103 #endif
104  }
105 
106  /** Insert an element at the end of the container (for 1D vectors/arrays) */
107  EIGEN_STRONG_INLINE void push_back(Scalar val)
108  {
109  const Index N = size();
110  derived().conservativeResize(N+1);
111  coeffRef(N) = val;
112  }
113 
114  EIGEN_STRONG_INLINE bool isSquare() const { return cols()==rows(); }
115  EIGEN_STRONG_INLINE bool isSingular(const Scalar absThreshold = 0) const { return std::abs(derived().determinant())<absThreshold; }
116 
117  /** Read a matrix from a string in Matlab-like format, for example "[1 0 2; 0 4 -1]"
118  * The string must start with '[' and end with ']'. Rows are separated by semicolons ';' and
119  * columns in each row by one or more whitespaces ' ', commas ',' or tabs '\t'.
120  *
121  * This format is also used for CConfigFile::read_matrix.
122  *
123  * This template method can be instantiated for matrices of the types: int, long, unsinged int, unsigned long, float, double, long double
124  *
125  * \return true on success. false if the string is malformed, and then the matrix will be resized to 0x0.
126  * \sa inMatlabFormat, CConfigFile::read_matrix
127  */
128  bool fromMatlabStringFormat(const std::string &s, std::ostream *dump_errors_here = NULL);
129  // Method implemented in eigen_plugins_impl.h
130 
131  /** Dump matrix in matlab format.
132  * This template method can be instantiated for matrices of the types: int, long, unsinged int, unsigned long, float, double, long double
133  * \sa fromMatlabStringFormat
134  */
135  std::string inMatlabFormat(const size_t decimal_digits = 6) const;
136  // Method implemented in eigen_plugins_impl.h
137 
138  /** Save matrix to a text file, compatible with MATLAB text format (see also the methods of matrix classes themselves).
139  * \param theMatrix It can be a CMatrixTemplate or a CMatrixFixedNumeric.
140  * \param file The target filename.
141  * \param fileFormat See TMatrixTextFileFormat. The format of the numbers in the text file.
142  * \param appendMRPTHeader Insert this header to the file "% File generated by MRPT. Load with MATLAB with: VAR=load(FILENAME);"
143  * \param userHeader Additional text to be written at the head of the file. Typically MALAB comments "% This file blah blah". Final end-of-line is not needed.
144  * \sa loadFromTextFile, CMatrixTemplate::inMatlabFormat, SAVE_MATRIX
145  */
146  void saveToTextFile(
147  const std::string &file,
149  bool appendMRPTHeader = false,
150  const std::string &userHeader = std::string()
151  ) const;
152  // Method implemented in eigen_plugins_impl.h
153 
154  /** Load matrix from a text file, compatible with MATLAB text format.
155  * Lines starting with '%' or '#' are interpreted as comments and ignored.
156  * \sa saveToTextFile, fromMatlabStringFormat
157  */
158  void loadFromTextFile(const std::string &file);
159  // Method implemented in eigen_plugins_impl.h
160 
161  //! \overload
162  void loadFromTextFile(std::istream &_input_text_stream);
163  // Method implemented in eigen_plugins_impl.h
164 
165  EIGEN_STRONG_INLINE void multiplyColumnByScalar(size_t c, Scalar s) { derived().col(c)*=s; }
166  EIGEN_STRONG_INLINE void multiplyRowByScalar(size_t r, Scalar s) { derived().row(r)*=s; }
167 
168  EIGEN_STRONG_INLINE void swapCols(size_t i1,size_t i2) { derived().col(i1).swap( derived().col(i2) ); }
169  EIGEN_STRONG_INLINE void swapRows(size_t i1,size_t i2) { derived().row(i1).swap( derived().row(i2) ); }
170 
171  EIGEN_STRONG_INLINE size_t countNonZero() const { return ((*static_cast<const Derived*>(this))!= 0).count(); }
172 
173  /** [VECTORS OR MATRICES] Finds the maximum value
174  * \exception std::exception On an empty input container
175  */
176  EIGEN_STRONG_INLINE Scalar maximum() const
177  {
178  if (size()==0) throw std::runtime_error("maximum: container is empty");
179  return derived().maxCoeff();
180  }
181 
182  /** [VECTORS OR MATRICES] Finds the minimum value
183  * \sa maximum, minimum_maximum
184  * \exception std::exception On an empty input container */
185  EIGEN_STRONG_INLINE Scalar minimum() const
186  {
187  if (size()==0) throw std::runtime_error("minimum: container is empty");
188  return derived().minCoeff();
189  }
190 
191  /** [VECTORS OR MATRICES] Compute the minimum and maximum of a container at once
192  * \sa maximum, minimum
193  * \exception std::exception On an empty input container */
194  EIGEN_STRONG_INLINE void minimum_maximum(
195  Scalar & out_min,
196  Scalar & out_max) const
197  {
198  out_min = minimum();
199  out_max = maximum();
200  }
201 
202 
203  /** [VECTORS ONLY] Finds the maximum value (and the corresponding zero-based index) from a given container.
204  * \exception std::exception On an empty input vector
205  */
206  EIGEN_STRONG_INLINE Scalar maximum(size_t *maxIndex) const
207  {
208  if (size()==0) throw std::runtime_error("maximum: container is empty");
209  Index idx;
210  const Scalar m = derived().maxCoeff(&idx);
211  if (maxIndex) *maxIndex = idx;
212  return m;
213  }
214 
215  /** [VECTORS OR MATRICES] Finds the maximum value (and the corresponding zero-based index) from a given container.
216  * \exception std::exception On an empty input vector
217  */
218  void find_index_max_value(size_t &u,size_t &v,Scalar &valMax) const
219  {
220  if (cols()==0 || rows()==0) throw std::runtime_error("find_index_max_value: container is empty");
221  Index idx1,idx2;
222  valMax = derived().maxCoeff(&idx1,&idx2);
223  u = idx1; v = idx2;
224  }
225 
226 
227  /** [VECTORS ONLY] Finds the minimum value (and the corresponding zero-based index) from a given container.
228  * \sa maximum, minimum_maximum
229  * \exception std::exception On an empty input vector */
230  EIGEN_STRONG_INLINE Scalar minimum(size_t *minIndex) const
231  {
232  if (size()==0) throw std::runtime_error("minimum: container is empty");
233  Index idx;
234  const Scalar m =derived().minCoeff(&idx);
235  if (minIndex) *minIndex = idx;
236  return m;
237  }
238 
239  /** [VECTORS ONLY] Compute the minimum and maximum of a container at once
240  * \sa maximum, minimum
241  * \exception std::exception On an empty input vector */
242  EIGEN_STRONG_INLINE void minimum_maximum(
243  Scalar & out_min,
244  Scalar & out_max,
245  size_t *minIndex,
246  size_t *maxIndex) const
247  {
248  out_min = minimum(minIndex);
249  out_max = maximum(maxIndex);
250  }
251 
252  /** Compute the norm-infinite of a vector ($f[ ||\mathbf{v}||_\infnty $f]), ie the maximum absolute value of the elements. */
253  EIGEN_STRONG_INLINE Scalar norm_inf() const { return lpNorm<Eigen::Infinity>(); }
254 
255  /** Compute the square norm of a vector/array/matrix (the Euclidean distance to the origin, taking all the elements as a single vector). \sa norm */
256  EIGEN_STRONG_INLINE Scalar squareNorm() const { return squaredNorm(); }
257 
258  /*! Sum all the elements, returning a value of the same type than the container */
259  EIGEN_STRONG_INLINE Scalar sumAll() const { return derived().sum(); }
260 
261  /** Computes the laplacian of this square graph weight matrix.
262  * The laplacian matrix is L = D - W, with D a diagonal matrix with the degree of each node, W the
263  */
264  template<typename OtherDerived>
265  EIGEN_STRONG_INLINE void laplacian(Eigen::MatrixBase <OtherDerived>& ret) const
266  {
267  if (rows()!=cols()) throw std::runtime_error("laplacian: Defined for square matrixes only");
268  const Index N = rows();
269  ret = -(*this);
270  for (Index i=0;i<N;i++)
271  {
272  Scalar deg = 0;
273  for (Index j=0;j<N;j++) deg+= derived().coeff(j,i);
274  ret.coeffRef(i,i)+=deg;
275  }
276  }
277 
278 
279  /** Changes the size of matrix, maintaining its previous content as possible and padding with zeros where applicable.
280  * **WARNING**: MRPT's add-on method \a setSize() pads with zeros, while Eigen's \a resize() does NOT (new elements are undefined).
281  */
282  EIGEN_STRONG_INLINE void setSize(size_t row, size_t col)
283  {
284 #ifdef _DEBUG
285  if ((Derived::RowsAtCompileTime!=Eigen::Dynamic && Derived::RowsAtCompileTime!=int(row)) || (Derived::ColsAtCompileTime!=Eigen::Dynamic && Derived::ColsAtCompileTime!=int(col))) {
286  std::stringstream ss; ss << "setSize: Trying to change a fixed sized matrix from " << rows() << "x" << cols() << " to " << row << "x" << col;
287  throw std::runtime_error(ss.str());
288  }
289 #endif
290  const Index oldCols = cols();
291  const Index oldRows = rows();
292  const int nNewCols = int(col) - int(cols());
293  const int nNewRows = int(row) - int(rows());
294  ::mrpt::math::detail::TAuxResizer<Eigen::MatrixBase<Derived>,SizeAtCompileTime>::internal_resize(*this,row,col);
295  if (nNewCols>0) derived().block(0,oldCols,row,nNewCols).setZero();
296  if (nNewRows>0) derived().block(oldRows,0,nNewRows,col).setZero();
297  }
298 
299  /** Efficiently computes only the biggest eigenvector of the matrix using the Power Method, and returns it in the passed vector "x". */
300  template <class OUTVECT>
302  OUTVECT &x,
303  Scalar resolution = Scalar(0.01),
304  size_t maxIterations = 6,
305  int *out_Iterations = NULL,
306  float *out_estimatedResolution = NULL ) const
307  {
308  // Apply the iterative Power Method:
309  size_t iter=0;
310  const Index n = rows();
311  x.resize(n);
312  x.setConstant(1); // Initially, set to all ones, for example...
313  Scalar dif;
314  do // Iterative loop:
315  {
316  Eigen::Matrix<Scalar,Derived::RowsAtCompileTime,1> xx = (*this) * x;
317  xx *= Scalar(1.0/xx.norm());
318  dif = (x-xx).array().abs().sum(); // Compute diference between iterations:
319  x = xx; // Set as current estimation:
320  iter++; // Iteration counter:
321  } while (iter<maxIterations && dif>resolution);
322  if (out_Iterations) *out_Iterations=static_cast<int>(iter);
323  if (out_estimatedResolution) *out_estimatedResolution=dif;
324  }
325 
326  /** Combined matrix power and assignment operator */
327  MatrixBase<Derived>& operator ^= (const unsigned int pow)
328  {
329  if (pow==0)
330  derived().setIdentity();
331  else
332  for (unsigned int i=1;i<pow;i++)
333  derived() *= derived();
334  return *this;
335  }
336 
337  /** Scalar power of all elements to a given power, this is diferent of ^ operator. */
338  EIGEN_STRONG_INLINE void scalarPow(const Scalar s) { (*this)=array().pow(s); }
339 
340  /** Checks for matrix type */
341  EIGEN_STRONG_INLINE bool isDiagonal() const
342  {
343  for (Index c=0;c<cols();c++)
344  for (Index r=0;r<rows();r++)
345  if (r!=c && coeff(r,c)!=0)
346  return false;
347  return true;
348  }
349 
350  /** Finds the maximum value in the diagonal of the matrix. */
351  EIGEN_STRONG_INLINE Scalar maximumDiagonal() const { return diagonal().maxCoeff(); }
352 
353  /** Computes the mean of the entire matrix
354  * \sa meanAndStdAll */
355  EIGEN_STRONG_INLINE double mean() const
356  {
357  if ( size()==0) throw std::runtime_error("mean: Empty container.");
358  return derived().sum()/static_cast<double>(size());
359  }
360 
361  /** Computes a row with the mean values of each column in the matrix and the associated vector with the standard deviation of each column.
362  * \sa mean,meanAndStdAll \exception std::exception If the matrix/vector is empty.
363  * \param unbiased_variance Standard deviation is sum(vals-mean)/K, with K=N-1 or N for unbiased_variance=true or false, respectively.
364  */
365  template <class VEC>
367  VEC &outMeanVector,
368  VEC &outStdVector,
369  const bool unbiased_variance = true ) const
370  {
371  const size_t N = rows();
372  if (N==0) throw std::runtime_error("meanAndStd: Empty container.");
373  const double N_inv = 1.0/N;
374  const double N_ = unbiased_variance ? (N>1 ? 1.0/(N-1) : 1.0) : 1.0/N;
375  outMeanVector.resize(cols());
376  outStdVector.resize(cols());
377  for (Index i=0;i<cols();i++)
378  {
379  outMeanVector[i]= this->col(i).array().sum() * N_inv;
380  outStdVector[i] = std::sqrt( (this->col(i).array()-outMeanVector[i]).square().sum() * N_ );
381  }
382  }
383 
384  /** Computes the mean and standard deviation of all the elements in the matrix as a whole.
385  * \sa mean,meanAndStd \exception std::exception If the matrix/vector is empty.
386  * \param unbiased_variance Standard deviation is sum(vals-mean)/K, with K=N-1 or N for unbiased_variance=true or false, respectively.
387  */
389  double &outMean,
390  double &outStd,
391  const bool unbiased_variance = true ) const
392  {
393  const size_t N = size();
394  if (N==0) throw std::runtime_error("meanAndStdAll: Empty container.");
395  const double N_ = unbiased_variance ? (N>1 ? 1.0/(N-1) : 1.0) : 1.0/N;
396  outMean = derived().array().sum()/static_cast<double>(size());
397  outStd = std::sqrt( (this->array() - outMean).square().sum()*N_);
398  }
399 
400  /** Insert matrix "m" into this matrix at indices (r,c), that is, (*this)(r,c)=m(0,0) and so on */
401  template <typename MAT>
402  EIGEN_STRONG_INLINE void insertMatrix(size_t r,size_t c, const MAT &m) { derived().block(r,c,m.rows(),m.cols())=m; }
403 
404  template <typename MAT>
405  EIGEN_STRONG_INLINE void insertMatrixTranspose(size_t r,size_t c, const MAT &m) { derived().block(r,c,m.cols(),m.rows())=m.adjoint(); }
406 
407  template <typename MAT> EIGEN_STRONG_INLINE void insertRow(size_t nRow, const MAT & aRow) { this->row(nRow) = aRow; }
408  template <typename MAT> EIGEN_STRONG_INLINE void insertCol(size_t nCol, const MAT & aCol) { this->col(nCol) = aCol; }
409 
410  template <typename R> void insertRow(size_t nRow, const std::vector<R> & aRow)
411  {
412  if (static_cast<Index>(aRow.size())!=cols()) throw std::runtime_error("insertRow: Row size doesn't fit the size of this matrix.");
413  for (Index j=0;j<cols();j++)
414  coeffRef(nRow,j) = aRow[j];
415  }
416  template <typename R> void insertCol(size_t nCol, const std::vector<R> & aCol)
417  {
418  if (static_cast<Index>(aCol.size())!=rows()) throw std::runtime_error("insertRow: Row size doesn't fit the size of this matrix.");
419  for (Index j=0;j<rows();j++)
420  coeffRef(j,nCol) = aCol[j];
421  }
422 
423  /** Remove columns of the matrix.*/
424  EIGEN_STRONG_INLINE void removeColumns(const std::vector<size_t> &idxsToRemove)
425  {
426  std::vector<size_t> idxs = idxsToRemove;
427  std::sort( idxs.begin(), idxs.end() );
428  std::vector<size_t>::iterator itEnd = std::unique( idxs.begin(), idxs.end() );
429  idxs.resize( itEnd - idxs.begin() );
430 
431  unsafeRemoveColumns( idxs );
432  }
433 
434  /** Remove columns of the matrix. The unsafe version assumes that, the indices are sorted in ascending order. */
435  EIGEN_STRONG_INLINE void unsafeRemoveColumns(const std::vector<size_t> &idxs)
436  {
437  size_t k = 1;
438  for (std::vector<size_t>::const_reverse_iterator it = idxs.rbegin(); it != idxs.rend(); ++it, ++k)
439  {
440  const size_t nC = cols() - *it - k;
441  if( nC > 0 )
442  derived().block(0,*it,rows(),nC) = derived().block(0,*it+1,rows(),nC).eval();
443  }
444  derived().conservativeResize(NoChange,cols()-idxs.size());
445  }
446 
447  /** Remove rows of the matrix. */
448  EIGEN_STRONG_INLINE void removeRows(const std::vector<size_t> &idxsToRemove)
449  {
450  std::vector<size_t> idxs = idxsToRemove;
451  std::sort( idxs.begin(), idxs.end() );
452  std::vector<size_t>::iterator itEnd = std::unique( idxs.begin(), idxs.end() );
453  idxs.resize( itEnd - idxs.begin() );
454 
455  unsafeRemoveRows( idxs );
456  }
457 
458  /** Remove rows of the matrix. The unsafe version assumes that, the indices are sorted in ascending order. */
459  EIGEN_STRONG_INLINE void unsafeRemoveRows(const std::vector<size_t> &idxs)
460  {
461  size_t k = 1;
462  for (std::vector<size_t>::reverse_iterator it = idxs.rbegin(); it != idxs.rend(); ++it, ++k)
463  {
464  const size_t nR = rows() - *it - k;
465  if( nR > 0 )
466  derived().block(*it,0,nR,cols()) = derived().block(*it+1,0,nR,cols()).eval();
467  }
468  derived().conservativeResize(rows()-idxs.size(),NoChange);
469  }
470 
471  /** Transpose */
472  EIGEN_STRONG_INLINE const AdjointReturnType t() const { return derived().adjoint(); }
473 
474  EIGEN_STRONG_INLINE PlainObject inv() const { PlainObject outMat = derived().inverse().eval(); return outMat; }
475  template <class MATRIX> EIGEN_STRONG_INLINE void inv(MATRIX &outMat) const { outMat = derived().inverse().eval(); }
476  template <class MATRIX> EIGEN_STRONG_INLINE void inv_fast(MATRIX &outMat) const { outMat = derived().inverse().eval(); }
477  EIGEN_STRONG_INLINE Scalar det() const { return derived().determinant(); }
478 
479  /** @} */ // end miscelaneous
480 
481 
482  /** @name MRPT plugin: Multiply and extra addition functions
483  @{ */
484 
485  EIGEN_STRONG_INLINE bool empty() const { return this->getColCount()==0 || this->getRowCount()==0; }
486 
487  /*! Add c (scalar) times A to this matrix: this += A * c */
488  template<typename OTHERMATRIX> EIGEN_STRONG_INLINE void add_Ac(const OTHERMATRIX &m,const Scalar c) { (*this)+=c*m; }
489  /*! Substract c (scalar) times A to this matrix: this -= A * c */
490  template<typename OTHERMATRIX> EIGEN_STRONG_INLINE void substract_Ac(const OTHERMATRIX &m,const Scalar c) { (*this) -= c*m; }
491 
492  /*! Substract A transposed to this matrix: this -= A.adjoint() */
493  template<typename OTHERMATRIX> EIGEN_STRONG_INLINE void substract_At(const OTHERMATRIX &m) { (*this) -= m.adjoint(); }
494 
495  /*! Substract n (integer) times A to this matrix: this -= A * n */
496  template<typename OTHERMATRIX> EIGEN_STRONG_INLINE void substract_An(const OTHERMATRIX& m, const size_t n) { this->noalias() -= n * m; }
497 
498  /*! this += A + A<sup>T</sup> */
499  template<typename OTHERMATRIX> EIGEN_STRONG_INLINE void add_AAt(const OTHERMATRIX &A) { this->noalias() += A; this->noalias() += A.adjoint(); }
500 
501  /*! this -= A + A<sup>T</sup> */ \
502  template<typename OTHERMATRIX> EIGEN_STRONG_INLINE void substract_AAt(const OTHERMATRIX &A) { this->noalias() -= A; this->noalias() -= A.adjoint(); }
503 
504 
505  template <class MATRIX1,class MATRIX2> EIGEN_STRONG_INLINE void multiply( const MATRIX1& A, const MATRIX2 &B ) /*!< this = A * B */ { (*this)= A*B; }
506 
507  template <class MATRIX1,class MATRIX2>
508  EIGEN_STRONG_INLINE void multiply_AB( const MATRIX1& A, const MATRIX2 &B ) /*!< this = A * B */ {
509  (*this)= A*B;
510  }
511 
512  template <typename MATRIX1,typename MATRIX2>
513  EIGEN_STRONG_INLINE void multiply_AtB(const MATRIX1 &A,const MATRIX2 &B) /*!< this=A^t * B */ {
514  *this = A.adjoint() * B;
515  }
516 
517  /*! Computes the vector vOut = this * vIn, where "vIn" is a column vector of the appropriate length. */
518  template<typename OTHERVECTOR1,typename OTHERVECTOR2>
519  EIGEN_STRONG_INLINE void multiply_Ab(const OTHERVECTOR1 &vIn,OTHERVECTOR2 &vOut,bool accumToOutput = false) const {
520  if (accumToOutput) vOut.noalias() += (*this) * vIn;
521  else vOut = (*this) * vIn;
522  }
523 
524  /*! Computes the vector vOut = this<sup>T</sup> * vIn, where "vIn" is a column vector of the appropriate length. */ \
525  template<typename OTHERVECTOR1,typename OTHERVECTOR2>
526  EIGEN_STRONG_INLINE void multiply_Atb(const OTHERVECTOR1 &vIn,OTHERVECTOR2 &vOut,bool accumToOutput = false) const {
527  if (accumToOutput) vOut.noalias() += this->adjoint() * vIn;
528  else vOut = this->adjoint() * vIn;
529  }
530 
531  template <typename MAT_C, typename MAT_R>
532  EIGEN_STRONG_INLINE void multiply_HCHt(const MAT_C &C,MAT_R &R,bool accumResultInOutput=false) const /*!< R = this * C * this<sup>T</sup> */ {
533  if (accumResultInOutput)
534  R.noalias() += (*this) * C * this->adjoint();
535  else R.noalias() = (*this) * C * this->adjoint();
536  }
537 
538  template <typename MAT_C, typename MAT_R>
539  EIGEN_STRONG_INLINE void multiply_HtCH(const MAT_C &C,MAT_R &R,bool accumResultInOutput=false) const /*!< R = this<sup>T</sup> * C * this */ {
540  if (accumResultInOutput)
541  R.noalias() += this->adjoint() * C * (*this);
542  else R.noalias() = this->adjoint() * C * (*this);
543  }
544 
545  /*! R = H * C * H<sup>T</sup> (with a vector H and a symmetric matrix C) In fact when H is a vector, multiply_HCHt_scalar and multiply_HtCH_scalar are exactly equivalent */
546  template <typename MAT_C>
547  EIGEN_STRONG_INLINE Scalar multiply_HCHt_scalar(const MAT_C &C) const {
548  return ( (*this) * C * this->adjoint() ).eval()(0,0);
549  }
550 
551  /*! R = H<sup>T</sup> * C * H (with a vector H and a symmetric matrix C) In fact when H is a vector, multiply_HCHt_scalar and multiply_HtCH_scalar are exactly equivalent */
552  template <typename MAT_C>
553  EIGEN_STRONG_INLINE Scalar multiply_HtCH_scalar(const MAT_C &C) const {
554  return ( this->adjoint() * C * (*this) ).eval()(0,0);
555  }
556 
557  /*! this = C * C<sup>T</sup> * f (with a matrix C and a scalar f). */
558  template<typename MAT_A>
559  EIGEN_STRONG_INLINE void multiply_AAt_scalar(const MAT_A &A,typename MAT_A::Scalar f) {
560  *this = (A * A.adjoint()) * f;
561  }
562 
563  /*! this = C<sup>T</sup> * C * f (with a matrix C and a scalar f). */
564  template<typename MAT_A> EIGEN_STRONG_INLINE void multiply_AtA_scalar(const MAT_A &A,typename MAT_A::Scalar f) {
565  *this = (A.adjoint() * A) * f;
566  }
567 
568  /*! this = A * skew(v), with \a v being a 3-vector (or 3-array) and skew(v) the skew symmetric matrix of v (see mrpt::math::skew_symmetric3) */
569  template <class MAT_A,class SKEW_3VECTOR> void multiply_A_skew3(const MAT_A &A,const SKEW_3VECTOR &v) {
570  mrpt::math::multiply_A_skew3(A,v,*this); }
571 
572  /*! this = skew(v)*A, with \a v being a 3-vector (or 3-array) and skew(v) the skew symmetric matrix of v (see mrpt::math::skew_symmetric3) */
573  template <class SKEW_3VECTOR,class MAT_A> void multiply_skew3_A(const SKEW_3VECTOR &v,const MAT_A &A) {
574  mrpt::math::multiply_skew3_A(v,A,*this); }
575 
576  /** outResult = this * A
577  */
578  template <class MAT_A,class MAT_OUT>
579  EIGEN_STRONG_INLINE void multiply_subMatrix(const MAT_A &A,MAT_OUT &outResult,const size_t A_cols_offset,const size_t A_rows_offset,const size_t A_col_count) const {
580  outResult = derived() * A.block(A_rows_offset,A_cols_offset,derived().cols(),A_col_count);
581  }
582 
583  template <class MAT_A,class MAT_B,class MAT_C>
584  void multiply_ABC(const MAT_A &A, const MAT_B &B, const MAT_C &C) /*!< this = A*B*C */ {
585  *this = A*B*C;
586  }
587 
588  template <class MAT_A,class MAT_B,class MAT_C>
589  void multiply_ABCt(const MAT_A &A, const MAT_B &B, const MAT_C &C) /*!< this = A*B*(C<sup>T</sup>) */ {
590  *this = A*B*C.adjoint();
591  }
592 
593  template <class MAT_A,class MAT_B,class MAT_C>
594  void multiply_AtBC(const MAT_A &A, const MAT_B &B, const MAT_C &C) /*!< this = A(<sup>T</sup>)*B*C */ {
595  *this = A.adjoint()*B*C;
596  }
597 
598  template <class MAT_A,class MAT_B>
599  EIGEN_STRONG_INLINE void multiply_ABt(const MAT_A &A,const MAT_B &B) /*!< this = A * B<sup>T</sup> */ {
600  *this = A*B.adjoint();
601  }
602 
603  template <class MAT_A>
604  EIGEN_STRONG_INLINE void multiply_AAt(const MAT_A &A) /*!< this = A * A<sup>T</sup> */ {
605  *this = A*A.adjoint();
606  }
607 
608  template <class MAT_A>
609  EIGEN_STRONG_INLINE void multiply_AtA(const MAT_A &A) /*!< this = A<sup>T</sup> * A */ {
610  *this = A.adjoint()*A;
611  }
612 
613  template <class MAT_A,class MAT_B>
614  EIGEN_STRONG_INLINE void multiply_result_is_symmetric(const MAT_A &A,const MAT_B &B) /*!< this = A * B (result is symmetric) */ {
615  *this = A*B;
616  }
617 
618  /** @} */ // end multiply functions
619 
620 
621  /** @name MRPT plugin: Eigenvalue / Eigenvectors
622  @{ */
623 
624  /** [For square matrices only] Compute the eigenvectors and eigenvalues (sorted), both returned as matrices: eigenvectors are the columns in "eVecs", and eigenvalues in ascending order as the diagonal of "eVals".
625  * \note Warning: Only the real part of complex eigenvectors and eigenvalues are returned.
626  * \sa eigenVectorsSymmetric, eigenVectorsVec
627  * \return false on error
628  */
629  template <class MATRIX1,class MATRIX2>
630  EIGEN_STRONG_INLINE bool eigenVectors( MATRIX1 & eVecs, MATRIX2 & eVals ) const;
631  // Implemented in eigen_plugins_impl.h (can't be here since Eigen::SelfAdjointEigenSolver isn't defined yet at this point.
632 
633  /** [For square matrices only] Compute the eigenvectors and eigenvalues (sorted), eigenvectors are the columns in "eVecs", and eigenvalues are returned in in ascending order in the vector "eVals".
634  * \note Warning: Only the real part of complex eigenvectors and eigenvalues are returned.
635  * \sa eigenVectorsSymmetric, eigenVectorsVec
636  * \return false on error
637  */
638  template <class MATRIX1,class VECTOR1>
639  EIGEN_STRONG_INLINE bool eigenVectorsVec( MATRIX1 & eVecs, VECTOR1 & eVals ) const;
640  // Implemented in eigen_plugins_impl.h
641 
642  /** [For square matrices only] Compute the eigenvectors and eigenvalues (sorted), and return only the eigenvalues in the vector "eVals".
643  * \note Warning: Only the real part of complex eigenvectors and eigenvalues are returned.
644  * \sa eigenVectorsSymmetric, eigenVectorsVec
645  */
646  template <class VECTOR>
647  EIGEN_STRONG_INLINE void eigenValues( VECTOR & eVals ) const
648  {
649  PlainObject eVecs;
650  eVecs.resizeLike(*this);
651  this->eigenVectorsVec(eVecs,eVals);
652  }
653 
654  /** [For symmetric matrices only] Compute the eigenvectors and eigenvalues (in no particular order), both returned as matrices: eigenvectors are the columns, and eigenvalues \sa eigenVectors
655  */
656  template <class MATRIX1,class MATRIX2>
657  EIGEN_STRONG_INLINE void eigenVectorsSymmetric( MATRIX1 & eVecs, MATRIX2 & eVals ) const;
658  // Implemented in eigen_plugins_impl.h (can't be here since Eigen::SelfAdjointEigenSolver isn't defined yet at this point.
659 
660  /** [For symmetric matrices only] Compute the eigenvectors and eigenvalues (in no particular order), both returned as matrices: eigenvectors are the columns, and eigenvalues \sa eigenVectorsVec
661  */
662  template <class MATRIX1,class VECTOR1>
663  EIGEN_STRONG_INLINE void eigenVectorsSymmetricVec( MATRIX1 & eVecs, VECTOR1 & eVals ) const;
664  // Implemented in eigen_plugins_impl.h
665 
666 
667  /** @} */ // end eigenvalues
668 
669 
670 
671  /** @name MRPT plugin: Linear algebra & decomposition-based methods
672  @{ */
673 
674  /** Cholesky M=U<sup>T</sup> * U decomposition for simetric matrix (upper-half of the matrix will be actually ignored) */
675  template <class MATRIX> EIGEN_STRONG_INLINE bool chol(MATRIX &U) const
676  {
677  Eigen::LLT<PlainObject> Chol = derived().template selfadjointView<Eigen::Lower>().llt();
678  if (Chol.info()==Eigen::NoConvergence)
679  return false;
680  U = PlainObject(Chol.matrixU());
681  return true;
682  }
683 
684  /** Gets the rank of the matrix via the Eigen::ColPivHouseholderQR method
685  * \param threshold If set to >0, it's used as threshold instead of Eigen's default one.
686  */
687  EIGEN_STRONG_INLINE size_t rank(double threshold=0) const
688  {
689  Eigen::ColPivHouseholderQR<PlainObject> QR = this->colPivHouseholderQr();
690  if (threshold>0) QR.setThreshold(threshold);
691  return QR.rank();
692  }
693 
694  /** @} */ // end linear algebra
695 
696 
697 
698  /** @name MRPT plugin: Scalar and element-wise extra operators
699  @{ */
700 
701  /** Scales all elements such as the minimum & maximum values are shifted to the given values */
702  void normalize(Scalar valMin, Scalar valMax)
703  {
704  if (size()==0) return;
705  Scalar curMin,curMax;
706  minimum_maximum(curMin,curMax);
707  Scalar minMaxDelta = curMax - curMin;
708  if (minMaxDelta==0) minMaxDelta = 1;
709  const Scalar minMaxDelta_ = (valMax-valMin)/minMaxDelta;
710  this->array() = (this->array()-curMin)*minMaxDelta_+valMin;
711  }
712  //! \overload
713  inline void adjustRange(Scalar valMin, Scalar valMax) { normalize(valMin,valMax); }
714 
715  /** @} */ // end Scalar
716 
717 
718  /** Extract one row from the matrix into a row vector */
719  template <class OtherDerived> EIGEN_STRONG_INLINE void extractRow(size_t nRow, Eigen::EigenBase<OtherDerived> &v, size_t startingCol = 0) const {
720  v = derived().block(nRow,startingCol,1,cols()-startingCol);
721  }
722  //! \overload
723  template <typename T> inline void extractRow(size_t nRow, std::vector<T> &v, size_t startingCol = 0) const {
724  const size_t N = cols()-startingCol;
725  v.resize(N);
726  for (size_t i=0;i<N;i++) v[i]=(*this)(nRow,startingCol+i);
727  }
728  /** Extract one row from the matrix into a column vector */
729  template <class VECTOR> EIGEN_STRONG_INLINE void extractRowAsCol(size_t nRow, VECTOR &v, size_t startingCol = 0) const
730  {
731  v = derived().adjoint().block(startingCol,nRow,cols()-startingCol,1);
732  }
733 
734 
735  /** Extract one column from the matrix into a column vector */
736  template <class VECTOR> EIGEN_STRONG_INLINE void extractCol(size_t nCol, VECTOR &v, size_t startingRow = 0) const {
737  v = derived().block(startingRow,nCol,rows()-startingRow,1);
738  }
739  //! \overload
740  template <typename T> inline void extractCol(size_t nCol, std::vector<T> &v, size_t startingRow = 0) const {
741  const size_t N = rows()-startingRow;
742  v.resize(N);
743  for (size_t i=0;i<N;i++) v[i]=(*this)(startingRow+i,nCol);
744  }
745 
746  template <class MATRIX> EIGEN_STRONG_INLINE void extractMatrix(const size_t firstRow, const size_t firstCol, MATRIX &m) const
747  {
748  m = derived().block(firstRow,firstCol,m.rows(),m.cols());
749  }
750  template <class MATRIX> EIGEN_STRONG_INLINE void extractMatrix(const size_t firstRow, const size_t firstCol, const size_t nRows, const size_t nCols, MATRIX &m) const
751  {
752  m.resize(nRows,nCols);
753  m = derived().block(firstRow,firstCol,nRows,nCols);
754  }
755 
756  /** Get a submatrix, given its bounds: first & last column and row (inclusive). */
757  template <class MATRIX>
758  EIGEN_STRONG_INLINE void extractSubmatrix(const size_t row_first,const size_t row_last,const size_t col_first,const size_t col_last,MATRIX &out) const
759  {
760  out.resize(row_last-row_first+1,col_last-col_first+1);
761  out = derived().block(row_first,col_first,row_last-row_first+1,col_last-col_first+1);
762  }
763 
764  /** Get a submatrix from a square matrix, by collecting the elements M(idxs,idxs), where idxs is a sequence {block_indices(i):block_indices(i)+block_size-1} for all "i" up to the size of block_indices.
765  * A perfect application of this method is in extracting covariance matrices of a subset of variables from the full covariance matrix.
766  * \sa extractSubmatrix, extractSubmatrixSymmetrical
767  */
768  template <class MATRIX>
770  const size_t block_size,
771  const std::vector<size_t> &block_indices,
772  MATRIX& out) const
773  {
774  if (block_size<1) throw std::runtime_error("extractSubmatrixSymmetricalBlocks: block_size must be >=1");
775  if (cols()!=rows()) throw std::runtime_error("extractSubmatrixSymmetricalBlocks: Matrix is not square.");
776 
777  const size_t N = block_indices.size();
778  const size_t nrows_out=N*block_size;
779  out.resize(nrows_out,nrows_out);
780  if (!N) return; // Done
781  for (size_t dst_row_blk=0;dst_row_blk<N; ++dst_row_blk )
782  {
783  for (size_t dst_col_blk=0;dst_col_blk<N; ++dst_col_blk )
784  {
785 #if defined(_DEBUG)
786  if (block_indices[dst_col_blk]*block_size + block_size-1>=size_t(cols())) throw std::runtime_error("extractSubmatrixSymmetricalBlocks: Indices out of range!");
787 #endif
788  out.block(dst_row_blk * block_size,dst_col_blk * block_size, block_size,block_size)
789  =
790  derived().block(block_indices[dst_row_blk] * block_size, block_indices[dst_col_blk] * block_size, block_size,block_size);
791  }
792  }
793  }
794 
795 
796  /** Get a submatrix from a square matrix, by collecting the elements M(idxs,idxs), where idxs is the sequence of indices passed as argument.
797  * A perfect application of this method is in extracting covariance matrices of a subset of variables from the full covariance matrix.
798  * \sa extractSubmatrix, extractSubmatrixSymmetricalBlocks
799  */
800  template <class MATRIX>
802  const std::vector<size_t> &indices,
803  MATRIX& out) const
804  {
805  if (cols()!=rows()) throw std::runtime_error("extractSubmatrixSymmetrical: Matrix is not square.");
806 
807  const size_t N = indices.size();
808  const size_t nrows_out=N;
809  out.resize(nrows_out,nrows_out);
810  if (!N) return; // Done
811  for (size_t dst_row_blk=0;dst_row_blk<N; ++dst_row_blk )
812  for (size_t dst_col_blk=0;dst_col_blk<N; ++dst_col_blk )
813  out.coeffRef(dst_row_blk,dst_col_blk) = this->coeff(indices[dst_row_blk],indices[dst_col_blk]);
814  }
815 
EIGEN_STRONG_INLINE void set_unsafe(const size_t row, const size_t col, const Scalar val)
Sets an element (Use with caution, bounds are not checked!)
Definition: eigen_plugins.h:98
EIGEN_STRONG_INLINE Scalar det() const
void meanAndStd(VEC &outMeanVector, VEC &outStdVector, const bool unbiased_variance=true) const
Computes a row with the mean values of each column in the matrix and the associated vector with the s...
void adjustRange(Scalar valMin, Scalar valMax)
EIGEN_STRONG_INLINE void extractRowAsCol(size_t nRow, VECTOR &v, size_t startingCol=0) const
Extract one row from the matrix into a column vector.
EIGEN_STRONG_INLINE void scalarPow(const Scalar s)
Scalar power of all elements to a given power, this is diferent of ^ operator.
EIGEN_STRONG_INLINE void multiply_Ab(const OTHERVECTOR1 &vIn, OTHERVECTOR2 &vOut, bool accumToOutput=false) const
EIGEN_STRONG_INLINE bool empty() const
EIGEN_STRONG_INLINE void removeRows(const std::vector< size_t > &idxsToRemove)
Remove rows of the matrix.
EIGEN_STRONG_INLINE Scalar multiply_HCHt_scalar(const MAT_C &C) const
std::string inMatlabFormat(const size_t decimal_digits=6) const
Dump matrix in matlab format.
EIGEN_STRONG_INLINE void fill(const Scalar v)
Definition: eigen_plugins.h:38
EIGEN_STRONG_INLINE bool isSquare() const
engineering format &#39;e&#39;
Definition: math_frwds.h:65
EIGEN_STRONG_INLINE bool eigenVectors(MATRIX1 &eVecs, MATRIX2 &eVals) const
[For square matrices only] Compute the eigenvectors and eigenvalues (sorted), both returned as matric...
EIGEN_STRONG_INLINE void add_AAt(const OTHERMATRIX &A)
EIGEN_STRONG_INLINE bool chol(MATRIX &U) const
Cholesky M=UT * U decomposition for simetric matrix (upper-half of the matrix will be actually ignore...
EIGEN_STRONG_INLINE size_t rank(double threshold=0) const
Gets the rank of the matrix via the Eigen::ColPivHouseholderQR method.
EIGEN_STRONG_INLINE void multiply(const MATRIX1 &A, const MATRIX2 &B)
void extractSubmatrixSymmetricalBlocks(const size_t block_size, const std::vector< size_t > &block_indices, MATRIX &out) const
Get a submatrix from a square matrix, by collecting the elements M(idxs,idxs), where idxs is a sequen...
void multiply_skew3_A(const SKEW_3VECTOR &v, const MAT_A &A)
EIGEN_STRONG_INLINE Scalar * get_unsafe_row(size_t row)
Fast but unsafe method to obtain a pointer to a given row of the matrix (Use only in time critical ap...
Definition: eigen_plugins.h:78
EIGEN_STRONG_INLINE Scalar maximumDiagonal() const
Finds the maximum value in the diagonal of the matrix.
EIGEN_STRONG_INLINE void substract_Ac(const OTHERMATRIX &m, const Scalar c)
GLenum GLsizei n
Definition: glext.h:4618
Scalar * iterator
Definition: eigen_plugins.h:23
void multiply_ABCt(const MAT_A &A, const MAT_B &B, const MAT_C &C)
void multiply_A_skew3(const MAT_A &A, const SKEW_3VECTOR &v, MAT_OUT &out)
Only for vectors/arrays "v" of length3, compute out = A * Skew(v), where Skew(v) is the skew symmetri...
Definition: ops_matrices.h:150
EIGEN_STRONG_INLINE void multiply_HCHt(const MAT_C &C, MAT_R &R, bool accumResultInOutput=false) const
EIGEN_STRONG_INLINE void swapRows(size_t i1, size_t i2)
GLuint GLuint GLsizei GLenum const GLvoid * indices
Definition: glext.h:3512
EIGEN_STRONG_INLINE iterator begin()
Definition: eigen_plugins.h:26
EIGEN_STRONG_INLINE void push_back(Scalar val)
Insert an element at the end of the container (for 1D vectors/arrays)
EIGEN_STRONG_INLINE bool isSingular(const Scalar absThreshold=0) const
EIGEN_STRONG_INLINE void extractCol(size_t nCol, VECTOR &v, size_t startingRow=0) const
Extract one column from the matrix into a column vector.
void find_index_max_value(size_t &u, size_t &v, Scalar &valMax) const
[VECTORS OR MATRICES] Finds the maximum value (and the corresponding zero-based index) from a given c...
EIGEN_STRONG_INLINE const AdjointReturnType t() const
Transpose.
void saveToTextFile(const std::string &file, mrpt::math::TMatrixTextFileFormat fileFormat=mrpt::math::MATRIX_FORMAT_ENG, bool appendMRPTHeader=false, const std::string &userHeader=std::string()) const
Save matrix to a text file, compatible with MATLAB text format (see also the methods of matrix classe...
EIGEN_STRONG_INLINE void eigenVectorsSymmetricVec(MATRIX1 &eVecs, VECTOR1 &eVals) const
[For symmetric matrices only] Compute the eigenvectors and eigenvalues (in no particular order)...
const Scalar * const_iterator
Definition: eigen_plugins.h:24
EIGEN_STRONG_INLINE void multiply_AtB(const MATRIX1 &A, const MATRIX2 &B)
EIGEN_STRONG_INLINE bool eigenVectorsVec(MATRIX1 &eVecs, VECTOR1 &eVals) const
[For square matrices only] Compute the eigenvectors and eigenvalues (sorted), eigenvectors are the co...
EIGEN_STRONG_INLINE Scalar get_unsafe(const size_t row, const size_t col) const
Read-only access to one element (Use with caution, bounds are not checked!)
Definition: eigen_plugins.h:82
EIGEN_STRONG_INLINE Scalar squareNorm() const
Compute the square norm of a vector/array/matrix (the Euclidean distance to the origin, taking all the elements as a single vector).
GLdouble s
Definition: glext.h:3602
EIGEN_STRONG_INLINE Scalar norm_inf() const
Compute the norm-infinite of a vector ($f[ ||{v}||_ $f]), ie the maximum absolute value of the elemen...
EIGEN_STRONG_INLINE void minimum_maximum(Scalar &out_min, Scalar &out_max) const
[VECTORS OR MATRICES] Compute the minimum and maximum of a container at once
T square(const T x)
Inline function for the square of a number.
Definition: bits.h:52
void multiply_A_skew3(const MAT_A &A, const SKEW_3VECTOR &v)
EIGEN_STRONG_INLINE void insertMatrixTranspose(size_t r, size_t c, const MAT &m)
EIGEN_STRONG_INLINE size_t getColCount() const
Get number of columns.
Definition: eigen_plugins.h:48
EIGEN_STRONG_INLINE void laplacian(Eigen::MatrixBase< OtherDerived > &ret) const
Computes the laplacian of this square graph weight matrix.
EIGEN_STRONG_INLINE void eigenVectorsSymmetric(MATRIX1 &eVecs, MATRIX2 &eVals) const
[For symmetric matrices only] Compute the eigenvectors and eigenvalues (in no particular order)...
EIGEN_STRONG_INLINE void swapCols(size_t i1, size_t i2)
EIGEN_STRONG_INLINE void multiply_Atb(const OTHERVECTOR1 &vIn, OTHERVECTOR2 &vOut, bool accumToOutput=false) const
const GLubyte * c
Definition: glext.h:5590
CONTAINER::Scalar sum(const CONTAINER &v)
Computes the sum of all the elements.
EIGEN_STRONG_INLINE size_t countNonZero() const
EIGEN_STRONG_INLINE void insertMatrix(size_t r, size_t c, const MAT &m)
Insert matrix "m" into this matrix at indices (r,c), that is, (*this)(r,c)=m(0,0) and so on...
EIGEN_STRONG_INLINE void substract_At(const OTHERMATRIX &m)
EIGEN_STRONG_INLINE void substract_AAt(const OTHERMATRIX &A)
int val
Definition: mrpt_jpeglib.h:953
EIGEN_STRONG_INLINE Scalar multiply_HtCH_scalar(const MAT_C &C) const
EIGEN_STRONG_INLINE void setSize(size_t row, size_t col)
Changes the size of matrix, maintaining its previous content as possible and padding with zeros where...
MatrixBase< Derived > & operator^=(const unsigned int pow)
Combined matrix power and assignment operator.
EIGEN_STRONG_INLINE void multiply_AB(const MATRIX1 &A, const MATRIX2 &B)
EIGEN_STRONG_INLINE void insertCol(size_t nCol, const MAT &aCol)
EIGEN_STRONG_INLINE void multiply_AAt_scalar(const MAT_A &A, typename MAT_A::Scalar f)
EIGEN_STRONG_INLINE void unit(const size_t nRows, const Scalar diag_vals)
Make the matrix an identity matrix (the diagonal values can be 1.0 or any other value) ...
Definition: eigen_plugins.h:51
EIGEN_STRONG_INLINE void removeColumns(const std::vector< size_t > &idxsToRemove)
Remove columns of the matrix.
EIGEN_STRONG_INLINE void unsafeRemoveRows(const std::vector< size_t > &idxs)
Remove rows of the matrix.
TMatrixTextFileFormat
Definition: math_frwds.h:63
EIGEN_STRONG_INLINE void assign(const Scalar v)
Definition: eigen_plugins.h:41
GLsizei const GLchar ** string
Definition: glext.h:3919
EIGEN_STRONG_INLINE void extractSubmatrix(const size_t row_first, const size_t row_last, const size_t col_first, const size_t col_last, MATRIX &out) const
Get a submatrix, given its bounds: first & last column and row (inclusive).
EIGEN_STRONG_INLINE void unsafeRemoveColumns(const std::vector< size_t > &idxs)
Remove columns of the matrix.
void multiply_AtBC(const MAT_A &A, const MAT_B &B, const MAT_C &C)
void meanAndStdAll(double &outMean, double &outStd, const bool unbiased_variance=true) const
Computes the mean and standard deviation of all the elements in the matrix as a whole.
EIGEN_STRONG_INLINE void insertRow(size_t nRow, const MAT &aRow)
EIGEN_STRONG_INLINE void substract_An(const OTHERMATRIX &m, const size_t n)
EIGEN_STRONG_INLINE iterator end()
Definition: eigen_plugins.h:27
EIGEN_STRONG_INLINE Scalar maximum() const
[VECTORS OR MATRICES] Finds the maximum value
void loadFromTextFile(const std::string &file)
Load matrix from a text file, compatible with MATLAB text format.
void normalize(Scalar valMin, Scalar valMax)
Scales all elements such as the minimum & maximum values are shifted to the given values...
Internal resize which compiles to nothing on fixed-size matrices.
Definition: math_frwds.h:49
const GLdouble * v
Definition: glext.h:3603
EIGEN_STRONG_INLINE Scalar sumAll() const
GLdouble GLdouble GLdouble r
Definition: glext.h:3618
void largestEigenvector(OUTVECT &x, Scalar resolution=Scalar(0.01), size_t maxIterations=6, int *out_Iterations=NULL, float *out_estimatedResolution=NULL) const
Efficiently computes only the biggest eigenvector of the matrix using the Power Method, and returns it in the passed vector "x".
EIGEN_STRONG_INLINE void multiply_ABt(const MAT_A &A, const MAT_B &B)
const float R
EIGEN_STRONG_INLINE bool isDiagonal() const
Checks for matrix type.
EIGEN_STRONG_INLINE void multiply_AtA(const MAT_A &A)
GLenum GLenum GLvoid * row
Definition: glext.h:3533
void extractSubmatrixSymmetrical(const std::vector< size_t > &indices, MATRIX &out) const
Get a submatrix from a square matrix, by collecting the elements M(idxs,idxs), where idxs is the sequ...
EIGEN_STRONG_INLINE void multiply_HtCH(const MAT_C &C, MAT_R &R, bool accumResultInOutput=false) const
EIGEN_STRONG_INLINE void add_Ac(const OTHERMATRIX &m, const Scalar c)
EIGEN_STRONG_INLINE void extractMatrix(const size_t firstRow, const size_t firstCol, MATRIX &m) const
EIGEN_STRONG_INLINE void multiply_subMatrix(const MAT_A &A, MAT_OUT &outResult, const size_t A_cols_offset, const size_t A_rows_offset, const size_t A_col_count) const
outResult = this * A
EIGEN_STRONG_INLINE void zeros()
Set all elements to zero.
Definition: eigen_plugins.h:66
bool fromMatlabStringFormat(const std::string &s, std::ostream *dump_errors_here=NULL)
Read a matrix from a string in Matlab-like format, for example "[1 0 2; 0 4 -1]" The string must star...
EIGEN_STRONG_INLINE void inv_fast(MATRIX &outMat) const
EIGEN_STRONG_INLINE void eigenValues(VECTOR &eVals) const
[For square matrices only] Compute the eigenvectors and eigenvalues (sorted), and return only the eig...
EIGEN_STRONG_INLINE size_t getRowCount() const
Get number of rows.
Definition: eigen_plugins.h:46
void multiply_skew3_A(const SKEW_3VECTOR &v, const MAT_A &A, MAT_OUT &out)
Only for vectors/arrays "v" of length3, compute out = Skew(v) * A, where Skew(v) is the skew symmetri...
Definition: ops_matrices.h:169
EIGEN_STRONG_INLINE void multiply_AtA_scalar(const MAT_A &A, typename MAT_A::Scalar f)
GLsizeiptr size
Definition: glext.h:3779
void multiply_ABC(const MAT_A &A, const MAT_B &B, const MAT_C &C)
EIGEN_STRONG_INLINE void multiply_result_is_symmetric(const MAT_A &A, const MAT_B &B)
GLenum GLint x
Definition: glext.h:3516
EIGEN_STRONG_INLINE void eye()
Make the matrix an identity matrix.
Definition: eigen_plugins.h:63
EIGEN_STRONG_INLINE void multiplyRowByScalar(size_t r, Scalar s)
EIGEN_STRONG_INLINE void ones(const size_t row, const size_t col)
Resize matrix and set all elements to one.
Definition: eigen_plugins.h:71
EIGEN_STRONG_INLINE void extractRow(size_t nRow, Eigen::EigenBase< OtherDerived > &v, size_t startingCol=0) const
Extract one row from the matrix into a row vector.
GLsizei GLsizei GLenum GLenum const GLvoid * data
Definition: glext.h:3520
double Scalar
Definition: KmUtils.h:41
EIGEN_STRONG_INLINE Scalar minimum() const
[VECTORS OR MATRICES] Finds the minimum value
EIGEN_STRONG_INLINE double mean() const
Computes the mean of the entire matrix.
EIGEN_STRONG_INLINE void multiplyColumnByScalar(size_t c, Scalar s)
EIGEN_STRONG_INLINE PlainObject inv() const
EIGEN_STRONG_INLINE void multiply_AAt(const MAT_A &A)



Page generated by Doxygen 1.8.14 for MRPT 1.5.6 Git: 4c65e8431 Tue Apr 24 08:18:17 2018 +0200 at lun oct 28 01:35:26 CET 2019