libMesh::LaspackLinearSolver< T > Class Template Reference

#include <laspack_linear_solver.h>

Inheritance diagram for libMesh::LaspackLinearSolver< T >:

List of all members.

Public Member Functions

 LaspackLinearSolver ()
 ~LaspackLinearSolver ()
void clear ()
void init ()
std::pair< unsigned int, Realsolve (SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const double tol, const unsigned int m_its)
std::pair< unsigned int, Realadjoint_solve (SparseMatrix< T > &matrix, NumericVector< T > &solution, NumericVector< T > &rhs, const double tol, const unsigned int m_its)
std::pair< unsigned int, Realsolve (SparseMatrix< T > &matrix, SparseMatrix< T > &pc, NumericVector< T > &solution, NumericVector< T > &rhs, const double tol, const unsigned int m_its)
std::pair< unsigned int, Realsolve (const ShellMatrix< T > &shell_matrix, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its)
virtual std::pair< unsigned
int, Real
solve (const ShellMatrix< T > &shell_matrix, const SparseMatrix< T > &precond_matrix, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its)
virtual void print_converged_reason ()
bool initialized () const
SolverType solver_type () const
void set_solver_type (const SolverType st)
PreconditionerType preconditioner_type () const
void set_preconditioner_type (const PreconditionerType pct)
void attach_preconditioner (Preconditioner< T > *preconditioner)
virtual void reuse_preconditioner (bool)
bool get_same_preconditioner ()
virtual void restrict_solve_to (const std::vector< unsigned int > *const dofs, const SubsetSolveMode subset_solve_mode=SUBSET_ZERO)
std::pair< unsigned int, Realsolve (SparseMatrix< T > &matrix, SparseMatrix< T > *precond_matrix, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)
std::pair< unsigned int, Realsolve (const ShellMatrix< T > &matrix, const SparseMatrix< T > *precond_matrix, NumericVector< T > &, NumericVector< T > &, const double, const unsigned int)

Static Public Member Functions

static AutoPtr< LinearSolver< T > > build (const SolverPackage solver_package=libMesh::default_solver_package())
static std::string get_info ()
static void print_info (std::ostream &out=libMesh::out)
static unsigned int n_objects ()
static void enable_print_counter_info ()
static void disable_print_counter_info ()

Protected Types

typedef std::map< std::string,
std::pair< unsigned int,
unsigned int > > 
Counts

Protected Member Functions

void increment_constructor_count (const std::string &name)
void increment_destructor_count (const std::string &name)

Protected Attributes

SolverType _solver_type
PreconditionerType _preconditioner_type
bool _is_initialized
Preconditioner< T > * _preconditioner
bool same_preconditioner

Static Protected Attributes

static Counts _counts
static Threads::atomic
< unsigned int > 
_n_objects
static Threads::spin_mutex _mutex
static bool _enable_print_counter = true

Private Member Functions

void set_laspack_preconditioner_type ()

Private Attributes

PrecondProcType _precond_type

Detailed Description

template<typename T>
class libMesh::LaspackLinearSolver< T >

This class provides an interface to Laspack iterative solvers that is compatible with the libMesh LinearSolver<>

Author:
Benjamin Kirk, 2002-2007

Definition at line 53 of file laspack_linear_solver.h.


Member Typedef Documentation

typedef std::map<std::string, std::pair<unsigned int, unsigned int> > libMesh::ReferenceCounter::Counts [protected, inherited]

Data structure to log the information. The log is identified by the class name.

Definition at line 113 of file reference_counter.h.


Constructor & Destructor Documentation

template<typename T >
libMesh::LaspackLinearSolver< T >::LaspackLinearSolver (  )  [inline]

Constructor. Initializes Laspack data structures

Definition at line 154 of file laspack_linear_solver.h.

00154                                              :
00155   _precond_type (ILUPrecond)
00156 {
00157 }

template<typename T >
libMesh::LaspackLinearSolver< T >::~LaspackLinearSolver (  )  [inline]

Destructor.

Definition at line 163 of file laspack_linear_solver.h.

References libMesh::LaspackLinearSolver< T >::clear().

00164 {
00165   this->clear ();
00166 }


Member Function Documentation

