libMesh
trilinos_aztec_linear_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2017 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 
20 #ifndef LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H
21 #define LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H
22 
23 
24 
25 // Local includes
26 #include "libmesh/libmesh_common.h"
27 #include "libmesh/linear_solver.h"
28 
29 #ifdef LIBMESH_TRILINOS_HAVE_AZTECOO
30 
31 // Trilinos include files. AztecOO requires Epetra, so there's no
32 // need to check for both.
33 #include <Epetra_LinearProblem.h>
34 #include <AztecOO.h>
35 
36 // C++ includes
37 #include <vector>
38 
39 namespace libMesh
40 {
41 
50 template <typename T>
52 {
53 public:
58  LIBMESH_CAN_DEFAULT_TO_COMMWORLD);
59 
64 
68  virtual void clear () libmesh_override;
69 
73  virtual void init (const char * name=libmesh_nullptr) libmesh_override;
74 
79  virtual std::pair<unsigned int, Real>
80  solve (SparseMatrix<T> & matrix_in,
81  NumericVector<T> & solution_in,
82  NumericVector<T> & rhs_in,
83  const double tol,
84  const unsigned int m_its) libmesh_override
85  {
86  return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
87  }
88 
96  virtual std::pair<unsigned int, Real>
97  solve (SparseMatrix<T> & matrix,
98  SparseMatrix<T> & preconditioner,
99  NumericVector<T> & solution,
100  NumericVector<T> & rhs,
101  const double tol,
102  const unsigned int m_its) libmesh_override;
103 
107  virtual std::pair<unsigned int, Real>
108  solve (const ShellMatrix<T> & shell_matrix,
109  NumericVector<T> & solution_in,
110  NumericVector<T> & rhs_in,
111  const double tol,
112  const unsigned int m_its) libmesh_override;
113 
119  virtual std::pair<unsigned int, Real>
120  solve (const ShellMatrix<T> & shell_matrix,
121  const SparseMatrix<T> & precond_matrix,
122  NumericVector<T> & solution_in,
123  NumericVector<T> & rhs_in,
124  const double tol,
125  const unsigned int m_its) libmesh_override;
126 
131  void get_residual_history(std::vector<double> & hist);
132 
141 
146  virtual void print_converged_reason() const libmesh_override;
147 
151  virtual LinearConvergenceReason get_converged_reason() const libmesh_override;
152 
153 private:
154 
159  void set_solver_type ();
160 
164  Epetra_LinearProblem * _linear_problem;
165 
169  AztecOO * _linear_solver;
170 };
171 
172 
173 /*----------------------- functions ----------------------------------*/
174 template <typename T>
175 inline
177  LinearSolver<T>(comm)
178 {
179  if (this->n_processors() == 1)
181  else
183 }
184 
185 
186 
187 template <typename T>
188 inline
190 {
191  this->clear ();
192 }
193 
194 } // namespace libMesh
195 
196 
197 
198 #endif // #ifdef LIBMESH_TRILINOS_HAVE_AZTECOO
199 #endif // LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H
std::string name(const ElemQuality q)
This function returns a string containing some name for q.
Definition: elem_quality.C:39
Encapsulates the MPI_Comm object.
Definition: parallel.h:657
virtual void init(const char *name=libmesh_nullptr) libmesh_override
Initialize data structures if not done so already.
virtual void print_converged_reason() const libmesh_override
Prints a useful message about why the latest linear solve con(di)verged.
Epetra_LinearProblem * _linear_problem
The Epetra linear problem object.
processor_id_type n_processors() const
virtual void clear() libmesh_override
Release all memory and clear data structures.
void get_residual_history(std::vector< double > &hist)
Fills the input vector with the sequence of residual norms from the latest iterative solve...
const class libmesh_nullptr_t libmesh_nullptr
Numeric vector.
Definition: dof_map.h:66
The libMesh namespace provides an interface to certain functionality in the library.
Generic sparse matrix.
Definition: dof_map.h:65
This base class can be inherited from to provide interfaces to linear solvers from different packages...
Definition: linear_solver.h:58
virtual LinearConvergenceReason get_converged_reason() const libmesh_override
This class provides an interface to AztecOO iterative solvers that is compatible with the libMesh Lin...
void set_solver_type()
Tells AztecOO to use the user-specified solver stored in _solver_type.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual std::pair< unsigned int, Real > solve(SparseMatrix< T > &matrix_in, NumericVector< T > &solution_in, NumericVector< T > &rhs_in, const double tol, const unsigned int m_its) libmesh_override
Call the Aztec solver.
const Parallel::Communicator & comm() const
PreconditionerType _preconditioner_type
Enum stating with type of preconditioner to use.
LinearConvergenceReason
Linear solver convergence flags (taken from the PETSc flags)
AztecLinearSolver(const libMesh::Parallel::Communicator &comm LIBMESH_CAN_DEFAULT_TO_COMMWORLD)
Constructor.
Generic shell matrix, i.e.
AztecOO * _linear_solver
The AztecOO solver object.