AMDiS  2.10
The Adaptive Multi-Dimensional Simulation Toolbox
Solvers.hpp
1 #pragma once
2 
3 #include <memory>
4 
5 #include <dune/common/classname.hh>
6 #include <dune/common/version.hh>
7 
8 #include <dune/common/ftraits.hh>
9 
10 #include <amdis/CreatorMap.hpp>
11 #include <amdis/Output.hpp>
12 
13 #include <amdis/linearalgebra/istl/ISTLSolverCreator.hpp>
14 #include <amdis/linearalgebra/istl/precompiled/Solvers.hpp>
15 
16 namespace AMDiS
17 {
19 
35  template <class Traits>
36  class DefaultCreators<tag::solver<Traits>>
37  {
38  using M = typename Traits::M;
39  using X = typename Traits::X;
40  using Y = typename Traits::Y;
41 
42  using FTraits = Dune::FieldTraits<typename M::field_type>;
43 
44  template <class S>
45  using Creator = ISTLSolverCreator<S,Traits>;
46 
48 
49  public:
50  static void init()
51  {
52  auto cg = new Creator<Dune::CGSolver<X>>;
53  Map::addCreator("cg", cg);
54 
55  // Generalized preconditioned conjugate gradient solver.
56  auto pcg = new Creator<Dune::GeneralizedPCGSolver<X>>;
57  Map::addCreator("pcg", pcg);
58 
59  auto fcg = new Creator<Dune::RestartedFCGSolver<X>>;
60  Map::addCreator("fcg", fcg);
61 
62  auto cfcg = new Creator<Dune::CompleteFCGSolver<X>>;
63  Map::addCreator("cfcg", cfcg);
64 
65  auto bicgstab = new Creator<Dune::BiCGSTABSolver<X>>;
66  Map::addCreator("bicgstab", bicgstab);
67  Map::addCreator("bcgs", bicgstab);
68  Map::addCreator("default", bicgstab);
69 
70  auto minres = new Creator<Dune::MINRESSolver<X>>;
71  Map::addCreator("minres", minres);
72 
73  auto gmres = new Creator<Dune::RestartedGMResSolver<X,Y>>;
74  Map::addCreator("gmres", gmres);
75 
76  auto fgmres = new Creator<Dune::RestartedFlexibleGMResSolver<X,Y>>;
77  Map::addCreator("fgmres", fgmres);
78 
79  init_direct(std::is_same<typename FTraits::real_type, double>{});
80  }
81 
82  static void init_direct(std::false_type)
83  {
84  warning("Direct solvers not created for the matrix with real_type = {}.",
85  Dune::className<typename FTraits::real_type>());
86  }
87 
88  static void init_direct(std::true_type)
89  {
90 #if HAVE_SUITESPARSE_UMFPACK
91  auto umfpack = new Creator<Dune::UMFPack<M>>;
92  Map::addCreator("umfpack", umfpack);
93 #endif
94 
95 #if HAVE_SUITESPARSE_LDL
96  auto ldl = new Creator<Dune::LDL<M>>;
97  Map::addCreator("ldl", ldl);
98 #endif
99 
100 #if HAVE_SUITESPARSE_SPQR
101  auto spqr = new Creator<Dune::SPQR<M>>;
102  Map::addCreator("spqr", spqr);
103 #endif
104 
105 #if HAVE_SUPERLU
106  auto superlu = new Creator<Dune::SuperLU<M>>;
107  Map::addCreator("superlu", superlu);
108 #endif
109 
110  // default direct solvers
111 #if HAVE_SUITESPARSE_UMFPACK
112  Map::addCreator("direct", umfpack);
113 #elif HAVE_SUPERLU
114  Map::addCreator("direct", superlu);
115 #endif
116  }
117  };
118 
119 
120  // extern template declarations:
122 
123 } // end namespace AMDiS
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
Definition: CreatorMap.hpp:16