DOLFIN
DOLFIN C++ interface
PETScKrylovSolver.h
1// Copyright (C) 2004-2015 Johan Jansson and Garth N. Wells
2//
3// This file is part of DOLFIN.
4//
5// DOLFIN is free software: you can redistribute it and/or modify
6// it under the terms of the GNU Lesser General Public License as published by
7// the Free Software Foundation, either version 3 of the License, or
8// (at your option) any later version.
9//
10// DOLFIN is distributed in the hope that it will be useful,
11// but WITHOUT ANY WARRANTY; without even the implied warranty of
12// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
13// GNU Lesser General Public License for more details.
14//
15// You should have received a copy of the GNU Lesser General Public License
16// along with DOLFIN. If not, see <http://www.gnu.org/licenses/>.
17//
18// Modified by Anders Logg 2005-2012
19// Modified by Johan Hoffman 2005
20// Modified by Andy R. Terrel 2005
21// Modified by Garth N. Wells 2005-2010
22
23#ifndef __DOLFIN_PETSC_KRYLOV_SOLVER_H
24#define __DOLFIN_PETSC_KRYLOV_SOLVER_H
25
26#ifdef HAS_PETSC
27
28#include <map>
29#include <memory>
30#include <string>
31#include <petscksp.h>
32#include <dolfin/common/types.h>
33#include "GenericLinearSolver.h"
34#include "PETScObject.h"
35
36namespace dolfin
37{
38
40 class GenericMatrix;
41 class GenericVector;
42 class PETScBaseMatrix;
43 class PETScMatrix;
44 class PETScVector;
45 class PETScPreconditioner;
46 class PETScSNESSolver;
47 class VectorSpaceBasis;
48
49 class PETScDMCollection;
50
53
55 {
56 public:
57
62 enum class norm_type {none, default_norm, preconditioned, unpreconditioned, natural};
63
66 PETScKrylovSolver(MPI_Comm comm,
67 std::string method="default",
68 std::string preconditioner="default");
69
72 PETScKrylovSolver(std::string method="default",
73 std::string preconditioner="default");
74
77 PETScKrylovSolver(MPI_Comm comm,
78 std::string method,
79 std::shared_ptr<PETScPreconditioner> preconditioner);
80
83 PETScKrylovSolver(std::string method,
84 std::shared_ptr<PETScPreconditioner> preconditioner);
85
87 explicit PETScKrylovSolver(KSP ksp);
88
90 virtual ~PETScKrylovSolver();
91
93 void set_operator(std::shared_ptr<const GenericLinearOperator> A);
94
97 void set_operator(const PETScBaseMatrix& A);
98
100 void set_operators(std::shared_ptr<const GenericLinearOperator> A,
101 std::shared_ptr<const GenericLinearOperator> P);
102
106 void set_operators(const PETScBaseMatrix& A, const PETScBaseMatrix& P);
107
109 std::size_t solve(GenericVector& x, const GenericVector& b);
110
113 std::size_t solve(GenericVector& x, const GenericVector& b,
114 bool transpose);
115
118 std::size_t solve(PETScVector& x, const PETScVector& b,
119 bool transpose=false);
120
122 std::size_t solve(const GenericLinearOperator& A, GenericVector& x,
123 const GenericVector& b);
124
128 void set_nonzero_guess(bool nonzero_guess);
129
133 void set_reuse_preconditioner(bool reuse_pc);
134
137 void set_tolerances(double relative, double absolute, double diverged,
138 int max_iter);
139
142 void set_norm_type(norm_type type);
143
145 norm_type get_norm_type() const;
146
148 void monitor(bool monitor_convergence);
149
152 void set_options_prefix(std::string options_prefix);
153
156 std::string get_options_prefix() const;
157
159 void set_from_options() const;
160
162 std::string str(bool verbose) const;
163
165 MPI_Comm mpi_comm() const;
166
168 KSP ksp() const;
169
171 static std::map<std::string, std::string> methods();
172
174 static std::map<std::string, std::string> preconditioners();
175
178
180 std::string parameter_type() const
181 { return "krylov_solver"; }
182
184 void set_dm(DM dm);
185
187 void set_dm_active(bool val);
188
189 friend class PETScSNESSolver;
190 friend class PETScTAOSolver;
191
192 private:
193
194 // Return norm_type enum for norm string
196
197 // Solve linear system Ax = b and return number of iterations
198 std::size_t _solve(const PETScBaseMatrix& A, PETScVector& x,
199 const PETScVector& b);
200
201 // Report the number of iterations
202 void write_report(int num_iterations, KSPConvergedReason reason);
203
204 void check_dimensions(const PETScBaseMatrix& A, const GenericVector& x,
205 const GenericVector& b) const;
206
207 // Available solvers
208 static const std::map<std::string, const KSPType> _methods;
209
210 // Available solvers descriptions
211 static const std::map<std::string, std::string> _methods_descr;
212
213 // PETSc solver pointer
214 KSP _ksp;
215
216 // Preconditioner
217 std::shared_ptr<PETScPreconditioner> _preconditioner;
218
219 bool preconditioner_set;
220
221 };
222
223}
224
225#endif
226
227#endif
Definition: GenericLinearOperator.h:43
This class provides a general solver for linear systems Ax = b.
Definition: GenericLinearSolver.h:38
This class defines a common interface for vectors.
Definition: GenericVector.h:48
Definition: PETScBaseMatrix.h:50
Definition: PETScKrylovSolver.h:55
void set_dm(DM dm)
Set the DM.
Definition: PETScKrylovSolver.cpp:446
static Parameters default_parameters()
Default parameter values.
Definition: PETScKrylovSolver.cpp:80
static std::map< std::string, std::string > methods()
Return a list of available solver methods.
Definition: PETScKrylovSolver.cpp:70
void set_operators(std::shared_ptr< const GenericLinearOperator > A, std::shared_ptr< const GenericLinearOperator > P)
Set operator (matrix) and preconditioner matrix.
Definition: PETScKrylovSolver.cpp:200
std::size_t solve(GenericVector &x, const GenericVector &b)
Solve linear system Ax = b and return number of iterations.
Definition: PETScKrylovSolver.cpp:220
void set_operator(std::shared_ptr< const GenericLinearOperator > A)
Set operator (matrix)
Definition: PETScKrylovSolver.cpp:190
void set_options_prefix(std::string options_prefix)
Definition: PETScKrylovSolver.cpp:512
void set_norm_type(norm_type type)
Definition: PETScKrylovSolver.cpp:416
virtual ~PETScKrylovSolver()
Destructor.
Definition: PETScKrylovSolver.cpp:181
KSP ksp() const
Return PETSc KSP pointer.
Definition: PETScKrylovSolver.cpp:561
norm_type get_norm_type() const
Get norm type used in convergence testing.
Definition: PETScKrylovSolver.cpp:461
void monitor(bool monitor_convergence)
Monitor residual at each iteration.
Definition: PETScKrylovSolver.cpp:489
std::string get_options_prefix() const
Definition: PETScKrylovSolver.cpp:520
void set_nonzero_guess(bool nonzero_guess)
Definition: PETScKrylovSolver.cpp:391
void set_reuse_preconditioner(bool reuse_pc)
Definition: PETScKrylovSolver.cpp:399
PETScKrylovSolver(MPI_Comm comm, std::string method="default", std::string preconditioner="default")
Definition: PETScKrylovSolver.cpp:92
static std::map< std::string, std::string > preconditioners()
Return a list of available named preconditioners.
Definition: PETScKrylovSolver.cpp:75
void set_tolerances(double relative, double absolute, double diverged, int max_iter)
Definition: PETScKrylovSolver.cpp:407
void set_from_options() const
Set options from PETSc options database.
Definition: PETScKrylovSolver.cpp:529
std::string str(bool verbose) const
Return informal string representation (pretty-print)
Definition: PETScKrylovSolver.cpp:536
void set_dm_active(bool val)
Activate/deactivate DM.
Definition: PETScKrylovSolver.cpp:452
norm_type
Definition: PETScKrylovSolver.h:62
std::string parameter_type() const
Return parameter type: "krylov_solver" or "lu_solver".
Definition: PETScKrylovSolver.h:180
MPI_Comm mpi_comm() const
Return MPI communicator.
Definition: PETScKrylovSolver.cpp:553
Definition: PETScObject.h:34
Definition: PETScSNESSolver.h:46
Definition: PETScTAOSolver.h:52
Definition: PETScVector.h:61
Definition: Parameters.h:95
Definition: adapt.h:30
double norm(const GenericVector &x, std::string norm_type="l2")
Definition: solve.cpp:173