template<typename T >
std::pair< unsigned int, Real > libMesh::LaspackLinearSolver< T >::adjoint_solve ( SparseMatrix< T > &  matrix,
NumericVector< T > &  solution,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  m_its 
) [inline, virtual]

Call the Laspack solver to solve A^T x = b

Reimplemented from libMesh::LinearSolver< T >.

Definition at line 284 of file laspack_linear_solver.C.

References libMesh::LaspackLinearSolver< T >::_precond_type, libMesh::LaspackMatrix< T >::_QMat, libMesh::LinearSolver< T >::_solver_type, libMeshEnums::BICG, libMeshEnums::BICGSTAB, libMeshEnums::CG, libMeshEnums::CGN, libMeshEnums::CGS, libMesh::LaspackMatrix< T >::close(), libMesh::err, libMeshEnums::GMRES, libMesh::LaspackLinearSolver< T >::init(), libMeshEnums::JACOBI, libMeshEnums::QMR, libMesh::LaspackLinearSolver< T >::set_laspack_preconditioner_type(), libMesh::LaspackLinearSolver< T >::solve(), and libMeshEnums::SSOR.

00289 {
00290   START_LOG("adjoint_solve()", "LaspackLinearSolver");
00291   this->init ();
00292 
00293   // Make sure the data passed in are really in Laspack types
00294   LaspackMatrix<T>* matrix   = libmesh_cast_ptr<LaspackMatrix<T>*>(&matrix_in);
00295   LaspackVector<T>* solution = libmesh_cast_ptr<LaspackVector<T>*>(&solution_in);
00296   LaspackVector<T>* rhs      = libmesh_cast_ptr<LaspackVector<T>*>(&rhs_in);
00297 
00298   // Zero-out the solution to prevent the solver from exiting in 0
00299   // iterations (?)
00300   //TODO:[BSK] Why does Laspack do this?  Comment out this and try ex13...
00301   solution->zero();
00302 
00303   // Close the matrix and vectors in case this wasn't already done.
00304   matrix->close ();
00305   solution->close ();
00306   rhs->close ();
00307 
00308   // Set the preconditioner type
00309   this->set_laspack_preconditioner_type ();
00310 
00311   // Set the solver tolerance
00312   SetRTCAccuracy (tol);
00313 
00314   // Solve the linear system
00315   switch (this->_solver_type)
00316     {
00317       // Conjugate-Gradient
00318     case CG:
00319       {
00320         CGIter (Transp_Q(&matrix->_QMat),
00321                 &solution->_vec,
00322                 &rhs->_vec,
00323                 m_its,
00324                 _precond_type,
00325                 1.);
00326         break;
00327       }
00328 
00329       // Conjugate-Gradient Normalized
00330     case CGN:
00331       {
00332         CGNIter (Transp_Q(&matrix->_QMat),
00333                  &solution->_vec,
00334                  &rhs->_vec,
00335                  m_its,
00336                  _precond_type,
00337                  1.);
00338         break;
00339       }
00340 
00341       // Conjugate-Gradient Squared
00342     case CGS:
00343       {
00344         CGSIter (Transp_Q(&matrix->_QMat),
00345                  &solution->_vec,
00346                  &rhs->_vec,
00347                  m_its,
00348                  _precond_type,
00349                  1.);
00350         break;
00351       }
00352 
00353       // Bi-Conjugate Gradient
00354     case BICG:
00355       {
00356         BiCGIter (Transp_Q(&matrix->_QMat),
00357                   &solution->_vec,
00358                   &rhs->_vec,
00359                   m_its,
00360                   _precond_type,
00361                   1.);
00362         break;
00363       }
00364 
00365       // Bi-Conjugate Gradient Stabilized
00366     case BICGSTAB:
00367       {
00368         BiCGSTABIter (Transp_Q(&matrix->_QMat),
00369                       &solution->_vec,
00370                       &rhs->_vec,
00371                       m_its,
00372                       _precond_type,
00373                       1.);
00374         break;
00375       }
00376 
00377       // Quasi-Minimum Residual
00378     case QMR:
00379       {
00380         QMRIter (Transp_Q(&matrix->_QMat),
00381                  &solution->_vec,
00382                  &rhs->_vec,
00383                  m_its,
00384                  _precond_type,
00385                  1.);
00386         break;
00387       }
00388 
00389       // Symmetric over-relaxation
00390     case SSOR:
00391       {
00392         SSORIter (Transp_Q(&matrix->_QMat),
00393                   &solution->_vec,
00394                   &rhs->_vec,
00395                   m_its,
00396                   _precond_type,
00397                   1.);
00398         break;
00399       }
00400 
00401       // Jacobi Relaxation
00402     case JACOBI:
00403       {
00404         JacobiIter (Transp_Q(&matrix->_QMat),
00405                     &solution->_vec,
00406                     &rhs->_vec,
00407                     m_its,
00408                     _precond_type,
00409                     1.);
00410         break;
00411       }
00412 
00413       // Generalized Minimum Residual
00414     case GMRES:
00415       {
00416         SetGMRESRestart (30);
00417         GMRESIter (Transp_Q(&matrix->_QMat),
00418                    &solution->_vec,
00419                    &rhs->_vec,
00420                    m_its,
00421                    _precond_type,
00422                    1.);
00423         break;
00424       }
00425 
00426       // Unknown solver, use GMRES
00427     default:
00428       {
00429         libMesh::err << "ERROR:  Unsupported LASPACK Solver: "
00430                       << this->_solver_type      << std::endl
00431                       << "Continuing with GMRES" << std::endl;
00432 
00433         this->_solver_type = GMRES;
00434 
00435         return this->solve (*matrix,
00436                             *solution,
00437                             *rhs,
00438                             tol,
00439                             m_its);
00440       }
00441     }
00442 
00443   // Check for an error
00444   if (LASResult() != LASOK)
00445     {
00446       libMesh::err << "ERROR:  LASPACK Error: " << std::endl;
00447       WriteLASErrDescr(stdout);
00448       libmesh_error();
00449     }
00450 
00451   STOP_LOG("adjoint_solve()", "LaspackLinearSolver");
00452   // Get the convergence step # and residual
00453   return std::make_pair(GetLastNoIter(), GetLastAccuracy());
00454 }

