g2o
Loading...
Searching...
No Matches
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
g2o::LinearSolverCSparse< MatrixType > Class Template Reference

linear solver which uses CSparse More...

#include <linear_solver_csparse.h>

Inheritance diagram for g2o::LinearSolverCSparse< MatrixType >:
Inheritance graph
[legend]
Collaboration diagram for g2o::LinearSolverCSparse< MatrixType >:
Collaboration graph
[legend]

Public Member Functions

 LinearSolverCSparse ()
 
 LinearSolverCSparse (LinearSolverCSparse< MatrixType > const &)=delete
 
LinearSolverCSparseoperator= (LinearSolverCSparse< MatrixType > const &)=delete
 
virtual ~LinearSolverCSparse ()=default
 
virtual bool init ()
 
bool solve (const SparseBlockMatrix< MatrixType > &A, double *x, double *b)
 
- Public Member Functions inherited from g2o::LinearSolverCCS< MatrixType >
 LinearSolverCCS ()
 
 ~LinearSolverCCS ()
 
virtual bool solveBlocks (double **&blocks, const SparseBlockMatrix< MatrixType > &A)
 
virtual bool solvePattern (SparseBlockMatrix< MatrixX > &spinv, const std::vector< std::pair< int, int > > &blockIndices, const SparseBlockMatrix< MatrixType > &A)
 
bool blockOrdering () const
 do the AMD ordering on the blocks or on the scalar matrix
 
void setBlockOrdering (bool blockOrdering)
 
- Public Member Functions inherited from g2o::LinearSolver< MatrixType >
 LinearSolver ()
 
virtual ~LinearSolver ()
 
bool writeDebug () const
 write a debug dump of the system matrix if it is not PSD in solve
 
void setWriteDebug (bool b)
 

Protected Member Functions

void prepareSolve (const SparseBlockMatrix< MatrixType > &A)
 
void computeSymbolicDecomposition (const SparseBlockMatrix< MatrixType > &A)
 
void fillCSparse (const SparseBlockMatrix< MatrixType > &A, bool onlyValues)
 
bool solveBlocks_impl (const SparseBlockMatrix< MatrixType > &A, std::function< void(MarginalCovarianceCholesky &)> compute)
 
- Protected Member Functions inherited from g2o::LinearSolverCCS< MatrixType >
void initMatrixStructure (const SparseBlockMatrix< MatrixType > &A)
 

Protected Attributes

csparse::CSparse csparse
 
MatrixStructure _matrixStructure
 
- Protected Attributes inherited from g2o::LinearSolverCCS< MatrixType >
SparseBlockMatrixCCS< MatrixType > * _ccsMatrix
 
bool _blockOrdering
 
- Protected Attributes inherited from g2o::LinearSolver< MatrixType >
bool _writeDebug
 

Additional Inherited Members

- Static Public Member Functions inherited from g2o::LinearSolver< MatrixType >
static void allocateBlocks (const SparseBlockMatrix< MatrixType > &A, double **&blocks)
 allocate block memory structure
 
static void deallocateBlocks (const SparseBlockMatrix< MatrixType > &A, double **&blocks)
 de-allocate the block structure
 
template<typename BlockDerived , typename ScalarDerived >
static void blockToScalarPermutation (const SparseBlockMatrix< MatrixType > &A, const Eigen::MatrixBase< BlockDerived > &p, const Eigen::MatrixBase< ScalarDerived > &scalar)
 

Detailed Description

template<typename MatrixType>
class g2o::LinearSolverCSparse< MatrixType >

linear solver which uses CSparse

Definition at line 47 of file linear_solver_csparse.h.

Constructor & Destructor Documentation

◆ LinearSolverCSparse() [1/2]

template<typename MatrixType >
g2o::LinearSolverCSparse< MatrixType >::LinearSolverCSparse ( )
inline

Definition at line 49 of file linear_solver_csparse.h.

49: LinearSolverCCS<MatrixType>() {}

◆ LinearSolverCSparse() [2/2]

template<typename MatrixType >
g2o::LinearSolverCSparse< MatrixType >::LinearSolverCSparse ( LinearSolverCSparse< MatrixType > const &  )
delete

◆ ~LinearSolverCSparse()

