MRPT  1.9.9
CSparseMatrix.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-2018, 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 #pragma once
10 
11 #include <mrpt/math/math_frwds.h>
15 #include <cstring> // memcpy
16 #include <stdexcept>
17 
18 // Include CSparse lib headers, either from the system or embedded:
19 extern "C" {
20 #if MRPT_HAS_SUITESPARSE
21 #define NCOMPLEX // In MRPT we don't need complex numbers, so avoid the
22 // annoying warning: 'cs_ci_house' has C-linkage specified,
23 // but returns UDT 'std::complex<double>' which is
24 // incompatible with C
25 #include "cs.h"
26 #else
27 #include <mrpt/otherlibs/CSparse/cs.h>
28 #endif
29 }
30 
31 namespace mrpt::math
32 {
33 /** Used in mrpt::math::CSparseMatrix */
34 struct CExceptionNotDefPos : public std::runtime_error
35 {
36  CExceptionNotDefPos(const char* s) : std::runtime_error(s) {}
37 };
38 
39 /** A sparse matrix structure, wrapping T. Davis' CSparse library (part of
40  *suitesparse)
41  * The type of the matrix entries is fixed to "double".
42  *
43  * There are two formats for the non-zero entries in this matrix:
44  * - A "triplet" matrix: a set of (r,c)=val triplet entries.
45  * - A column-compressed sparse (CCS) matrix.
46  *
47  * The latter is the "normal" format, which is expected by all mathematical
48  *operations defined
49  * in this class. There're three ways of initializing and populating a sparse
50  *matrix:
51  *
52  * <ol>
53  * <li> <b>As a triplet (empty), then add entries, then compress:</b>
54  * \code
55  * CSparseMatrix SM(100,100);
56  * SM.insert_entry(i,j, val); // or
57  * SM.insert_submatrix(i,j, MAT); // ...
58  * SM.compressFromTriplet();
59  * \endcode
60  * </li>
61  * <li> <b>As a triplet from a CSparseMatrixTemplate<double>:</b>
62  * \code
63  * CSparseMatrixTemplate<double> data;
64  * data(row,col) = val;
65  * ...
66  * CSparseMatrix SM(data);
67  * \endcode
68  * </li>
69  * <li> <b>From an existing dense matrix:</b>
70  * \code
71  * CMatrixDouble data(100,100); // or
72  * CMatrixFloat data(100,100); // or
73  * CMatrixFixedNumeric<double,4,6> data; // etc...
74  * CSparseMatrix SM(data);
75  * \endcode
76  * </li>
77  *
78  * </ol>
79  *
80  * Due to its practical utility, there is a special inner class
81  *CSparseMatrix::CholeskyDecomp to handle Cholesky-related methods and data.
82  *
83  * \note This class was initially adapted from "robotvision", by Hauke
84  *Strasdat, Steven Lovegrove and Andrew J. Davison. See
85  *http://www.openslam.org/robotvision.html
86  * \note CSparse is maintained by Timothy Davis:
87  *http://people.sc.fsu.edu/~jburkardt/c_src/csparse/csparse.html .
88  * \note See also his book "Direct methods for sparse linear systems".
89  *http://books.google.es/books?id=TvwiyF8vy3EC&pg=PA12&lpg=PA12&dq=cs_compress&source=bl&ots=od9uGJ793j&sig=Wa-fBk4sZkZv3Y0Op8FNH8PvCUs&hl=es&ei=UjA0TJf-EoSmsQay3aXPAw&sa=X&oi=book_result&ct=result&resnum=8&ved=0CEQQ6AEwBw#v=onepage&q&f=false
90  *
91  * \sa mrpt::math::MatrixBlockSparseCols, mrpt::math::CMatrixFixedNumeric,
92  *mrpt::math::CMatrixTemplateNumeric, etc.
93  * \ingroup mrpt_math_grp
94  */
96 {
97  private:
99 
100  /** Initialization from a dense matrix of any kind existing in MRPT. */
101  template <class MATRIX>
102  void construct_from_mrpt_mat(const MATRIX& C)
103  {
104  std::vector<int> row_list, col_list; // Use "int" for convenience so we
105  // can do a memcpy below...
106  std::vector<double> content_list;
107  const int nCol = C.cols();
108  const int nRow = C.rows();
109  for (int c = 0; c < nCol; ++c)
110  {
111  col_list.push_back(row_list.size());
112  for (int r = 0; r < nRow; ++r)
113  if (C.get_unsafe(r, c) != 0)
114  {
115  row_list.push_back(r);
116  content_list.push_back(C(r, c));
117  }
118  }
119  col_list.push_back(row_list.size());
120 
121  sparse_matrix.m = nRow;
122  sparse_matrix.n = nCol;
123  sparse_matrix.nzmax = content_list.size();
124  sparse_matrix.i = (int*)malloc(sizeof(int) * row_list.size());
125  sparse_matrix.p = (int*)malloc(sizeof(int) * col_list.size());
126  sparse_matrix.x = (double*)malloc(sizeof(double) * content_list.size());
127 
128  ::memcpy(
129  sparse_matrix.i, &row_list[0],
130  sizeof(row_list[0]) * row_list.size());
131  ::memcpy(
132  sparse_matrix.p, &col_list[0],
133  sizeof(col_list[0]) * col_list.size());
134  ::memcpy(
135  sparse_matrix.x, &content_list[0],
136  sizeof(content_list[0]) * content_list.size());
137 
138  sparse_matrix.nz =
139  -1; // <0 means "column compressed", ">=0" means triplet.
140  }
141 
142  /** Initialization from a triplet "cs", which is first compressed */
143  void construct_from_triplet(const cs& triplet);
144 
145  /** To be called by constructors only, assume previous pointers are trash
146  * and overwrite them */
147  void construct_from_existing_cs(const cs& sm);
148 
149  /** free buffers (deallocate the memory of the i,p,x buffers) */
150  void internal_free_mem();
151 
152  /** Copy the data from an existing "cs" CSparse data structure */
153  void copy(const cs* const sm);
154 
155  /** Fast copy the data from an existing "cs" CSparse data structure, copying
156  * the pointers and leaving NULLs in the source structure. */
157  void copy_fast(cs* const sm);
158 
159  public:
160  /** @name Constructors, destructor & copy operations
161  @{ */
162 
163  /** Create an initially empty sparse matrix, in the "triplet" form.
164  * Notice that you must call "compressFromTriplet" after populating the
165  * matrix and before using the math operatons on this matrix.
166  * The initial size can be later on extended with insert_entry() or
167  * setRowCount() & setColCount().
168  * \sa insert_entry, setRowCount, setColCount
169  */
170  CSparseMatrix(const size_t nRows = 0, const size_t nCols = 0);
171 
172  /** A good way to initialize a sparse matrix from a list of non nullptr
173  * elements.
174  * This constructor takes all the non-zero entries in "data" and builds a
175  * column-compressed sparse representation.
176  */
177  template <typename T>
179  {
180  ASSERTMSG_(
181  !data.empty(),
182  "Input data must contain at least one non-zero element.");
183  sparse_matrix.i =
184  nullptr; // This is to know they shouldn't be tried to free()
185  sparse_matrix.p = nullptr;
186  sparse_matrix.x = nullptr;
187 
188  // 1) Create triplet matrix
189  CSparseMatrix triplet(data.rows(), data.cols());
190  // 2) Put data in:
192  data.begin();
193  it != data.end(); ++it)
194  triplet.insert_entry_fast(
195  it->first.first, it->first.second, it->second);
196 
197  // 3) Compress:
198  construct_from_triplet(triplet.sparse_matrix);
199  }
200 
201  // We can't do a simple "template <class ANYMATRIX>" since it would be tried
202  // to match against "cs*"...
203 
204  /** Constructor from a dense matrix of any kind existing in MRPT, creating a
205  * "column-compressed" sparse matrix. */
206  template <typename T, size_t N, size_t M>
207  inline explicit CSparseMatrix(const CMatrixFixedNumeric<T, N, M>& MAT)
208  {
210  }
211 
212  /** Constructor from a dense matrix of any kind existing in MRPT, creating a
213  * "column-compressed" sparse matrix. */
214  template <typename T>
215  inline explicit CSparseMatrix(const CMatrixTemplateNumeric<T>& MAT)
216  {
218  }
219 
220  /** Copy constructor */
221  CSparseMatrix(const CSparseMatrix& other);
222 
223  /** Copy constructor from an existing "cs" CSparse data structure */
224  explicit CSparseMatrix(const cs* const sm);
225 
226  /** Destructor */
227  virtual ~CSparseMatrix();
228 
229  /** Copy operator from another existing object */
230  void operator=(const CSparseMatrix& other);
231 
232  /** Fast swap contents with another sparse matrix */
233  void swap(CSparseMatrix& other);
234 
235  /** Erase all previous contents and leave the matrix as a "triplet" ROWS x
236  * COLS matrix without any nonzero entry. */
237  void clear(const size_t nRows = 1, const size_t nCols = 1);
238 
239  /** @} */
240 
241  /** @name Math operations (the interesting stuff...)
242  @{ */
243 
244  /** this = A+B */
245  void add_AB(const CSparseMatrix& A, const CSparseMatrix& B);
246  /** this = A*B */
247  void multiply_AB(const CSparseMatrix& A, const CSparseMatrix& B);
248  /** out_res = this * b */
249  void multiply_Ab(
251  mrpt::math::CVectorDouble& out_res) const;
252 
253  inline CSparseMatrix operator+(const CSparseMatrix& other) const
254  {
255  CSparseMatrix RES;
256  RES.add_AB(*this, other);
257  return RES;
258  }
259  inline CSparseMatrix operator*(const CSparseMatrix& other) const
260  {
261  CSparseMatrix RES;
262  RES.multiply_AB(*this, other);
263  return RES;
264  }
266  const mrpt::math::CVectorDouble& other) const
267  {
269  multiply_Ab(other, res);
270  return res;
271  }
272  inline void operator+=(const CSparseMatrix& other)
273  {
274  this->add_AB(*this, other);
275  }
276  inline void operator*=(const CSparseMatrix& other)
277  {
278  this->multiply_AB(*this, other);
279  }
280  CSparseMatrix transpose() const;
281 
282  /** @} */
283 
284  /** @ Access the matrix, get, set elements, etc.
285  @{ */
286 
287  /** ONLY for TRIPLET matrices: insert a new non-zero entry in the matrix.
288  * This method cannot be used once the matrix is in column-compressed
289  * form.
290  * The size of the matrix will be automatically extended if the indices
291  * are out of the current limits.
292  * \sa isTriplet, compressFromTriplet
293  */
294  void insert_entry(const size_t row, const size_t col, const double val);
295 
296  /** This was an optimized version, but is now equivalent to insert_entry()
297  * due to the need to be compatible with unmodified CSparse system
298  * libraries. */
299  inline void insert_entry_fast(
300  const size_t row, const size_t col, const double val)
301  {
302  insert_entry(row, col, val);
303  }
304 
305  /** ONLY for TRIPLET matrices: insert a given matrix (in any of the MRPT
306  * formats) at a given location of the sparse matrix.
307  * This method cannot be used once the matrix is in column-compressed
308  * form.
309  * The size of the matrix will be automatically extended if the indices
310  * are out of the current limits.
311  * Since this is inline, it can be very efficient for fixed-size MRPT
312  * matrices.
313  * \sa isTriplet, compressFromTriplet, insert_entry
314  */
315  template <class MATRIX>
316  inline void insert_submatrix(
317  const size_t row, const size_t col, const MATRIX& M)
318  {
319  if (!isTriplet())
321  "insert_entry() is only available for sparse matrix in "
322  "'triplet' format.")
323  const size_t nR = M.rows();
324  const size_t nC = M.cols();
325  for (size_t r = 0; r < nR; r++)
326  for (size_t c = 0; c < nC; c++)
327  insert_entry_fast(row + r, col + c, M.get_unsafe(r, c));
328  // If needed, extend the size of the matrix:
329  sparse_matrix.m = std::max(sparse_matrix.m, int(row + nR));
330  sparse_matrix.n = std::max(sparse_matrix.n, int(col + nC));
331  }
332 
333  /** ONLY for TRIPLET matrices: convert the matrix in a column-compressed
334  * form.
335  * \sa insert_entry
336  */
337  void compressFromTriplet();
338 
339  /** Return a dense representation of the sparse matrix.
340  * \sa saveToTextFile_dense
341  */
342  void get_dense(CMatrixDouble& outMat) const;
343 
344  /** Static method to convert a "cs" structure into a dense representation of
345  * the sparse matrix.
346  */
347  static void cs2dense(const cs& SM, CMatrixDouble& outMat);
348 
349  /** save as a dense matrix to a text file \return False on any error.
350  */
351  bool saveToTextFile_dense(const std::string& filName);
352 
353  /** Save sparse structure to a text file loadable from MATLAB (can be called
354  * on triplet or CCS matrices).
355  *
356  * The format of the text file is:
357  * \code
358  * NUM_ROWS NUM_COLS NUM_NON_ZERO_MAX
359  * row_1 col_1 value_1
360  * row_2 col_2 value_2
361  * ...
362  * \endcode
363  *
364  * Instructions for loading from MATLAB in triplet form will be
365  * automatically writen to the
366  * output file as comments in th first lines:
367  *
368  * \code
369  * D=load('file.txt');
370  * SM=spconvert(D(2:end,:));
371  * or, to always preserve the actual matrix size m x n:
372  * m=D(1,1); n=D(1,2); nzmax=D(1,3);
373  * Di=D(2:end,1); Dj=D(2:end,2); Ds=D(2:end,3);
374  * M=sparse(Di,Dj,Ds, m,n, nzmax);
375  * \endcode
376  * \return False on any error.
377  */
378  bool saveToTextFile_sparse(const std::string& filName);
379 
380  // Very basic, standard methods that MRPT methods expect for any matrix:
381  inline size_t rows() const { return sparse_matrix.m; }
382  inline size_t cols() const { return sparse_matrix.n; }
383  /** Change the number of rows in the matrix (can't be lower than current
384  * size) */
385  inline void setRowCount(const size_t nRows)
386  {
387  ASSERT_(nRows >= (size_t)sparse_matrix.m);
388  sparse_matrix.m = nRows;
389  }
390  inline void setColCount(const size_t nCols)
391  {
392  ASSERT_(nCols >= (size_t)sparse_matrix.n);
393  sparse_matrix.n = nCols;
394  }
395 
396  /** Returns true if this sparse matrix is in "triplet" form. \sa
397  * isColumnCompressed */
398  inline bool isTriplet() const
399  {
400  return sparse_matrix.nz >= 0;
401  } // <0 means "column compressed", ">=0" means triplet.
402 
403  /** Returns true if this sparse matrix is in "column compressed" form. \sa
404  * isTriplet */
405  inline bool isColumnCompressed() const
406  {
407  return sparse_matrix.nz < 0;
408  } // <0 means "column compressed", ">=0" means triplet.
409 
410  /** @} */
411 
412  /** @name Cholesky factorization
413  @{ */
414 
415  /** Auxiliary class to hold the results of a Cholesky factorization of a
416  *sparse matrix.
417  * This implementation does not allow updating/downdating.
418  *
419  * Usage example:
420  * \code
421  * CSparseMatrix SM(100,100);
422  * SM.insert_entry(i,j, val); ...
423  * SM.compressFromTriplet();
424  *
425  * // Do Cholesky decomposition:
426  * CSparseMatrix::CholeskyDecomp CD(SM);
427  * CD.get_inverse();
428  * ...
429  * \endcode
430  *
431  * \note Only the upper triangular part of the input matrix is accessed.
432  * \note This class was initially adapted from "robotvision", by Hauke
433  *Strasdat, Steven Lovegrove and Andrew J. Davison. See
434  *http://www.openslam.org/robotvision.html
435  * \note This class designed to be "uncopiable".
436  * \sa The main class: CSparseMatrix
437  */
439  {
440  private:
443  /** A const reference to the original matrix used to build this
444  * decomposition. */
446 
447  public:
448  /** Constructor from a square definite-positive sparse matrix A, which
449  * can be use to solve Ax=b
450  * The actual Cholesky decomposition takes places in this
451  * constructor.
452  * \note Only the upper triangular part of the matrix is accessed.
453  * \exception std::runtime_error On non-square input matrix.
454  * \exception mrpt::math::CExceptionNotDefPos On non-definite-positive
455  * matrix as input.
456  */
457  CholeskyDecomp(const CSparseMatrix& A);
458  CholeskyDecomp(const CholeskyDecomp& A) = delete;
459 
460  CholeskyDecomp& operator=(const CholeskyDecomp&) = delete;
461 
462  /** Destructor */
463  virtual ~CholeskyDecomp();
464 
465  /** Return the L matrix (L*L' = M), as a dense matrix. */
466  inline CMatrixDouble get_L() const
467  {
468  CMatrixDouble L;
469  get_L(L);
470  return L;
471  }
472 
473  /** Return the L matrix (L*L' = M), as a dense matrix. */
474  void get_L(CMatrixDouble& out_L) const;
475 
476  /** Return the vector from a back-substitution step that solves: Ux=b */
477  template <class VECTOR>
478  inline VECTOR backsub(const VECTOR& b) const
479  {
480  VECTOR res;
481  backsub(b, res);
482  return res;
483  }
484 
485  /** Return the vector from a back-substitution step that solves: Ux=b.
486  * Vectors can be Eigen::VectorXd or mrpt::math::CVectorDouble */
487  void backsub(const Eigen::VectorXd& b, Eigen::VectorXd& result_x) const;
488 
489  /** overload for double pointers which assume the user has reserved the
490  * output memory for \a result */
491  void backsub(const double* b, double* result, const size_t N) const;
492 
493  /** Update the Cholesky factorization from an updated vesion of the
494  * original input, square definite-positive sparse matrix.
495  * NOTE: This new matrix MUST HAVE exactly the same sparse structure
496  * than the original one.
497  */
498  void update(const CSparseMatrix& new_SM);
499  };
500  static_assert(
503  "Copy Check");
504 
505  /** @} */
506 
507 }; // end class CSparseMatrix
508 }
509 
void setColCount(const size_t nCols)
#define THROW_EXCEPTION(msg)
Definition: exceptions.h:41
const CSparseMatrix * m_originalSM
A const reference to the original matrix used to build this decomposition.
Used in mrpt::math::CSparseMatrix.
Definition: CSparseMatrix.h:34
void insert_entry_fast(const size_t row, const size_t col, const double val)
This was an optimized version, but is now equivalent to insert_entry() due to the need to be compatib...
void construct_from_triplet(const cs &triplet)
Initialization from a triplet "cs", which is first compressed.
Column vector, like Eigen::MatrixX*, but automatically initialized to zeros since construction...
Definition: eigen_frwds.h:44
Auxiliary class to hold the results of a Cholesky factorization of a sparse matrix.
A sparse matrix structure, wrapping T.
Definition: CSparseMatrix.h:95
STL namespace.
void compressFromTriplet()
ONLY for TRIPLET matrices: convert the matrix in a column-compressed form.
void construct_from_mrpt_mat(const MATRIX &C)
Initialization from a dense matrix of any kind existing in MRPT.
GLdouble s
Definition: glext.h:3676
void swap(CSparseMatrix &other)
Fast swap contents with another sparse matrix.
void insert_entry(const size_t row, const size_t col, const double val)
@ Access the matrix, get, set elements, etc.
A sparse matrix container (with cells of any type), with iterators.
CMatrixDouble get_L() const
Return the L matrix (L*L&#39; = M), as a dense matrix.
void get_dense(CMatrixDouble &outMat) const
Return a dense representation of the sparse matrix.
#define ASSERT_(f)
Defines an assertion mechanism.
Definition: exceptions.h:113
A numeric matrix of compile-time fixed size.
CSparseMatrix operator*(const CSparseMatrix &other) const
This base provides a set of functions for maths stuff.
void multiply_AB(const CSparseMatrix &A, const CSparseMatrix &B)
this = A*B
CSparseMatrix(const CMatrixTemplateNumeric< T > &MAT)
Constructor from a dense matrix of any kind existing in MRPT, creating a "column-compressed" sparse m...
void construct_from_existing_cs(const cs &sm)
To be called by constructors only, assume previous pointers are trash and overwrite them...
CSparseMatrix(const size_t nRows=0, const size_t nCols=0)
Create an initially empty sparse matrix, in the "triplet" form.
const GLubyte * c
Definition: glext.h:6313
void add_AB(const CSparseMatrix &A, const CSparseMatrix &B)
this = A+B
int val
Definition: mrpt_jpeglib.h:955
void operator=(const CSparseMatrix &other)
Copy operator from another existing object.
GLubyte GLubyte b
Definition: glext.h:6279
#define ASSERTMSG_(f, __ERROR_MSG)
Defines an assertion mechanism.
Definition: exceptions.h:101
void multiply_Ab(const mrpt::math::CVectorDouble &b, mrpt::math::CVectorDouble &out_res) const
out_res = this * b
void setRowCount(const size_t nRows)
Change the number of rows in the matrix (can&#39;t be lower than current size)
GLsizei const GLchar ** string
Definition: glext.h:4101
CSparseMatrix operator+(const CSparseMatrix &other) const
void clear(const size_t nRows=1, const size_t nCols=1)
Erase all previous contents and leave the matrix as a "triplet" ROWS x COLS matrix without any nonzer...
bool isTriplet() const
Returns true if this sparse matrix is in "triplet" form.
CholeskyDecomp & operator=(const CholeskyDecomp &)=delete
CSparseMatrix(const CMatrixFixedNumeric< T, N, M > &MAT)
Constructor from a dense matrix of any kind existing in MRPT, creating a "column-compressed" sparse m...
GLdouble GLdouble GLdouble r
Definition: glext.h:3705
mrpt::math::CVectorDouble operator*(const mrpt::math::CVectorDouble &other) const
void operator*=(const CSparseMatrix &other)
CholeskyDecomp(const CSparseMatrix &A)
Constructor from a square definite-positive sparse matrix A, which can be use to solve Ax=b The actua...
bool saveToTextFile_dense(const std::string &filName)
save as a dense matrix to a text file
void update(const CSparseMatrix &new_SM)
Update the Cholesky factorization from an updated vesion of the original input, square definite-posit...
GLenum GLenum GLvoid * row
Definition: glext.h:3576
void operator+=(const CSparseMatrix &other)
virtual ~CSparseMatrix()
Destructor.
void insert_submatrix(const size_t row, const size_t col, const MATRIX &M)
ONLY for TRIPLET matrices: insert a given matrix (in any of the MRPT formats) at a given location of ...
CSparseMatrix transpose() const
bool saveToTextFile_sparse(const std::string &filName)
Save sparse structure to a text file loadable from MATLAB (can be called on triplet or CCS matrices)...
static void cs2dense(const cs &SM, CMatrixDouble &outMat)
Static method to convert a "cs" structure into a dense representation of the sparse matrix...
GLsizei const GLfloat * value
Definition: glext.h:4117
GLuint res
Definition: glext.h:7268
VECTOR backsub(const VECTOR &b) const
Return the vector from a back-substitution step that solves: Ux=b.
CSparseMatrix(const CSparseMatrixTemplate< T > &data)
A good way to initialize a sparse matrix from a list of non nullptr elements.
bool isColumnCompressed() const
Returns true if this sparse matrix is in "column compressed" form.
void copy(const cs *const sm)
Copy the data from an existing "cs" CSparse data structure.
GLsizei GLsizei GLenum GLenum const GLvoid * data
Definition: glext.h:3546
void copy_fast(cs *const sm)
Fast copy the data from an existing "cs" CSparse data structure, copying the pointers and leaving NUL...
typename SparseMatrixMap::const_iterator const_iterator
Const iterator to move through the matrix.
void internal_free_mem()
free buffers (deallocate the memory of the i,p,x buffers)
void memcpy(void *dest, size_t destSize, const void *src, size_t copyCount) noexcept
An OS and compiler independent version of "memcpy".
Definition: os.cpp:356



Page generated by Doxygen 1.8.14 for MRPT 1.9.9 Git: 7d5e6d718 Fri Aug 24 01:51:28 2018 +0200 at lun nov 2 08:35:50 CET 2020