trilinos_aztec_linear_solver.h
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2014 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 
32 #ifdef LIBMESH_HAVE_AZTECOO
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  void clear ();
69 
73  void init ();
74 
79  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)
85  {
86  return this->solve(matrix_in, matrix_in, solution_in, rhs_in, tol, m_its);
87  }
88 
95  std::pair<unsigned int, Real>
96  solve (SparseMatrix<T> &matrix,
97  SparseMatrix<T> &preconditioner,
98  NumericVector<T> &solution,
99  NumericVector<T> &rhs,
100  const double tol,
101  const unsigned int m_its);
102 
106  std::pair<unsigned int, Real>
107  solve (const ShellMatrix<T>& shell_matrix,
108  NumericVector<T>& solution_in,
109  NumericVector<T>& rhs_in,
110  const double tol,
111  const unsigned int m_its);
112 
118  virtual std::pair<unsigned int, Real>
119  solve (const ShellMatrix<T>& shell_matrix,
120  const SparseMatrix<T>& precond_matrix,
121  NumericVector<T>& solution_in,
122  NumericVector<T>& rhs_in,
123  const double tol,
124  const unsigned int m_its);
125 
130  void get_residual_history(std::vector<double>& hist);
131 
139 
144  virtual void print_converged_reason();
145 
146 private:
147 
152  void set_solver_type ();
153 
157  Epetra_LinearProblem * _linear_problem;
158 
162  AztecOO * _linear_solver;
163 };
164 
165 
166 /*----------------------- functions ----------------------------------*/
167 template <typename T>
168 inline
170  LinearSolver<T>(comm)
171 {
172  if (this->n_processors() == 1)
174  else
176 }
177 
178 
179 
180 template <typename T>
181 inline
183 {
184  this->clear ();
185 }
186 
187 } // namespace libMesh
188 
189 
190 
191 #endif // #ifdef LIBMESH_HAVE_AZTECOO
192 #endif // LIBMESH_TRILINOS_AZTEC_LINEAR_SOLVER_H

Site Created By: libMesh Developers
Last modified: February 07 2014 16:57:07 UTC

Hosted By:
SourceForge.net Logo