10 #ifndef EIGEN_ITERATIVE_SOLVER_BASE_H
11 #define EIGEN_ITERATIVE_SOLVER_BASE_H
20 template<
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); }
52 template<
typename InputDerived>
66 template<
typename InputDerived>
70 m_preconditioner.analyzePattern(*mp_matrix);
71 m_isInitialized =
true;
72 m_analysisIsOk =
true;
86 template<
typename InputDerived>
90 eigen_assert(m_analysisIsOk &&
"You must first call analyzePattern()");
91 m_preconditioner.factorize(*mp_matrix);
92 m_factorizationIsOk =
true;
107 template<
typename InputDerived>
111 m_preconditioner.compute(*mp_matrix);
112 m_isInitialized =
true;
113 m_analysisIsOk =
true;
114 m_factorizationIsOk =
true;
120 Index rows()
const {
return mp_matrix ? mp_matrix->rows() : 0; }
122 Index cols()
const {
return mp_matrix ? mp_matrix->cols() : 0; }
143 return (mp_matrix && m_maxIterations<0) ? mp_matrix->cols() : m_maxIterations;
149 m_maxIterations = maxIters;
156 eigen_assert(m_isInitialized &&
"ConjugateGradient is not initialized.");
163 eigen_assert(m_isInitialized &&
"ConjugateGradient is not initialized.");
171 template<
typename Rhs>
inline const internal::solve_retval<Derived, Rhs>
174 eigen_assert(m_isInitialized &&
"IterativeSolverBase is not initialized.");
175 eigen_assert(rows()==b.rows()
176 &&
"IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
177 return internal::solve_retval<Derived, Rhs>(derived(), b.derived());
184 template<
typename Rhs>
185 inline const internal::sparse_solve_retval<IterativeSolverBase, Rhs>
188 eigen_assert(m_isInitialized &&
"IterativeSolverBase is not initialized.");
189 eigen_assert(rows()==b.
rows()
190 &&
"IterativeSolverBase::solve(): invalid number of rows of the right hand side matrix b");
191 return internal::sparse_solve_retval<IterativeSolverBase, Rhs>(*
this, b.
derived());
197 eigen_assert(m_isInitialized &&
"IterativeSolverBase is not initialized.");
202 template<
typename Rhs,
typename DestScalar,
int DestOptions,
typename DestIndex>
205 eigen_assert(rows()==b.rows());
207 int rhsCols = b.cols();
211 for(
int k=0; k<rhsCols; ++k)
214 tx = derived().solve(tb);
215 dest.
col(k) = tx.sparseView(0);
221 template<
typename InputDerived>
222 void grabInput(
const EigenBase<InputDerived>& A)
225 grabInput_impl(A.const_cast_derived());
228 template<
typename InputDerived>
229 void grabInput_impl(
const EigenBase<InputDerived>& A)
232 mp_matrix = &m_copyMatrix;
235 void grabInput_impl(MatrixType& A)
237 if(MatrixType::RowsAtCompileTime==
Dynamic && MatrixType::ColsAtCompileTime==
Dynamic)
238 m_copyMatrix.resize(0,0);
244 m_isInitialized =
false;
245 m_analysisIsOk =
false;
246 m_factorizationIsOk =
false;
247 m_maxIterations = -1;
248 m_tolerance = NumTraits<Scalar>::epsilon();
250 MatrixType m_copyMatrix;
251 const MatrixType* mp_matrix;
252 Preconditioner m_preconditioner;
255 RealScalar m_tolerance;
257 mutable RealScalar m_error;
258 mutable int m_iterations;
260 mutable bool m_isInitialized, m_analysisIsOk, m_factorizationIsOk;
265 template<
typename Derived,
typename Rhs>
266 struct sparse_solve_retval<IterativeSolverBase<Derived>, Rhs>
267 : sparse_solve_retval_base<IterativeSolverBase<Derived>, Rhs>
269 typedef IterativeSolverBase<Derived> Dec;
270 EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
272 template<typename Dest>
void evalTo(Dest& dst)
const
274 dec().derived()._solve_sparse(rhs(),dst);
282 #endif // EIGEN_ITERATIVE_SOLVER_BASE_H
ComputationInfo info() const
Definition: IterativeSolverBase.h:195
Preconditioner & preconditioner()
Definition: IterativeSolverBase.h:135
A versatible sparse matrix representation.
Definition: SparseMatrix.h:85
Derived & analyzePattern(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:67
const int Dynamic
Definition: Constants.h:21
Definition: EigenBase.h:26
const internal::sparse_solve_retval< IterativeSolverBase, Rhs > solve(const SparseMatrixBase< Rhs > &b) const
Definition: IterativeSolverBase.h:186
IterativeSolverBase()
Definition: IterativeSolverBase.h:36
Base class of any sparse matrices or sparse expressions.
Definition: ForwardDeclarations.h:239
Derived & derived()
Definition: EigenBase.h:34
Derived & factorize(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:87
Derived & setTolerance(const RealScalar &tolerance)
Definition: IterativeSolverBase.h:128
RealScalar error() const
Definition: IterativeSolverBase.h:161
IterativeSolverBase(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:53
RealScalar tolerance() const
Definition: IterativeSolverBase.h:125
const Preconditioner & preconditioner() const
Definition: IterativeSolverBase.h:138
Definition: Constants.h:376
const internal::solve_retval< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: IterativeSolverBase.h:172
The matrix class, also used for vectors and row-vectors.
Definition: Matrix.h:127
Index rows() const
Definition: SparseMatrixBase.h:159
ComputationInfo
Definition: Constants.h:374
Derived & setMaxIterations(int maxIters)
Definition: IterativeSolverBase.h:147
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:21
int iterations() const
Definition: IterativeSolverBase.h:154
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:48
ColXpr col(Index i)
Definition: SparseMatrixBase.h:733
int maxIterations() const
Definition: IterativeSolverBase.h:141
Derived & compute(const EigenBase< InputDerived > &A)
Definition: IterativeSolverBase.h:108