template<typename T>
void libMesh::LinearSolver< T >::attach_preconditioner ( Preconditioner< T > *  preconditioner  )  [inline, inherited]

Attaches a Preconditioner object to be used

Definition at line 105 of file linear_solver.C.

References libMesh::LinearSolver< T >::_is_initialized, libMesh::LinearSolver< T >::_preconditioner, libMesh::LinearSolver< T >::_preconditioner_type, libMesh::err, and libMeshEnums::SHELL_PRECOND.

00106 {
00107   if(this->_is_initialized)
00108   {
00109     libMesh::err<<"Preconditioner must be attached before the solver is initialized!"<<std::endl;
00110     libmesh_error();
00111   }
00112 
00113   _preconditioner_type = SHELL_PRECOND;
00114   _preconditioner = preconditioner;
00115 }

template<typename T >
AutoPtr< LinearSolver< T > > libMesh::LinearSolver< T >::build ( const SolverPackage  solver_package = libMesh::default_solver_package()  )  [inline, static, inherited]

Builds a LinearSolver using the linear solver package specified by solver_package

Definition at line 39 of file linear_solver.C.

References libMesh::err, libMesh::LASPACK_SOLVERS, libMeshEnums::PETSC_SOLVERS, and libMesh::TRILINOS_SOLVERS.

00040 {
00041   // Build the appropriate solver
00042   switch (solver_package)
00043     {
00044 
00045 
00046 #ifdef LIBMESH_HAVE_LASPACK
00047     case LASPACK_SOLVERS:
00048       {
00049         AutoPtr<LinearSolver<T> > ap(new LaspackLinearSolver<T>);
00050         return ap;
00051       }
00052 #endif
00053 
00054 
00055 #ifdef LIBMESH_HAVE_PETSC
00056     case PETSC_SOLVERS:
00057       {
00058         AutoPtr<LinearSolver<T> > ap(new PetscLinearSolver<T>);
00059         return ap;
00060       }
00061 #endif
00062 
00063 
00064 #ifdef LIBMESH_HAVE_TRILINOS
00065     case TRILINOS_SOLVERS:
00066       {
00067         AutoPtr<LinearSolver<T> > ap(new AztecLinearSolver<T>);
00068         return ap;
00069       }
00070 #endif
00071 
00072     default:
00073       libMesh::err << "ERROR:  Unrecognized solver package: "
00074                     << solver_package
00075                     << std::endl;
00076       libmesh_error();
00077     }
00078 
00079   AutoPtr<LinearSolver<T> > ap(NULL);
00080   return ap;
00081 }

