20 #include "libmesh/libmesh_common.h"    22 #if defined(LIBMESH_HAVE_LASPACK)    26 #include "libmesh/laspack_linear_solver.h"    27 #include "libmesh/libmesh_logging.h"    28 #include "libmesh/enum_to_string.h"    29 #include "libmesh/enum_solver_type.h"    30 #include "libmesh/enum_preconditioner_type.h"    31 #include "libmesh/enum_convergence_flags.h"    87       this->_solver_type         = 
GMRES;
   108 template <
typename T>
   109 std::pair<unsigned int, Real>
   113                                const std::optional<double> tol,
   114                                const std::optional<unsigned int> m_its)
   116   LOG_SCOPE(
"solve()", 
"LaspackLinearSolver");
   124   const int max_its = this->get_int_solver_setting(
"max_its", m_its);
   125   const double abs_tol = this->get_real_solver_setting(
"abs_tol", tol);
   138   this->set_laspack_preconditioner_type ();
   141   SetRTCAccuracy (abs_tol);
   144   switch (this->_solver_type)
   149         CGIter (&matrix->
_QMat,
   161         CGNIter (&matrix->
_QMat,
   173         CGSIter (&matrix->
_QMat,
   185         BiCGIter (&matrix->
_QMat,
   197         BiCGSTABIter (&matrix->
_QMat,
   209         QMRIter (&matrix->
_QMat,
   221         SSORIter (&matrix->
_QMat,
   233         JacobiIter (&matrix->
_QMat,
   245         SetGMRESRestart (30);
   246         GMRESIter (&matrix->
_QMat,
   260                      << 
"Continuing with GMRES" << std::endl;
   262         this->_solver_type = 
GMRES;
   264         return this->solve (*matrix,
   273   if (LASResult() != LASOK)
   275       WriteLASErrDescr(stdout);
   276       libmesh_error_msg(
"Exiting after LASPACK Error!");
   280   return std::make_pair(GetLastNoIter(), GetLastAccuracy());
   285 template <
typename T>
   286 std::pair<unsigned int, Real>
   290                                        const std::optional<double> tol,
   291                                        const std::optional<unsigned int> m_its)
   293   LOG_SCOPE(
"adjoint_solve()", 
"LaspackLinearSolver");
   301   const int max_its = this->get_int_solver_setting(
"max_its", m_its);
   302   const double abs_tol = this->get_real_solver_setting(
"abs_tol", tol);
   315   this->set_laspack_preconditioner_type ();
   318   SetRTCAccuracy (abs_tol);
   321   switch (this->_solver_type)
   326         CGIter (Transp_Q(&matrix->
_QMat),
   338         CGNIter (Transp_Q(&matrix->
_QMat),
   350         CGSIter (Transp_Q(&matrix->
_QMat),
   362         BiCGIter (Transp_Q(&matrix->
_QMat),
   374         BiCGSTABIter (Transp_Q(&matrix->
_QMat),
   386         QMRIter (Transp_Q(&matrix->
_QMat),
   398         SSORIter (Transp_Q(&matrix->
_QMat),
   410         JacobiIter (Transp_Q(&matrix->
_QMat),
   422         SetGMRESRestart (30);
   423         GMRESIter (Transp_Q(&matrix->
_QMat),
   437                      << 
"Continuing with GMRES" << std::endl;
   439         this->_solver_type = 
GMRES;
   441         return this->solve (*matrix,
   450   if (LASResult() != LASOK)
   452       WriteLASErrDescr(stdout);
   453       libmesh_error_msg(
"Exiting after LASPACK Error!");
   457   return std::make_pair(GetLastNoIter(), GetLastAccuracy());
   463 template <
typename T>
   464 std::pair<unsigned int, Real>
   468                                const std::optional<double> ,
   469                                const std::optional<unsigned int> )
   471   libmesh_not_implemented();
   472   return std::make_pair(0,0.0);
   477 template <
typename T>
   478 std::pair<unsigned int, Real>
   483                                const std::optional<double> ,
   484                                const std::optional<unsigned int> )
   486   libmesh_not_implemented();
   487   return std::make_pair(0,0.0);
   492 template <
typename T>
   495   switch (this->_preconditioner_type)
   498       _precond_type = 
nullptr; 
return;
   501       _precond_type = ILUPrecond; 
return;
   504       _precond_type = JacobiPrecond; 
return;
   507       _precond_type = SSORPrecond; 
return;
   511       libMesh::err << 
"ERROR:  Unsupported LASPACK Preconditioner: "   512                    << this->_preconditioner_type << std::endl
   513                    << 
"Continuing with ILU"      << std::endl;
   515       this->set_laspack_preconditioner_type();
   521 template <
typename T>
   532   libMesh::out << 
"Detailed reporting is currently only supported"   533                << 
"with Petsc 2.3.1 and later." << std::endl;
   538 template <
typename T>
   557 #endif // #ifdef LIBMESH_HAVE_LASPACK 
virtual void print_converged_reason() const override
Prints a useful message about why the latest linear solve con(di)verged. 
 
void set_laspack_preconditioner_type()
Tells LASPACK to use the user-specified preconditioner stored in _preconditioner_type. 
 
LinearConvergenceReason
Linear solver convergence flags (taken from the PETSc flags). 
 
virtual void init(const char *name=nullptr) override
Initialize data structures if not done so already. 
 
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
 
QMatrix _QMat
The Laspack sparse matrix pointer. 
 
The libMesh namespace provides an interface to certain functionality in the library. 
 
This class provides an interface to Laspack iterative solvers that is compatible with the libMesh Lin...
 
virtual void close() override
Calls the SparseMatrix's internal assembly routines, ensuring that the values are consistent across p...
 
bool _is_initialized
Flag that tells if init() has been called. 
 
virtual void clear() override
Release all memory and clear data structures. 
 
void init(triangulateio &t)
Initializes the fields of t to nullptr/0 as necessary. 
 
virtual LinearConvergenceReason get_converged_reason() const override
 
std::string enum_to_string(const T e)
 
This class provides a nice interface to the Laspack C-based data structures for serial vectors...
 
bool initialized()
Checks that library initialization has been done. 
 
virtual std::pair< unsigned int, Real > adjoint_solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Laspack solver to solve A^T x = b. 
 
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const std::optional< double > tol=std::nullopt, const std::optional< unsigned int > m_its=std::nullopt) override
Call the Laspack solver. 
 
The LaspackMatrix class wraps a QMatrix object from the Laspack library. 
 
Generic shell matrix, i.e.