AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
IterativeRunner.hpp
1 #pragma once
2 
3 #include <memory>
4 
5 #include <amdis/CreatorInterface.hpp>
6 #include <amdis/Initfile.hpp>
7 #include <amdis/linearalgebra/LinearSolverInterface.hpp>
8 #include <amdis/linearalgebra/eigen/PreconConfig.hpp>
9 #include <amdis/linearalgebra/eigen/SolverConfig.hpp>
10 
11 namespace AMDiS
12 {
13  template <class M, class X, class Y, class IterativeSolver>
15  : public LinearSolverInterface<M,X,Y>
16  {
19 
20  public:
21  IterativeRunner(std::string const& prefix)
22  : solver_{}
23  {
24  SolverCfg::init(prefix, solver_);
25  PreconCfg::init(prefix + "->precon", solver_.preconditioner());
26  Parameters::get(prefix + "->reuse pattern", reusePattern_);
27  }
28 
30  void init(M const& A)
31  {
32  if (!reusePattern_ || !initialized_) {
33  solver_.analyzePattern(A);
34  initialized_ = true;
35  }
36  solver_.factorize(A);
37 
38  test_exit(solver_.info() == Eigen::Success, "Error in solver.compute(matrix)");
39  }
40 
42  void finish()
43  {
44  initialized_ = false;
45  }
46 
48  void apply(X& x, Y const& b, Dune::InverseOperatorResult& stat) override
49  {
50  Dune::Timer t;
51  x = solver_.solveWithGuess(b, x);
52 
53  stat.reduction = solver_.error();
54  stat.iterations = solver_.iterations();
55  stat.converged = (solver_.info() == Eigen::Success);
56  stat.elapsed = t.elapsed();
57  }
58 
59  private:
60  IterativeSolver solver_;
61  bool reusePattern_ = false;
62  bool initialized_ = false;
63  };
64 
65 } // end namespace AMDiS
void apply(X &x, Y const &b, Dune::InverseOperatorResult &stat) override
Implements LinearSolverInterface::apply()
Definition: IterativeRunner.hpp:48
void finish()
Implements LinearSolverInterface::finish()
Definition: IterativeRunner.hpp:42
Definition: PreconConfig.hpp:12
Definition: AdaptBase.hpp:6
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
Definition: IterativeRunner.hpp:14
Definition: LinearSolverInterface.hpp:8
void init(M const &A)
initialize the matrix and maybe compute factorization
Definition: IterativeRunner.hpp:30
Definition: SolverConfig.hpp:13