template<typename T >
void libMesh::LaspackLinearSolver< T >::clear (  )  [inline, virtual]
void libMesh::ReferenceCounter::disable_print_counter_info (  )  [static, inherited]

Definition at line 106 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

00107 {
00108   _enable_print_counter = false;
00109   return;
00110 }

void libMesh::ReferenceCounter::enable_print_counter_info (  )  [static, inherited]

Methods to enable/disable the reference counter output from print_info()

Definition at line 100 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter.

00101 {
00102   _enable_print_counter = true;
00103   return;
00104 }

std::string libMesh::ReferenceCounter::get_info (  )  [static, inherited]

Gets a string containing the reference information.

Definition at line 47 of file reference_counter.C.

References libMesh::ReferenceCounter::_counts, and libMesh::Quality::name().

Referenced by libMesh::ReferenceCounter::print_info().

00048 {
00049 #if defined(LIBMESH_ENABLE_REFERENCE_COUNTING) && defined(DEBUG)
00050 
00051   std::ostringstream oss;
00052 
00053   oss << '\n'
00054       << " ---------------------------------------------------------------------------- \n"
00055       << "| Reference count information                                                |\n"
00056       << " ---------------------------------------------------------------------------- \n";
00057 
00058   for (Counts::iterator it = _counts.begin();
00059        it != _counts.end(); ++it)
00060     {
00061       const std::string name(it->first);
00062       const unsigned int creations    = it->second.first;
00063       const unsigned int destructions = it->second.second;
00064 
00065       oss << "| " << name << " reference count information:\n"
00066           << "|  Creations:    " << creations    << '\n'
00067           << "|  Destructions: " << destructions << '\n';
00068     }
00069 
00070   oss << " ---------------------------------------------------------------------------- \n";
00071 
00072   return oss.str();
00073 
00074 #else
00075 
00076   return "";
00077 
00078 #endif
00079 }

template<typename T >
bool libMesh::LinearSolver< T >::get_same_preconditioner (  )  [inline, inherited]

Definition at line 285 of file linear_solver.h.

References libMesh::LinearSolver< T >::same_preconditioner.

00286 {
00287   return same_preconditioner;
00288 }

void libMesh::ReferenceCounter::increment_constructor_count ( const std::string &  name  )  [inline, protected, inherited]

Increments the construction counter. Should be called in the constructor of any derived class that will be reference counted.

Definition at line 163 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::ReferenceCountedObject().

00164 {
00165   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00166   std::pair<unsigned int, unsigned int>& p = _counts[name];
00167 
00168   p.first++;
00169 }

void libMesh::ReferenceCounter::increment_destructor_count ( const std::string &  name  )  [inline, protected, inherited]

Increments the destruction counter. Should be called in the destructor of any derived class that will be reference counted.

Definition at line 176 of file reference_counter.h.

References libMesh::ReferenceCounter::_counts, and libMesh::Threads::spin_mtx.

Referenced by libMesh::ReferenceCountedObject< RBParametrized >::~ReferenceCountedObject().

00177 {
00178   Threads::spin_mutex::scoped_lock lock(Threads::spin_mtx);
00179   std::pair<unsigned int, unsigned int>& p = _counts[name];
00180 
00181   p.second++;
00182 }

template<typename T >
void libMesh::LaspackLinearSolver< T >::init (  )  [inline, virtual]

Initialize data structures if not done so already.

Implements libMesh::LinearSolver< T >.

Definition at line 93 of file laspack_linear_solver.C.

References libMesh::LinearSolver< T >::_is_initialized, and libMesh::LinearSolver< T >::initialized().

Referenced by libMesh::LaspackLinearSolver< T >::adjoint_solve(), and libMesh::LaspackLinearSolver< T >::solve().

