22 #ifdef LIBMESH_HAVE_EIGEN 35 #include <unsupported/Eigen/IterativeSolvers> 45 _comp_info(Eigen::Success)
82 std::pair<unsigned int, Real>
87 const unsigned int m_its)
89 LOG_SCOPE(
"solve()",
"EigenSparseLinearSolver");
102 std::pair<unsigned int, Real> retval(0,0.);
105 switch (this->_solver_type)
110 Eigen::ConjugateGradient<EigenSM> solver (matrix.
_mat);
111 solver.setMaxIterations(m_its);
112 solver.setTolerance(tol);
113 solution._vec = solver.solveWithGuess(rhs._vec,solution._vec);
114 libMesh::out <<
"#iterations: " << solver.iterations() << std::endl;
115 libMesh::out <<
"estimated error: " << solver.error() << std::endl;
116 retval = std::make_pair(solver.iterations(), solver.error());
117 _comp_info = solver.info();
124 Eigen::BiCGSTAB<EigenSM> solver (matrix.
_mat);
125 solver.setMaxIterations(m_its);
126 solver.setTolerance(tol);
127 solution._vec = solver.solveWithGuess(rhs._vec,solution._vec);
128 libMesh::out <<
"#iterations: " << solver.iterations() << std::endl;
129 libMesh::out <<
"estimated error: " << solver.error() << std::endl;
130 retval = std::make_pair(solver.iterations(), solver.error());
131 _comp_info = solver.info();
138 Eigen::GMRES<EigenSM> solver (matrix.
_mat);
139 solver.setMaxIterations(m_its);
140 solver.setTolerance(tol);
145 if (this->_solver_configuration)
147 auto it = this->_solver_configuration->int_valued_data.find(
"gmres_restart");
149 if (it != this->_solver_configuration->int_valued_data.end())
150 solver.set_restart(it->second);
153 libMesh::out <<
"Eigen GMRES solver, restart = " << solver.get_restart() << std::endl;
154 solution._vec = solver.solveWithGuess(rhs._vec, solution._vec);
155 libMesh::out <<
"#iterations: " << solver.iterations() << std::endl;
156 libMesh::out <<
"estimated error: " << solver.error() << std::endl;
157 retval = std::make_pair(solver.iterations(), solver.error());
158 _comp_info = solver.info();
177 matrix.
_mat.makeCompressed();
188 Eigen::SparseLU<EigenSM> solver;
191 solver.analyzePattern(matrix.
_mat);
194 solver.factorize(matrix.
_mat);
197 solution._vec = solver.solve(rhs._vec);
202 retval = std::make_pair(1, 0);
205 _comp_info = solver.info();
214 <<
"Continuing with BICGSTAB" << std::endl;
218 return this->solve (matrix,
231 template <
typename T>
232 std::pair<unsigned int, Real>
237 const unsigned int m_its)
239 LOG_SCOPE(
"adjoint_solve()",
"EigenSparseLinearSolver");
241 libmesh_experimental();
245 std::pair<unsigned int, Real> retval = this->solve (mat_trans,
257 template <
typename T>
258 std::pair<unsigned int, Real>
265 libmesh_not_implemented();
266 return std::make_pair(0,0.0);
271 template <
typename T>
272 std::pair<unsigned int, Real>
280 libmesh_not_implemented();
281 return std::make_pair(0,0.0);
286 template <
typename T>
289 libmesh_not_implemented();
317 template <
typename T>
320 auto it = _convergence_reasons.find(_comp_info);
324 if (it == _convergence_reasons.end())
326 libmesh_warning(
"Warning: unknown Eigen::ComputationInfo: " \
328 <<
" returning CONVERGED_ITS." \
345 #endif // #ifdef LIBMESH_HAVE_EIGEN
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const double tol, const unsigned int m_its) override
virtual void close() override
EigenSparseLinearSolver(const libMesh::Parallel::Communicator &comm_in)
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const double tol, const unsigned int m_its) override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
virtual void init(const char *name=nullptr) override
void init(triangulateio &t)
virtual LinearConvergenceReason get_converged_reason() const override
OStreamProxy err(std::cerr)
std::string enum_to_string(const T e)
virtual void clear() override
void set_eigen_preconditioner_type()
virtual void get_transpose(SparseMatrix< T > &dest) const =0
OStreamProxy out(std::cout)