/home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/timestepping/linear_solver/GMRES.hh Source File#

DiFfRG: /home/runner/work/DiFfRG_current/DiFfRG_current/DiFfRG/include/DiFfRG/timestepping/linear_solver/GMRES.hh Source File
DiFfRG
GMRES.hh
Go to the documentation of this file.
1#pragma once
2
3// external libraries
4#include <deal.II/lac/precondition.h>
5#include <deal.II/lac/solver_gmres.h>
6
7// DiFfRG
9
10namespace DiFfRG
11{
12 template <typename SparseMatrixType, typename VectorType,
13 typename PreconditionerType = dealii::PreconditionJacobi<SparseMatrixType>>
14 class GMRES : public AbstractLinearSolver<SparseMatrixType, VectorType>
15 {
16 public:
17 GMRES() : matrix(nullptr) {}
18
19 void init(const SparseMatrixType &matrix)
20 {
21 this->matrix = &matrix;
22 preconditioner.initialize(matrix, 1.0);
23 }
24
25 bool invert() { return false; }
26
27 int solve(const VectorType &src, VectorType &dst, const double tol)
28 {
29 if (!matrix) throw std::runtime_error("GMRES::solve: matrix not initialized");
30 dealii::SolverControl solver_control(std::max<std::size_t>(1000, src.size() / 10), tol);
31 dealii::SolverGMRES<VectorType> solver(solver_control);
32
33 try {
34 solver.solve(*matrix, dst, src, preconditioner);
35 } catch (std::exception &e) {
36 std::cerr << "GMRES linear solver failed: " << e.what() << std::endl;
37 throw;
38 }
39
40 int steps = solver_control.last_step();
41 return steps;
42 }
43
44 private:
45 const SparseMatrixType *matrix;
46 PreconditionerType preconditioner;
47 };
48} // namespace DiFfRG
Definition abstract_linear_solver.hh:9
Definition GMRES.hh:15
GMRES()
Definition GMRES.hh:17
bool invert()
Definition GMRES.hh:25
void init(const SparseMatrixType &matrix)
Definition GMRES.hh:19
int solve(const VectorType &src, VectorType &dst, const double tol)
Definition GMRES.hh:27
PreconditionerType preconditioner
Definition GMRES.hh:46
const SparseMatrixType * matrix
Definition GMRES.hh:45
Definition complex_math.hh:10