00094 {
00095   // Initialize the data structures if not done so already.
00096   if (!this->initialized())
00097     {
00098       this->_is_initialized = true;
00099     }
00100 
00101  // SetRTCAuxProc (print_iter_accuracy);
00102 }

template<typename T>
bool libMesh::LinearSolver< T >::initialized (  )  const [inline, inherited]
static unsigned int libMesh::ReferenceCounter::n_objects (  )  [inline, static, inherited]

Prints the number of outstanding (created, but not yet destroyed) objects.

Definition at line 79 of file reference_counter.h.

References libMesh::ReferenceCounter::_n_objects.

00080   { return _n_objects; }

template<typename T >
PreconditionerType libMesh::LinearSolver< T >::preconditioner_type (  )  const [inline, inherited]

Returns the type of preconditioner to use.

Definition at line 85 of file linear_solver.C.

References libMesh::LinearSolver< T >::_preconditioner, and libMesh::LinearSolver< T >::_preconditioner_type.

00086 {
00087   if(_preconditioner)
00088     return _preconditioner->type();
00089 
00090   return _preconditioner_type;
00091 }

template<typename T >
void libMesh::LaspackLinearSolver< T >::print_converged_reason (  )  [inline, virtual]

Prints a useful message about why the latest linear solve con(di)verged.

Implements libMesh::LinearSolver< T >.

Definition at line 518 of file laspack_linear_solver.C.

References libMesh::out.

00519 {
00520   libMesh::out << "print_converged_reason() is currently only supported"
00521                 << "with Petsc 2.3.1 and later." << std::endl;
00522 }

void libMesh::ReferenceCounter::print_info ( std::ostream &  out = libMesh::out  )  [static, inherited]

Prints the reference information, by default to libMesh::out.

Definition at line 88 of file reference_counter.C.

References libMesh::ReferenceCounter::_enable_print_counter, and libMesh::ReferenceCounter::get_info().

00089 {
00090   if( _enable_print_counter ) out_stream << ReferenceCounter::get_info();
00091 }

template<typename T >
void libMesh::LinearSolver< T >::restrict_solve_to ( const std::vector< unsigned int > *const   dofs,
const SubsetSolveMode  subset_solve_mode = SUBSET_ZERO 
) [inline, virtual, inherited]

After calling this method, all successive solves will be restricted to the given set of dofs, which must contain local dofs on each processor only and not contain any duplicates. This mode can be disabled by calling this method with dofs being a NULL pointer.

Reimplemented in libMesh::PetscLinearSolver< T >.

Definition at line 126 of file linear_solver.C.

00128 {
00129   if(dofs!=NULL)
00130     {
00131       libmesh_not_implemented();
00132     }
00133 }

template<typename T >
void libMesh::LinearSolver< T >::reuse_preconditioner ( bool  reuse_flag  )  [inline, virtual, inherited]

Definition at line 119 of file linear_solver.C.

References libMesh::LinearSolver< T >::same_preconditioner.

Referenced by libMesh::ImplicitSystem::disable_cache().

00120   {
00121     same_preconditioner = reuse_flag;
00122   }

template<typename T >
void libMesh::LaspackLinearSolver< T >::set_laspack_preconditioner_type (  )  [inline, private]

Tells LASPACK to use the user-specified preconditioner stored in _preconditioner_type

Definition at line 489 of file laspack_linear_solver.C.

References libMesh::LaspackLinearSolver< T >::_precond_type, libMesh::LinearSolver< T >::_preconditioner_type, libMesh::err, libMeshEnums::IDENTITY_PRECOND, libMeshEnums::ILU_PRECOND, libMeshEnums::JACOBI_PRECOND, and libMeshEnums::SSOR_PRECOND.

Referenced by libMesh::LaspackLinearSolver< T >::adjoint_solve(), and libMesh::LaspackLinearSolver< T >::solve().

