AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
Solvers.hpp
1 #pragma once
2 
3 #include <string>
4 
5 // MTL4 includes
6 #include <boost/numeric/itl/itl.hpp>
7 
8 // more solvers defined in AMDiS
9 #include <amdis/linearalgebra/mtl/itl/minres.hpp>
10 #include <amdis/linearalgebra/mtl/itl/gcr.hpp>
11 #include <amdis/linearalgebra/mtl/itl/fgmres.hpp>
12 #include <amdis/linearalgebra/mtl/itl/fgmres_householder.hpp>
13 #include <amdis/linearalgebra/mtl/itl/gmres2.hpp>
14 #include <amdis/linearalgebra/mtl/itl/gmres_householder.hpp>
15 #include <amdis/linearalgebra/mtl/itl/preonly.hpp>
16 
17 #include <amdis/CreatorMap.hpp>
18 #include <amdis/Initfile.hpp>
19 #include <amdis/linearalgebra/mtl/KrylovRunner.hpp>
20 #include <amdis/linearalgebra/mtl/Traits.hpp>
21 #include <amdis/linearalgebra/mtl/UmfpackRunner.hpp>
22 
23 namespace AMDiS
24 {
35  {
36  public:
37  explicit cg_solver_type(std::string const& /*prefix*/) {}
38  template <class LinOp, class X, class B, class P, class I>
39  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
40  {
41  return itl::cg(A, x, b, precon, iter);
42  }
43  };
44 
45 
55  {
56  public:
57  explicit cgs_solver_type(std::string const& /*prefix*/) {}
58  template <class LinOp, class X, class B, class P, class I>
59  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
60  {
61  return itl::cgs(A, x, b, precon, iter);
62  }
63  };
64 
65 
75  {
76  public:
77  explicit bicg_solver_type(std::string const& /*prefix*/) {}
78  template <class LinOp, class X, class B, class P, class I>
79  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
80  {
81  return itl::bicg(A, x, b, precon, iter);
82  }
83  };
84 
85 
95  {
96  public:
97  explicit bicgstab_type(std::string const& /*prefix*/) {}
98  template <class LinOp, class X, class B, class P, class I>
99  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
100  {
101  return itl::bicgstab(A, x, b, precon, iter);
102  }
103  };
104 
105 
115  {
116  public:
117  explicit bicgstab2_type(std::string const& /*prefix*/) {}
118  template <class LinOp, class X, class B, class P, class I>
119  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
120  {
121  return itl::bicgstab_2(A, x, b, precon, iter);
122  }
123  };
124 
125 
134  {
135  public:
136  explicit qmr_solver_type(std::string const& /*prefix*/) {}
137  template <class LinOp, class X, class B, class P, class I>
138  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
139  {
140  itl::pc::identity<LinOp> id(A);
141  return itl::qmr(A, x, b, precon, id, iter);
142  }
143  };
144 
145 
155  {
156  public:
157  explicit tfqmr_solver_type(std::string const& /*prefix*/) {}
158  template <class LinOp, class X, class B, class P, class I>
159  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
160  {
161  itl::pc::identity<LinOp> id(A);
162  return itl::tfqmr(A, x, b, precon, id, iter);
163  }
164  };
165 
166 
176  {
177  int ell;
178  public:
179  explicit bicgstab_ell_type(std::string const& prefix) : ell(3)
180  {
181  Parameters::get(prefix + "->ell", ell);
182  }
183  template <class LinOp, class X, class B, class P, class I>
184  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
185  {
186  itl::pc::identity<LinOp> id(A);
187  return itl::bicgstab_ell(A, x, b, precon, id, iter, ell);
188  }
189  };
190 
191 
202  {
203 
204  public:
205  explicit gmres_type(std::string const& prefix) : restart(30), ortho(GRAM_SCHMIDT)
206  {
207  Parameters::get(prefix + "->restart", restart);
208  Parameters::get(prefix + "->orthogonalization", ortho);
209  }
210  template <class LinOp, class X, class B, class P, class I>
211  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
212  {
213  itl::pc::identity<LinOp> id(A);
214  switch (ortho)
215  {
216  default:
217  case GRAM_SCHMIDT:
218  return itl::gmres2(A, x, b, id, precon, iter, restart);
219  break;
220  case HOUSEHOLDER:
221  return itl::gmres_householder(A, x, b, precon, iter, restart);
222  break;
223  }
224  }
225 
226  private:
227  static constexpr int GRAM_SCHMIDT = 1;
228  static constexpr int HOUSEHOLDER = 2;
229 
230  int restart;
231  int ortho;
232  };
233 
234 
250  {
251  int s;
252  public:
253  explicit idr_s_type(std::string const& prefix) : s(30)
254  {
255  Parameters::get(prefix + "->s", s);
256  }
257  template <class LinOp, class X, class B, class P, class I>
258  int operator()(LinOp const& A, X& x, B const& b, P const&, I& iter)
259  {
260  itl::pc::identity<LinOp> id(A);
261  return itl::idr_s(A, x, b, id, id, iter, s);
262  }
263  };
264 
265 
275  {
276  public:
277  explicit minres_solver_type(std::string const& /*prefix*/) {}
278  template <class LinOp, class X, class B, class P, class I>
279  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
280  {
281  return itl::minres(A, x, b, precon, iter);
282  }
283  };
284 
285 
295  class gcr_type
296  {
297  int restart;
298 
299  public:
300  explicit gcr_type(std::string const& prefix) : restart(30)
301  {
302  Parameters::get(prefix + "->restart", restart);
303  }
304  template <class LinOp, class X, class B, class P, class I>
305  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
306  {
307  itl::pc::identity<LinOp> id(A);
308  return itl::gcr(A, x, b, id, precon, iter, restart);
309  }
310  };
311 
312 
324  {
325  public:
326  explicit fgmres_type(std::string const& prefix)
327  : restart(30), orthogonalization(GRAM_SCHMIDT)
328  {
329  Parameters::get(prefix + "->restart", restart);
330  Parameters::get(prefix + "->orthogonalization", orthogonalization);
331  }
332  template <class LinOp, class X, class B, class P, class I>
333  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
334  {
335  itl::pc::identity<LinOp> id(A);
336  switch (orthogonalization)
337  {
338  default:
339  case GRAM_SCHMIDT:
340  return itl::fgmres(A, x, b, id, precon, iter, restart);
341  break;
342  case HOUSEHOLDER:
343  return itl::fgmres_householder(A, x, b, precon, iter, restart);
344  break;
345  }
346  }
347 
348  private:
349  static constexpr int GRAM_SCHMIDT = 1;
350  static constexpr int HOUSEHOLDER = 2;
351 
352  int restart;
353  int orthogonalization;
354  };
355 
356 
365  {
366  public:
367  explicit preonly_type(std::string const& /*prefix*/) {}
368  template <class LinOp, class X, class B, class P, class I>
369  int operator()(LinOp const& A, X& x, B const& b, P const& precon, I& iter)
370  {
371  return itl::preonly(A, x, b, precon, iter);
372  }
373  };
374 
375 
376  // ===========================================================================
377 
379 
397  template <class M, class X, class Y>
398  class DefaultCreators< LinearSolverInterface<M,X,Y> >
399  {
400  using SolverBase = LinearSolverInterface<M,X,Y>;
401 
402  template <class ITLSolver>
403  using SolverCreator = typename KrylovRunner<M,X,Y,ITLSolver>::Creator;
404 
405 #ifdef HAVE_UMFPACK
406  using UmfpackSolverCreator = typename UmfpackRunner<M,X,Y>::Creator;
407 #endif
408 
409  using Map = CreatorMap<SolverBase>;
410 
411  public:
412  static void init()
413  {
414  auto cg = new SolverCreator<cg_solver_type>;
415  Map::addCreator("cg", cg);
416 
417  auto cgs = new SolverCreator<cgs_solver_type>;
418  Map::addCreator("cgs", cgs);
419 
420  auto bicg = new SolverCreator<bicg_solver_type>;
421  Map::addCreator("bicg", bicg);
422 
423  auto bcgs = new SolverCreator<bicgstab_type>;
424  Map::addCreator("bicgstab", bcgs);
425  Map::addCreator("bcgs", bcgs);
426 
427  auto bcgs2 = new SolverCreator<bicgstab2_type>;
428  Map::addCreator("bicgstab2", bcgs2);
429  Map::addCreator("bcgs2", bcgs2);
430 
431  auto qmr = new SolverCreator<qmr_solver_type>;
432  Map::addCreator("qmr", qmr);
433 
434  auto tfqmr = new SolverCreator<tfqmr_solver_type>;
435  Map::addCreator("tfqmr", tfqmr);
436 
437  auto bcgsl = new SolverCreator<bicgstab_ell_type>;
438  Map::addCreator("bicgstab_ell", bcgsl);
439  Map::addCreator("bcgsl", bcgsl);
440 
441  auto gmres = new SolverCreator<gmres_type>;
442  Map::addCreator("gmres", gmres);
443 
444  auto fgmres = new SolverCreator<fgmres_type>;
445  Map::addCreator("fgmres", fgmres);
446 
447  auto minres = new SolverCreator<minres_solver_type>;
448  Map::addCreator("minres", minres);
449 
450  auto idrs = new SolverCreator<idr_s_type>;
451  Map::addCreator("idrs", idrs);
452 
453  auto gcr = new SolverCreator<gcr_type>;
454  Map::addCreator("gcr", gcr);
455 
456  auto preonly = new SolverCreator<preonly_type>;
457  Map::addCreator("preonly", preonly);
458 
459  // default iterative solver
460  Map::addCreator("default", gmres);
461 
462  init_direct(std::is_same<typename Dune::FieldTraits<typename M::value_type>::real_type, double>{});
463  }
464 
465  static void init_direct(std::false_type) {}
466  static void init_direct(std::true_type)
467  {
468 #ifdef HAVE_UMFPACK
469  auto umfpack = new UmfpackSolverCreator;
470  Map::addCreator("umfpack", umfpack);
471 
472  // default direct solvers
473  Map::addCreator("direct", umfpack);
474 #endif
475  }
476  };
477 
478 } // end namespace AMDiS
ITL_Solver <fgmres_type> implementation of flexible GMRes method.
Definition: Solvers.hpp:323
ITL_Solver <qmr_solver_type> implementation of Quasi-Minimal Residual method.
Definition: Solvers.hpp:133
ITL_Solver <bicgstab_type> implementation of stabilized bi-conjugate gradient method.
Definition: Solvers.hpp:94
ITL_Solver <tfqmr_solver_type> implementation of Transposed-Free Quasi-Minimal Residual method...
Definition: Solvers.hpp:154
ITL_Solver <gmres_type> implementation of generalized minimal residual method.
Definition: Solvers.hpp:201
ITL_Solver <bicgstab2_type> implementation of BiCGStab(l) method with l=2.
Definition: Solvers.hpp:114
Definition: UmfpackRunner.hpp:36
A CreatorMap is used to construct objects, which types depends on key words determined at run time...
Definition: CreatorMap.hpp:29
Definition: AdaptBase.hpp:6
static std::optional< T > get(std::string const &key)
Get parameter-values from parameter-tree.
Definition: Initfile.hpp:25
Definition: LinearSolverInterface.hpp:8
ITL_Solver <gcr_type> implementation of generalized conjugate residual method.
Definition: Solvers.hpp:295
Definition: CreatorMap.hpp:16
ITL_Solver <cg_solver_type> implementation of conjugate gradient method.
Definition: Solvers.hpp:34
ITL_Solver <minres_solver_type> implementation of minimal residual method.
Definition: Solvers.hpp:274
ITL_Solver <bicgstab_ell_type> implementation of stabilized BiCG(ell) method.
Definition: Solvers.hpp:175
ITL_Solver <idr_s_type> implementation of Induced Dimension Reduction method.
Definition: Solvers.hpp:249
ITL_Solver <bicg_solver_type> implementation of bi-conjugate gradient method.
Definition: Solvers.hpp:74
ITL_Solver <preonly_type> implementation of preconditioner as.
Definition: Solvers.hpp:364
ITL_Solver <cgs_solver_type> implementation of squared conjugate gradient method. ...
Definition: Solvers.hpp:54
Definition: KrylovRunner.hpp:32