trilinos_aztec_linear_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 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.
34 #include <Epetra_LinearProblem.h>
35 #include <AztecOO.h>
37 
38 // C++ includes
39 #include <vector>
40 
41 namespace libMesh
42 {
43 
52 template <typename T>
54 {
55 public:
60 
65 
69  virtual void clear () override;
70 
74  virtual void init (const char * name=nullptr) override;
75 
80  virtual std::pair<unsigned int, Real>
81  solve (SparseMatrix<T> & matrix_in,
82  NumericVector<T> & solution_in,
83  NumericVector<T> & rhs_in,
84  const double tol,
85  const unsigned int m_its) override
86  {
87  return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
88  }
89 
97  virtual std::pair<unsigned int, Real>
98  solve (SparseMatrix<T> & matrix,
99  SparseMatrix<T> & preconditioner,
100  NumericVector<T> & solution,
101  NumericVector<T> & rhs,
102  const double tol,
103  const unsigned int m_its) override;
104 
108  virtual std::pair<unsigned int, Real>
109  solve (const ShellMatrix<T> & shell_matrix,
110  NumericVector<T> & solution_in,
111  NumericVector<T> & rhs_in,
112  const double tol,
113  const unsigned int m_its) override;
114 
120  virtual std::pair<unsigned int, Real>
121  solve (const ShellMatrix<T> & shell_matrix,
122  const SparseMatrix<T> & precond_matrix,
123  NumericVector<T> & solution_in,
124  NumericVector<T> & rhs_in,
125  const double tol,
126  const unsigned int m_its) override;
127 
132  void get_residual_history(std::vector<double> & hist);
133 
142 
147  virtual void print_converged_reason() const override;
148 
152  virtual LinearConvergenceReason get_converged_reason() const override;
153 
154 private:
155 
160  void set_solver_type ();
161 
165  Epetra_LinearProblem * _linear_problem;
166 
170  AztecOO * _linear_solver;
171 };
172 
173 
174 /*----------------------- functions ----------------------------------*/
175 template <typename T>
176 inline
178 {
179  this->clear ();
180 }
181 
182 } // namespace libMesh
183 
184 
185 
186 #endif // #ifdef LIBMESH_TRILINOS_HAVE_AZTECOO
187 #endif // LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H
std::string name(const ElemQuality q)
Definition: elem_quality.C:42
void get_residual_history(std::vector< double > &hist)
Provides a uniform interface to vector storage schemes for different linear algebra libraries...
Definition: diff_context.h:40
const Parallel::Communicator & comm() const
AztecLinearSolver(const libMesh::Parallel::Communicator &comm)
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) override
virtual void init(const char *name=nullptr) override
virtual void print_converged_reason() const override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual LinearConvergenceReason get_converged_reason() const override