|
| | LeastSquaresConjugateGradient () |
| |
| template<typename MatrixDerived > |
| | LeastSquaresConjugateGradient (const EigenBase< MatrixDerived > &A) |
| |
| | ~LeastSquaresConjugateGradient () |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_vector_with_guess_impl (const Rhs &b, Dest &x) const |
| |
| | IterativeSolverBase () |
| |
| | IterativeSolverBase (const EigenBase< MatrixDerived > &A) |
| |
| | ~IterativeSolverBase () |
| |
| LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | analyzePattern (const EigenBase< MatrixDerived > &A) |
| |
| LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | factorize (const EigenBase< MatrixDerived > &A) |
| |
| LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | compute (const EigenBase< MatrixDerived > &A) |
| |
| EIGEN_CONSTEXPR Index | rows () const EIGEN_NOEXCEPT |
| |
| EIGEN_CONSTEXPR Index | cols () const EIGEN_NOEXCEPT |
| |
| RealScalar | tolerance () const |
| |
| LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | setTolerance (const RealScalar &tolerance) |
| |
| Preconditioner & | preconditioner () |
| |
| const Preconditioner & | preconditioner () const |
| |
| Index | maxIterations () const |
| |
| LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | setMaxIterations (Index maxIters) |
| |
| Index | iterations () const |
| |
| RealScalar | error () const |
| |
| const SolveWithGuess< LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >, Rhs, Guess > | solveWithGuess (const MatrixBase< Rhs > &b, const Guess &x0) const |
| |
| ComputationInfo | info () const |
| |
| void | _solve_with_guess_impl (const Rhs &b, SparseMatrixBase< DestDerived > &aDest) const |
| |
| internal::enable_if< Rhs::ColsAtCompileTime!=1 &&DestDerived::ColsAtCompileTime!=1 >::type | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &aDest) const |
| |
| internal::enable_if< Rhs::ColsAtCompileTime==1||DestDerived::ColsAtCompileTime==1 >::type | _solve_with_guess_impl (const Rhs &b, MatrixBase< DestDerived > &dest) const |
| |
| void | _solve_impl (const Rhs &b, Dest &x) const |
| |
| LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | derived () |
| |
| const LeastSquaresConjugateGradient< _MatrixType, _Preconditioner > & | derived () const |
| |
| | SparseSolverBase () |
| |
| | ~SparseSolverBase () |
| |
| Derived & | derived () |
| |
| const Derived & | derived () const |
| |
| template<typename Rhs > |
| const Solve< Derived, Rhs > | solve (const MatrixBase< Rhs > &b) const |
| |
| template<typename Rhs > |
| const Solve< Derived, Rhs > | solve (const SparseMatrixBase< Rhs > &b) const |
| |
| template<typename Rhs , typename Dest > |
| void | _solve_impl (const SparseMatrixBase< Rhs > &b, SparseMatrixBase< Dest > &dest) const |
| |
template<typename _MatrixType, typename _Preconditioner>
class Eigen::LeastSquaresConjugateGradient< _MatrixType, _Preconditioner >
A conjugate gradient solver for sparse (or dense) least-square problems.
This class allows to solve for A x = b linear problems using an iterative conjugate gradient algorithm. The matrix A can be non symmetric and rectangular, but the matrix A' A should be positive-definite to guaranty stability. Otherwise, the SparseLU or SparseQR classes might be preferable. The matrix A and the vectors x and b can be either dense or sparse.
- Template Parameters
-
\implsparsesolverconcept
The maximal number of iterations and tolerance value can be controlled via the setMaxIterations() and setTolerance() methods. The defaults are the size of the problem for the maximal number of iterations and NumTraits<Scalar>::epsilon() for the tolerance.
This class can be used as the direct solver classes. Here is a typical usage example:
int m=1000000,
n = 10000;
std::cout <<
"#iterations: " << lscg.
iterations() << std::endl;
std::cout <<
"estimated error: " << lscg.
error() << std::endl;
Matrix3f m
Definition AngleAxis_mimic_euler.cpp:1
int n
Definition BiCGSTAB_simple.cpp:1
Scalar * b
Definition benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition bench_gemm.cpp:48
RealScalar error() const
Definition IterativeSolverBase.h:305
Derived & compute(const EigenBase< MatrixDerived > &A)
Definition IterativeSolverBase.h:238
Index iterations() const
Definition IterativeSolverBase.h:296
A conjugate gradient solver for sparse (or dense) least-square problems.
Definition LeastSquareConjugateGradient.h:150
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition SparseSolverBase.h:88
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Definition gnuplot_common_settings.hh:12
By default the iterations start with x=0 as an initial guess of the solution. One can control the start using the solveWithGuess() method.
- See also
- class ConjugateGradient, SparseLU, SparseQR