22 #if defined(LIBMESH_HAVE_LASPACK) 87 this->_solver_type =
GMRES;
108 template <
typename T>
109 std::pair<unsigned int, Real>
114 const unsigned int m_its)
116 LOG_SCOPE(
"solve()",
"LaspackLinearSolver");
135 this->set_laspack_preconditioner_type ();
138 SetRTCAccuracy (tol);
141 switch (this->_solver_type)
146 CGIter (&matrix->
_QMat,
158 CGNIter (&matrix->
_QMat,
170 CGSIter (&matrix->
_QMat,
182 BiCGIter (&matrix->
_QMat,
194 BiCGSTABIter (&matrix->
_QMat,
206 QMRIter (&matrix->
_QMat,
218 SSORIter (&matrix->
_QMat,
230 JacobiIter (&matrix->
_QMat,
242 SetGMRESRestart (30);
243 GMRESIter (&matrix->
_QMat,
257 <<
"Continuing with GMRES" << std::endl;
259 this->_solver_type =
GMRES;
261 return this->solve (*matrix,
270 if (LASResult() != LASOK)
272 WriteLASErrDescr(stdout);
273 libmesh_error_msg(
"Exiting after LASPACK Error!");
277 return std::make_pair(GetLastNoIter(), GetLastAccuracy());
282 template <
typename T>
283 std::pair<unsigned int, Real>
288 const unsigned int m_its)
290 LOG_SCOPE(
"adjoint_solve()",
"LaspackLinearSolver");
309 this->set_laspack_preconditioner_type ();
312 SetRTCAccuracy (tol);
315 switch (this->_solver_type)
320 CGIter (Transp_Q(&matrix->
_QMat),
332 CGNIter (Transp_Q(&matrix->
_QMat),
344 CGSIter (Transp_Q(&matrix->
_QMat),
356 BiCGIter (Transp_Q(&matrix->
_QMat),
368 BiCGSTABIter (Transp_Q(&matrix->
_QMat),
380 QMRIter (Transp_Q(&matrix->
_QMat),
392 SSORIter (Transp_Q(&matrix->
_QMat),
404 JacobiIter (Transp_Q(&matrix->
_QMat),
416 SetGMRESRestart (30);
417 GMRESIter (Transp_Q(&matrix->
_QMat),
431 <<
"Continuing with GMRES" << std::endl;
433 this->_solver_type =
GMRES;
435 return this->solve (*matrix,
444 if (LASResult() != LASOK)
446 WriteLASErrDescr(stdout);
447 libmesh_error_msg(
"Exiting after LASPACK Error!");
451 return std::make_pair(GetLastNoIter(), GetLastAccuracy());
457 template <
typename T>
458 std::pair<unsigned int, Real>
465 libmesh_not_implemented();
466 return std::make_pair(0,0.0);
471 template <
typename T>
472 std::pair<unsigned int, Real>
480 libmesh_not_implemented();
481 return std::make_pair(0,0.0);
486 template <
typename T>
489 switch (this->_preconditioner_type)
492 _precond_type =
nullptr;
return;
495 _precond_type = ILUPrecond;
return;
498 _precond_type = JacobiPrecond;
return;
501 _precond_type = SSORPrecond;
return;
505 libMesh::err <<
"ERROR: Unsupported LASPACK Preconditioner: " 506 << this->_preconditioner_type << std::endl
507 <<
"Continuing with ILU" << std::endl;
509 this->set_laspack_preconditioner_type();
515 template <
typename T>
526 libMesh::out <<
"Detailed reporting is currently only supported" 527 <<
"with Petsc 2.3.1 and later." << std::endl;
532 template <
typename T>
551 #endif // #ifdef LIBMESH_HAVE_LASPACK virtual void print_converged_reason() const override
void set_laspack_preconditioner_type()
virtual void init(const char *name=nullptr) override
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
virtual void close() override
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 clear() 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 zero() override
OStreamProxy out(std::cout)
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