$extrastylesheet
trilinos_aztec_linear_solver.h
Go to the documentation of this file.
00001 // The libMesh Finite Element Library.
00002 // Copyright (C) 2002-2014 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
00003 
00004 // This library is free software; you can redistribute it and/or
00005 // modify it under the terms of the GNU Lesser General Public
00006 // License as published by the Free Software Foundation; either
00007 // version 2.1 of the License, or (at your option) any later version.
00008 
00009 // This library is distributed in the hope that it will be useful,
00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of
00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00012 // Lesser General Public License for more details.
00013 
00014 // You should have received a copy of the GNU Lesser General Public
00015 // License along with this library; if not, write to the Free Software
00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00017 
00018 
00019 
00020 #ifndef LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H
00021 #define LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H
00022 
00023 
00024 
00025 // Local includes
00026 #include "libmesh/libmesh_common.h"
00027 #include "libmesh/linear_solver.h"
00028 
00032 #ifdef LIBMESH_HAVE_AZTECOO
00033 #include <Epetra_LinearProblem.h>
00034 #include <AztecOO.h>
00035 
00036 // C++ includes
00037 #include <vector>
00038 
00039 namespace libMesh
00040 {
00041 
00050 template <typename T>
00051 class AztecLinearSolver : public LinearSolver<T>
00052 {
00053 public:
00057   AztecLinearSolver (const libMesh::Parallel::Communicator &comm
00058                      LIBMESH_CAN_DEFAULT_TO_COMMWORLD);
00059 
00063   ~AztecLinearSolver ();
00064 
00068   void clear ();
00069 
00073   void init (const char *name=NULL);
00074 
00079   std::pair<unsigned int, Real>
00080   solve (SparseMatrix<T>  &matrix_in,
00081          NumericVector<T> &solution_in,
00082          NumericVector<T> &rhs_in,
00083          const double tol,
00084          const unsigned int m_its)
00085   {
00086     return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
00087   }
00088 
00095   std::pair<unsigned int, Real>
00096   solve (SparseMatrix<T>  &matrix,
00097          SparseMatrix<T>  &preconditioner,
00098          NumericVector<T> &solution,
00099          NumericVector<T> &rhs,
00100          const double tol,
00101          const unsigned int m_its);
00102 
00106   std::pair<unsigned int, Real>
00107   solve (const ShellMatrix<T>& shell_matrix,
00108          NumericVector<T>& solution_in,
00109          NumericVector<T>& rhs_in,
00110          const double tol,
00111          const unsigned int m_its);
00112 
00118   virtual std::pair<unsigned int, Real>
00119   solve (const ShellMatrix<T>& shell_matrix,
00120          const SparseMatrix<T>& precond_matrix,
00121          NumericVector<T>& solution_in,
00122          NumericVector<T>& rhs_in,
00123          const double tol,
00124          const unsigned int m_its);
00125 
00130   void get_residual_history(std::vector<double>& hist);
00131 
00138   Real get_initial_residual();
00139 
00143   virtual LinearConvergenceReason get_converged_reason() const;
00144 
00145 private:
00146 
00151   void set_solver_type ();
00152 
00156   Epetra_LinearProblem * _linear_problem;
00157 
00161   AztecOO * _linear_solver;
00162 };
00163 
00164 
00165 /*----------------------- functions ----------------------------------*/
00166 template <typename T>
00167 inline
00168 AztecLinearSolver<T>::AztecLinearSolver (const libMesh::Parallel::Communicator &comm) :
00169   LinearSolver<T>(comm)
00170 {
00171   if (this->n_processors() == 1)
00172     this->_preconditioner_type = ILU_PRECOND;
00173   else
00174     this->_preconditioner_type = BLOCK_JACOBI_PRECOND;
00175 }
00176 
00177 
00178 
00179 template <typename T>
00180 inline
00181 AztecLinearSolver<T>::~AztecLinearSolver ()
00182 {
00183   this->clear ();
00184 }
00185 
00186 } // namespace libMesh
00187 
00188 
00189 
00190 #endif // #ifdef LIBMESH_HAVE_AZTECOO
00191 #endif // LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H