00490 {
00491   switch (this->_preconditioner_type)
00492     {
00493     case IDENTITY_PRECOND:
00494       _precond_type = NULL; return;
00495 
00496     case ILU_PRECOND:
00497       _precond_type = ILUPrecond; return;
00498 
00499     case JACOBI_PRECOND:
00500       _precond_type = JacobiPrecond; return;
00501 
00502     case SSOR_PRECOND:
00503       _precond_type = SSORPrecond; return;
00504 
00505 
00506     default:
00507       libMesh::err << "ERROR:  Unsupported LASPACK Preconditioner: "
00508                     << this->_preconditioner_type << std::endl
00509                     << "Continuing with ILU"      << std::endl;
00510       this->_preconditioner_type = ILU_PRECOND;
00511       this->set_laspack_preconditioner_type();
00512     }
00513 }

template<typename T >
void libMesh::LinearSolver< T >::set_preconditioner_type ( const PreconditionerType  pct  )  [inline, inherited]

Sets the type of preconditioner to use.

Definition at line 95 of file linear_solver.C.

References libMesh::LinearSolver< T >::_preconditioner, and libMesh::LinearSolver< T >::_preconditioner_type.

00096 {
00097   if(_preconditioner)
00098     _preconditioner->set_type(pct);
00099   else
00100     _preconditioner_type = pct;
00101 }

template<typename T>
void libMesh::LinearSolver< T >::set_solver_type ( const SolverType  st  )  [inline, inherited]

Sets the type of solver to use.

Definition at line 101 of file linear_solver.h.

00102   { _solver_type = st; }

template<typename T>
std::pair< unsigned int, Real > libMesh::LinearSolver< T >::solve ( const ShellMatrix< T > &  matrix,
const SparseMatrix< T > *  precond_matrix,
NumericVector< T > &  sol,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  n_iter 
) [inline, inherited]

This function solves a system whose matrix is a shell matrix, but an optional sparse matrix may be used as preconditioning matrix.

Definition at line 310 of file linear_solver.h.

References libMesh::LinearSolver< T >::solve().

00316 {
00317   if (pc_mat)
00318     return this->solve(mat, *pc_mat, sol, rhs, tol, n_iter);
00319   else
00320     return this->solve(mat, sol, rhs, tol, n_iter);
00321 }

template<typename T>
std::pair< unsigned int, Real > libMesh::LinearSolver< T >::solve ( SparseMatrix< T > &  matrix,
SparseMatrix< T > *  precond_matrix,
NumericVector< T > &  sol,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  n_iter 
) [inline, inherited]

This function calls the solver "_solver_type" preconditioned with the "_preconditioner_type" preconditioner. The preconditioning matrix is used if it is provided, or the system matrix is used if precond_matrix is null

Definition at line 293 of file linear_solver.h.

References libMesh::LinearSolver< T >::solve().

00299 {
00300   if (pc_mat)
00301     return this->solve(mat, *pc_mat, sol, rhs, tol, n_iter);
00302   else
00303     return this->solve(mat, sol, rhs, tol, n_iter);
00304 }

template<typename T >
std::pair< unsigned int, Real > libMesh::LaspackLinearSolver< T >::solve ( const ShellMatrix< T > &  shell_matrix,
const SparseMatrix< T > &  precond_matrix,
NumericVector< T > &  solution_in,
NumericVector< T > &  rhs_in,
const double  tol,
const unsigned int  m_its 
) [inline, virtual]

This function solves a system whose matrix is a shell matrix, but a sparse matrix is used as preconditioning matrix, this allowing other preconditioners than JACOBI.

Implements libMesh::LinearSolver< T >.

Definition at line 475 of file laspack_linear_solver.C.

00481 {
00482   libmesh_not_implemented();
00483   return std::make_pair(0,0.0);
00484 }

template<typename T >
std::pair< unsigned int, Real > libMesh::LaspackLinearSolver< T >::solve ( const ShellMatrix< T > &  shell_matrix,
NumericVector< T > &  solution_in,
NumericVector< T > &  rhs_in,
const double  tol,
const unsigned int  m_its 
) [inline, virtual]

This function solves a system whose matrix is a shell matrix.

Implements libMesh::LinearSolver< T >.

Definition at line 461 of file laspack_linear_solver.C.

00466 {
00467   libmesh_not_implemented();
00468   return std::make_pair(0,0.0);
00469 }

template<typename T >
std::pair< unsigned int, Real > libMesh::LaspackLinearSolver< T >::solve ( SparseMatrix< T > &  matrix,
SparseMatrix< T > &  pc,
NumericVector< T > &  solution,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  m_its 
) [inline, virtual]