template<typename MatrixType >
virtual g2o::LinearSolverCSparse< MatrixType >::~LinearSolverCSparse ( )
virtualdefault

Member Function Documentation

◆ computeSymbolicDecomposition()

template<typename MatrixType >
void g2o::LinearSolverCSparse< MatrixType >::computeSymbolicDecomposition ( const SparseBlockMatrix< MatrixType > &  A)
inlineprotected

Definition at line 98 of file linear_solver_csparse.h.

98 {
99 double t = get_monotonic_time();
100 if (!this->blockOrdering()) {
102 } else {
103 A.fillBlockStructure(_matrixStructure);
104
105 // prepare block structure for the CSparse call
106 double* structureX = nullptr;
107 int structureNz = _matrixStructure.nzMax();
108 int structureAllocated = _matrixStructure.n;
109 csparse::CSparse::SparseView auxBlock(
110 _matrixStructure.n, _matrixStructure.n, structureNz,
112 structureAllocated);
113
114 // AMD ordering on the block structure
115 VectorXI blockPermutation;
116 csparse.amd(auxBlock, blockPermutation);
117
118 // blow up the permutation to the scalar matrix
119 VectorXI scalarPermutation;
120 this->blockToScalarPermutation(A, blockPermutation, scalarPermutation);
121 csparse.analyze_p(scalarPermutation.data());
122 }
123 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
124 if (globalStats) {
125 globalStats->timeSymbolicDecomposition = get_monotonic_time() - t;
126 }
127 }
bool blockOrdering() const
do the AMD ordering on the blocks or on the scalar matrix
static void blockToScalarPermutation(const SparseBlockMatrix< MatrixType > &A, const Eigen::MatrixBase< BlockDerived > &p, const Eigen::MatrixBase< ScalarDerived > &scalar)
int * Aii
row indices of A, of size nz = Ap [n]
int n
A is m-by-n. n must be >= 0.
int nzMax() const
max number of non-zeros blocks
int * Ap
column pointers for A, of size n+1
bool amd(const SparseView &sparseView, VectorXI &result)
compute AMD ordering on the given SparseView, store into result
bool analyze_p(int *permutation)
Eigen::Matrix< int, Eigen::Dynamic, 1, Eigen::ColMajor > VectorXI
Definition eigen_types.h:41
double get_monotonic_time()
Definition timeutil.cpp:43
static G2OBatchStatistics * globalStats()
Definition batch_stats.h:77

References g2o::LinearSolverCSparse< MatrixType >::_matrixStructure, g2o::MatrixStructure::Aii, g2o::csparse::CSparse::amd(), g2o::csparse::CSparse::analyze(), g2o::csparse::CSparse::analyze_p(), g2o::MatrixStructure::Ap, g2o::LinearSolverCCS< MatrixType >::blockOrdering(), g2o::LinearSolver< MatrixType >::blockToScalarPermutation(), g2o::LinearSolverCSparse< MatrixType >::csparse, g2o::SparseBlockMatrix< MatrixType >::fillBlockStructure(), g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), g2o::MatrixStructure::n, g2o::MatrixStructure::nzMax(), and g2o::G2OBatchStatistics::timeSymbolicDecomposition.

Referenced by g2o::LinearSolverCSparse< MatrixType >::prepareSolve().

◆ fillCSparse()

template<typename MatrixType >
void g2o::LinearSolverCSparse< MatrixType >::fillCSparse ( const SparseBlockMatrix< MatrixType > &  A,
bool  onlyValues 
)
inlineprotected

Definition at line 129 of file linear_solver_csparse.h.

