10#ifndef EIGEN_ITERATIVE_SOLVER_BASE_H
11#define EIGEN_ITERATIVE_SOLVER_BASE_H
20template<
typename Derived>
24 typedef typename internal::traits<Derived>::MatrixType MatrixType;
25 typedef typename internal::traits<Derived>::Preconditioner Preconditioner;
26 typedef typename MatrixType::Scalar Scalar;
27 typedef typename MatrixType::Index Index;
28 typedef typename MatrixType::RealScalar RealScalar;
32 Derived& derived() {
return *
static_cast<Derived*
>(
this); }
33 const Derived& derived()
const {
return *
static_cast<const Derived*
>(
this); }
68 m_isInitialized =
true;
69 m_analysisIsOk =
true;
85 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
87 m_preconditioner.factorize(A);
88 m_factorizationIsOk =
true;
107 m_isInitialized =
true;
108 m_analysisIsOk =
true;
109 m_factorizationIsOk =
true;
115 Index rows()
const {
return mp_matrix ? mp_matrix->rows() : 0; }
117 Index cols()
const {
return mp_matrix ? mp_matrix->cols() : 0; }
138 return (mp_matrix && m_maxIterations<0) ? mp_matrix->cols() : m_maxIterations;
144 m_maxIterations = maxIters;
151 eigen_assert(m_isInitialized &&
"ConjugateGradient is not initialized.");
158 eigen_assert(m_isInitialized &&
"ConjugateGradient is not initialized.");
166 template<
typename Rhs>
inline const internal::solve_retval<Derived, Rhs>
169 eigen_assert(m_isInitialized &&
"IterativeSolverBase is not initialized.");
170 eigen_assert(rows()==b.rows()
171 &&
"IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
172 return internal::solve_retval<Derived, Rhs>(derived(), b.derived());
179 template<
typename Rhs>
180 inline const internal::sparse_solve_retval<IterativeSolverBase, Rhs>
183 eigen_assert(m_isInitialized &&
"IterativeSolverBase is not initialized.");
184 eigen_assert(rows()==b.
rows()
185 &&
"IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
186 return internal::sparse_solve_retval<IterativeSolverBase, Rhs>(*
this, b.
derived());
192 eigen_assert(m_isInitialized &&
"IterativeSolverBase is not initialized.");
197 template<
typename Rhs,
typename DestScalar,
int DestOptions,
typename DestIndex>
200 eigen_assert(rows()==b.rows());
202 int rhsCols = b.cols();
206 for(
int k=0; k<rhsCols; ++k)
209 tx = derived().solve(tb);
210 dest.
col(k) = tx.sparseView(0);
217 m_isInitialized =
false;
218 m_analysisIsOk =
false;
219 m_factorizationIsOk =
false;
220 m_maxIterations = -1;
221 m_tolerance = NumTraits<Scalar>::epsilon();
223 const MatrixType* mp_matrix;
224 Preconditioner m_preconditioner;
227 RealScalar m_tolerance;
229 mutable RealScalar m_error;
230 mutable int m_iterations;
232 mutable bool m_isInitialized, m_analysisIsOk, m_factorizationIsOk;
237template<
typename Derived,
typename Rhs>
238struct sparse_solve_retval<IterativeSolverBase<Derived>, Rhs>
239 : sparse_solve_retval_base<IterativeSolverBase<Derived>, Rhs>
241 typedef IterativeSolverBase<Derived> Dec;
242 EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
244 template<
typename Dest>
void evalTo(Dest& dst)
const
246 dec().derived()._solve_sparse(rhs(),dst);
Base class for linear iterative solvers.
Definition IterativeSolverBase.h:22
Preconditioner & preconditioner()
Definition IterativeSolverBase.h:130
Derived & setMaxIterations(int maxIters)
Definition IterativeSolverBase.h:142
Derived & compute(const MatrixType &A)
Definition IterativeSolverBase.h:103
int iterations() const
Definition IterativeSolverBase.h:149
RealScalar error() const
Definition IterativeSolverBase.h:156
IterativeSolverBase()
Definition IterativeSolverBase.h:36
const internal::solve_retval< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition IterativeSolverBase.h:167
RealScalar tolerance() const
Definition IterativeSolverBase.h:120
IterativeSolverBase(const MatrixType &A)
Definition IterativeSolverBase.h:52
int maxIterations() const
Definition IterativeSolverBase.h:136
const Preconditioner & preconditioner() const
Definition IterativeSolverBase.h:133
Derived & setTolerance(RealScalar tolerance)
Definition IterativeSolverBase.h:123
const internal::sparse_solve_retval< IterativeSolverBase, Rhs > solve(const SparseMatrixBase< Rhs > &b) const
Definition IterativeSolverBase.h:181
Derived & analyzePattern(const MatrixType &A)
Definition IterativeSolverBase.h:65
ComputationInfo info() const
Definition IterativeSolverBase.h:190
Derived & factorize(const MatrixType &A)
Definition IterativeSolverBase.h:83
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:129
Base class of any sparse matrices or sparse expressions.
Definition SparseMatrixBase.h:27
SparseInnerVectorSet< Derived, 1 > col(Index j)
Definition SparseBlock.h:306
A versatible sparse matrix representation.
Definition SparseMatrix.h:87
ComputationInfo
Definition Constants.h:367
@ Success
Definition Constants.h:369
Index rows() const
Definition EigenBase.h:44
Derived & derived()
Definition EigenBase.h:34