Call the Laspack solver

Implements libMesh::LinearSolver< T >.

Definition at line 173 of file laspack_linear_solver.h.

References libMesh::err.

00179 {
00180   libMesh::err << "ERROR: LASPACK does not support a user-supplied preconditioner!"
00181                 << std::endl;
00182   libmesh_error();
00183 
00184   std::pair<unsigned int, Real> p;
00185   return p;
00186 }

template<typename T >
std::pair< unsigned int, Real > libMesh::LaspackLinearSolver< T >::solve ( SparseMatrix< T > &  matrix,
NumericVector< T > &  solution,
NumericVector< T > &  rhs,
const double  tol,
const unsigned int  m_its 
) [inline, virtual]

Call the Laspack solver

Implements libMesh::LinearSolver< T >.

Definition at line 108 of file laspack_linear_solver.C.

References libMesh::LaspackLinearSolver< T >::_precond_type, libMesh::LaspackMatrix< T >::_QMat, libMesh::LinearSolver< T >::_solver_type, libMeshEnums::BICG, libMeshEnums::BICGSTAB, libMeshEnums::CG, libMeshEnums::CGN, libMeshEnums::CGS, libMesh::LaspackMatrix< T >::close(), libMesh::err, libMeshEnums::GMRES, libMesh::LaspackLinearSolver< T >::init(), libMeshEnums::JACOBI, libMeshEnums::QMR, libMesh::LaspackLinearSolver< T >::set_laspack_preconditioner_type(), and libMeshEnums::SSOR.

Referenced by libMesh::LaspackLinearSolver< T >::adjoint_solve().

00113 {
00114   START_LOG("solve()", "LaspackLinearSolver");
00115   this->init ();
00116 
00117   // Make sure the data passed in are really in Laspack types
00118   LaspackMatrix<T>* matrix   = libmesh_cast_ptr<LaspackMatrix<T>*>(&matrix_in);
00119   LaspackVector<T>* solution = libmesh_cast_ptr<LaspackVector<T>*>(&solution_in);
00120   LaspackVector<T>* rhs      = libmesh_cast_ptr<LaspackVector<T>*>(&rhs_in);
00121 
00122   // Zero-out the solution to prevent the solver from exiting in 0
00123   // iterations (?)
00124   //TODO:[BSK] Why does Laspack do this?  Comment out this and try ex13...
00125   solution->zero();
00126 
00127   // Close the matrix and vectors in case this wasn't already done.
00128   matrix->close ();
00129   solution->close ();
00130   rhs->close ();
00131 
00132   // Set the preconditioner type
00133   this->set_laspack_preconditioner_type ();
00134 
00135   // Set the solver tolerance
00136   SetRTCAccuracy (tol);
00137 
00138   // Solve the linear system
00139   switch (this->_solver_type)
00140     {
00141       // Conjugate-Gradient
00142     case CG:
00143       {
00144         CGIter (&matrix->_QMat,
00145                 &solution->_vec,
00146                 &rhs->_vec,
00147                 m_its,
00148                 _precond_type,
00149                 1.);
00150         break;
00151       }
00152 
00153       // Conjugate-Gradient Normalized
00154     case CGN:
00155       {
00156         CGNIter (&matrix->_QMat,
00157                  &solution->_vec,
00158                  &rhs->_vec,
00159                  m_its,
00160                  _precond_type,
00161                  1.);
00162         break;
00163       }
00164 
00165       // Conjugate-Gradient Squared
00166     case CGS:
00167       {
00168         CGSIter (&matrix->_QMat,
00169                  &solution->_vec,
00170                  &rhs->_vec,
00171                  m_its,
00172                  _precond_type,
00173                  1.);
00174         break;
00175       }
00176 
00177       // Bi-Conjugate Gradient
00178     case BICG:
00179       {
00180         BiCGIter (&matrix->_QMat,
00181                   &solution->_vec,
00182                   &rhs->_vec,
00183                   m_its,
00184                   _precond_type,
00185                   1.);
00186         break;
00187       }
00188 
00189       // Bi-Conjugate Gradient Stabilized
00190     case BICGSTAB:
00191       {
00192         BiCGSTABIter (&matrix->_QMat,
00193                       &solution->_vec,
00194                       &rhs->_vec,
00195                       m_its,
00196                       _precond_type,
00197                       1.);
00198         break;
00199       }
00200 
00201       // Quasi-Minimum Residual
00202     case QMR:
00203       {
00204         QMRIter (&matrix->_QMat,
00205                  &solution->_vec,
00206                  &rhs->_vec,
00207                  m_its,
00208                  _precond_type,
00209                  1.);
00210         break;
00211       }
00212 
00213       // Symmetric over-relaxation
00214     case SSOR:
00215       {
00216         SSORIter (&matrix->_QMat,
00217                   &solution->_vec,
00218                   &rhs->_vec,
00219                   m_its,
00220                   _precond_type,
00221                   1.);
00222         break;
00223       }
00224 
00225       // Jacobi Relaxation
00226     case JACOBI:
00227       {
00228         JacobiIter (&matrix->_QMat,
00229                     &solution->_vec,
00230                     &rhs->_vec,
00231                     m_its,
00232                     _precond_type,
00233                     1.);
00234         break;
00235       }
00236 
00237       // Generalized Minimum Residual
00238     case GMRES:
00239       {
00240         SetGMRESRestart (30);
00241         GMRESIter (&matrix->_QMat,
00242                    &solution->_vec,
00243                    &rhs->_vec,
00244                    m_its,
00245                    _precond_type,
00246                    1.);
00247         break;
00248       }
00249 
00250       // Unknown solver, use GMRES
00251     default:
00252       {
00253         libMesh::err << "ERROR:  Unsupported LASPACK Solver: "
00254                       << this->_solver_type      << std::endl
00255                       << "Continuing with GMRES" << std::endl;
00256 
00257         this->_solver_type = GMRES;
00258 
00259         return this->solve (*matrix,
00260                             *solution,
00261                             *rhs,
00262                             tol,
00263                             m_its);
00264       }
00265     }
00266 
00267   // Check for an error
00268   if (LASResult() != LASOK)
00269     {
00270       libMesh::err << "ERROR:  LASPACK Error: " << std::endl;
00271       WriteLASErrDescr(stdout);
00272       libmesh_error();
00273     }
00274 
00275   STOP_LOG("solve()", "LaspackLinearSolver");
00276   // Get the convergence step # and residual
00277   return std::make_pair(GetLastNoIter(), GetLastAccuracy());
00278 }