129 {
130 if (!onlyValues) this->initMatrixStructure(A);
131 int m = A.rows();
132 int n = A.cols();
133 assert(m > 0 && n > 0 && "Hessian has 0 rows/cols");
134
135 csparse::CSparse::SparseView ccsA = csparse.sparseView();
136
137 if (ccsA.columnsAllocated < n) {
138 // pre-allocate more space if re-allocating
139 ccsA.columnsAllocated = ccsA.columnsAllocated == 0 ? n : 2 * n;
140 delete[] ccsA.p;
141 ccsA.p = new int[ccsA.columnsAllocated + 1];
142 }
143
144 if (!onlyValues) {
145 int nzmax = A.nonZeros();
146 if (ccsA.nzmax < nzmax) {
147 // pre-allocate more space if re-allocating
148 ccsA.nzmax = ccsA.nzmax == 0 ? nzmax : 2 * nzmax;
149 delete[] ccsA.x;
150 delete[] ccsA.i;
151 ccsA.i = new int[ccsA.nzmax];
152 ccsA.x = new double[ccsA.nzmax];
153 }
154 }
155 ccsA.m = m;
156 ccsA.n = n;
157
158 if (onlyValues) {
159 this->_ccsMatrix->fillCCS(ccsA.x, true);
160 } else {
161 int nz = this->_ccsMatrix->fillCCS(ccsA.p, ccsA.i, ccsA.x, true);
162 (void)nz;
163 assert(nz <= ccsA.nzmax);
164 }
165 }
void initMatrixStructure(const SparseBlockMatrix< MatrixType > &A)
SparseBlockMatrixCCS< MatrixType > * _ccsMatrix

References g2o::LinearSolverCCS< MatrixType >::_ccsMatrix, g2o::SparseBlockMatrix< MatrixType >::cols(), g2o::csparse::CSparse::SparseView::columnsAllocated, g2o::LinearSolverCSparse< MatrixType >::csparse, g2o::csparse::CSparse::SparseView::i, g2o::LinearSolverCCS< MatrixType >::initMatrixStructure(), g2o::csparse::CSparse::SparseView::m, g2o::csparse::CSparse::SparseView::n, g2o::SparseBlockMatrix< MatrixType >::nonZeros(), g2o::csparse::CSparse::SparseView::nzmax, g2o::csparse::CSparse::SparseView::p, g2o::SparseBlockMatrix< MatrixType >::rows(), g2o::csparse::CSparse::sparseView(), and g2o::csparse::CSparse::SparseView::x.

Referenced by g2o::LinearSolverCSparse< MatrixType >::prepareSolve().

◆ init()

template<typename MatrixType >
virtual bool g2o::LinearSolverCSparse< MatrixType >::init ( )
inlinevirtual

init for operating on matrices with a different non-zero pattern like before

Implements g2o::LinearSolver< MatrixType >.

Definition at line 57 of file linear_solver_csparse.h.

57 {
59 return true;
60 }

References g2o::LinearSolverCSparse< MatrixType >::csparse, and g2o::csparse::CSparse::freeSymbolic().

◆ operator=()

template<typename MatrixType >
LinearSolverCSparse & g2o::LinearSolverCSparse< MatrixType >::operator= ( LinearSolverCSparse< MatrixType > const &  )
delete

◆ prepareSolve()

template<typename MatrixType >
void g2o::LinearSolverCSparse< MatrixType >::prepareSolve ( const SparseBlockMatrix< MatrixType > &  A)
inlineprotected

Definition at line 88 of file linear_solver_csparse.h.

88 {
89 bool hasSymbolic = csparse.hasSymbolic();
90 fillCSparse(A, hasSymbolic);
91 // perform symbolic cholesky once
92 if (!hasSymbolic) {
94 assert(csparse.hasSymbolic() && "Symbolic cholesky failed");
95 }
96 }
void computeSymbolicDecomposition(const SparseBlockMatrix< MatrixType > &A)
void fillCSparse(const SparseBlockMatrix< MatrixType > &A, bool onlyValues)

References g2o::LinearSolverCSparse< MatrixType >::computeSymbolicDecomposition(), g2o::LinearSolverCSparse< MatrixType >::csparse, g2o::LinearSolverCSparse< MatrixType >::fillCSparse(), and g2o::csparse::CSparse::hasSymbolic().

Referenced by g2o::LinearSolverCSparse< MatrixType >::solve(), and g2o::LinearSolverCSparse< MatrixType >::solveBlocks_impl().

◆ solve()

template<typename MatrixType >
bool g2o::LinearSolverCSparse< MatrixType >::solve ( const SparseBlockMatrix< MatrixType > &  A,
double *  x,
double *  b 
)
inlinevirtual

Assumes that A is the same matrix for several calls. Among other assumptions, the non-zero pattern does not change! If the matrix changes call init() before. solve system Ax = b, x and b have to allocated beforehand!!