template<typename T>
SolverType libMesh::LinearSolver< T >::solver_type (  )  const [inline, inherited]

Returns the type of solver to use.

Definition at line 96 of file linear_solver.h.

00096 { return _solver_type; }


Member Data Documentation

bool libMesh::ReferenceCounter::_enable_print_counter = true [static, protected, inherited]

Flag to control whether reference count information is printed when print_info is called.

Definition at line 137 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::disable_print_counter_info(), libMesh::ReferenceCounter::enable_print_counter_info(), and libMesh::ReferenceCounter::print_info().

Mutual exclusion object to enable thread-safe reference counting.

Definition at line 131 of file reference_counter.h.

Threads::atomic< unsigned int > libMesh::ReferenceCounter::_n_objects [static, protected, inherited]

The number of objects. Print the reference count information when the number returns to 0.

Definition at line 126 of file reference_counter.h.

Referenced by libMesh::ReferenceCounter::n_objects(), libMesh::ReferenceCounter::ReferenceCounter(), and libMesh::ReferenceCounter::~ReferenceCounter().

template<typename T>
bool libMesh::LinearSolver< T >::same_preconditioner [protected, inherited]

Boolean flag to indicate whether we want to use an identical preconditioner to the previous solve. This can save substantial work in the cases where the system matrix is the same for successive solves.

Definition at line 255 of file linear_solver.h.

Referenced by libMesh::PetscLinearSolver< T >::adjoint_solve(), libMesh::LinearSolver< T >::get_same_preconditioner(), libMesh::LinearSolver< T >::reuse_preconditioner(), and libMesh::PetscLinearSolver< T >::solve().


The documentation for this class was generated from the following files:

Site Created By: libMesh Developers
Last modified: February 05 2013 19:55:27 UTC

Hosted By:
SourceForge.net Logo