Implements g2o::LinearSolver< MatrixType >.

Definition at line 62 of file linear_solver_csparse.h.

62 {
63 prepareSolve(A);
64
65 double t = get_monotonic_time();
66 bool ok = csparse.solve(x, b);
67 if (!ok && this->writeDebug()) {
69 "Cholesky failure, writing debug.txt (Hessian loadable by Octave)");
70 csparse.writeSparse("debug.txt");
71 } else {
72 G2O_DEBUG("Cholesky failure");
73 }
74
75 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
76 if (globalStats) {
77 globalStats->timeNumericDecomposition = get_monotonic_time() - t;
78 globalStats->choleskyNNZ = static_cast<size_t>(csparse.choleskyNz());
79 }
80
81 return ok;
82 }
void prepareSolve(const SparseBlockMatrix< MatrixType > &A)
bool writeDebug() const
write a debug dump of the system matrix if it is not PSD in solve
bool solve(double *x, double *b) const
bool writeSparse(const std::string &filename) const
#define G2O_ERROR(...)
Definition logger.h:89
#define G2O_DEBUG(...)
Definition logger.h:90

References g2o::G2OBatchStatistics::choleskyNNZ, g2o::csparse::CSparse::choleskyNz(), g2o::LinearSolverCSparse< MatrixType >::csparse, G2O_DEBUG, G2O_ERROR, g2o::get_monotonic_time(), g2o::G2OBatchStatistics::globalStats(), g2o::LinearSolverCSparse< MatrixType >::prepareSolve(), g2o::csparse::CSparse::solve(), g2o::G2OBatchStatistics::timeNumericDecomposition, g2o::LinearSolver< MatrixType >::writeDebug(), and g2o::csparse::CSparse::writeSparse().

◆ solveBlocks_impl()

template<typename MatrixType >
bool g2o::LinearSolverCSparse< MatrixType >::solveBlocks_impl ( const SparseBlockMatrix< MatrixType > &  A,
std::function< void(MarginalCovarianceCholesky &)>  compute 
)
inlineprotectedvirtual

Implementation of the general parts for computing the inverse blocks of the linear system matrix. Here we call a function to do the underlying computation.

Implements g2o::LinearSolverCCS< MatrixType >.

Definition at line 172 of file linear_solver_csparse.h.

174 {
175 prepareSolve(A);
176 bool ok = csparse.factorize();
177 if (ok) {
178 csparse::CSparse::FactorView factor = csparse.factor();
179 MarginalCovarianceCholesky mcc;
180 mcc.setCholeskyFactor(factor.n, factor.p, factor.i, factor.x,
181 factor.pinv);
182 compute(mcc);
183 } else {
184 G2O_WARN("inverse fail (numeric decomposition)");
185 }
187 G2OBatchStatistics* globalStats = G2OBatchStatistics::globalStats();
188 if (globalStats) {
189 globalStats->choleskyNNZ = static_cast<size_t>(csparse.choleskyNz());
190 }
191 return ok;
192 }
#define G2O_WARN(...)
Definition logger.h:88

References g2o::G2OBatchStatistics::choleskyNNZ, g2o::csparse::CSparse::choleskyNz(), g2o::LinearSolverCSparse< MatrixType >::csparse, g2o::csparse::CSparse::factor(), g2o::csparse::CSparse::factorize(), g2o::csparse::CSparse::freeFactor(), G2O_WARN, g2o::G2OBatchStatistics::globalStats(), g2o::csparse::CSparse::FactorView::i, g2o::csparse::CSparse::FactorView::n, g2o::csparse::CSparse::FactorView::p, g2o::csparse::CSparse::FactorView::pinv, g2o::LinearSolverCSparse< MatrixType >::prepareSolve(), g2o::MarginalCovarianceCholesky::setCholeskyFactor(), and g2o::csparse::CSparse::FactorView::x.

Member Data Documentation

◆ _matrixStructure

template<typename MatrixType >
MatrixStructure g2o::LinearSolverCSparse< MatrixType >::_matrixStructure
protected

◆ csparse

template<typename MatrixType >
csparse::CSparse g2o::LinearSolverCSparse< MatrixType >::csparse
protected

The documentation for this class was generated from the following file: