AMDiS 2.11-git
The Adaptive Multi-Dimensional Simulation Toolbox
 
Loading...
Searching...
No Matches
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#include <dune/common/version.hh>
10
11#include <dune/istl/solvers.hh>
12
13#if HAVE_SUITESPARSE_UMFPACK
14#include <dune/istl/umfpack.hh>
15#endif
16#if HAVE_SUITESPARSE_LDL
17#include <dune/istl/ldl.hh>
18#endif
19#if HAVE_SUITESPARSE_SPQR
20#include <dune/istl/spqr.hh>
21#endif
22#if HAVE_SUPERLU
23#include <dune/istl/superlu.hh>
24#endif
25
26#include <amdis/CreatorMap.hpp>
27#include <amdis/Output.hpp>
28
29#include <amdis/linearalgebra/istl/ISTLSolverCreator.hpp>
30
31namespace AMDiS
32{
34
50 template <class Traits>
51 class DefaultCreators<tag::solver<Traits>>
52 {
53 using M = typename Traits::M;
54 using X = typename Traits::X;
55 using Y = typename Traits::Y;
56
57 using FTraits = Dune::FieldTraits<typename M::field_type>;
58
59 template <class S>
60 using Creator = ISTLSolverCreator<S,Traits>;
61
63
64 public:
65 static void init()
66 {
67 auto cg = new Creator<Dune::CGSolver<X>>;
68 Map::addCreator("cg", cg);
69
70 // Generalized preconditioned conjugate gradient solver.
71 auto pcg = new Creator<Dune::GeneralizedPCGSolver<X>>;
72 Map::addCreator("pcg", pcg);
73
74 auto fcg = new Creator<Dune::RestartedFCGSolver<X>>;
75 Map::addCreator("fcg", fcg);
76
77 auto cfcg = new Creator<Dune::CompleteFCGSolver<X>>;
78 Map::addCreator("cfcg", cfcg);
79
80 auto bicgstab = new Creator<Dune::BiCGSTABSolver<X>>;
81 Map::addCreator("bicgstab", bicgstab);
82 Map::addCreator("bcgs", bicgstab);
83 Map::addCreator("default", bicgstab);
84
85 auto minres = new Creator<Dune::MINRESSolver<X>>;
86 Map::addCreator("minres", minres);
87
88 auto gmres = new Creator<Dune::RestartedGMResSolver<X,Y>>;
89 Map::addCreator("gmres", gmres);
90
91 auto fgmres = new Creator<Dune::RestartedFlexibleGMResSolver<X,Y>>;
92 Map::addCreator("fgmres", fgmres);
93
94 init_direct(std::is_same<typename FTraits::real_type, double>{});
95 }
96
97 static void init_direct(std::false_type)
98 {
99 warning("Direct solvers not created for the matrix with real_type = {}.",
100 Dune::className<typename FTraits::real_type>());
101 }
102
103 static void init_direct(std::true_type)
104 {
105#if HAVE_SUITESPARSE_UMFPACK
106 auto umfpack = new Creator<Dune::UMFPack<M>>;
107 Map::addCreator("umfpack", umfpack);
108#endif
109
110#if HAVE_SUITESPARSE_LDL
111 auto ldl = new Creator<Dune::LDL<M>>;
112 Map::addCreator("ldl", ldl);
113#endif
114
115#if HAVE_SUITESPARSE_SPQR
116 auto spqr = new Creator<Dune::SPQR<M>>;
117 Map::addCreator("spqr", spqr);
118#endif
119
120#if HAVE_SUPERLU
121 auto superlu = new Creator<Dune::SuperLU<M>>;
122 Map::addCreator("superlu", superlu);
123#endif
124
125 // default direct solvers
126#if HAVE_SUITESPARSE_UMFPACK
127 Map::addCreator("direct", umfpack);
128#elif HAVE_SUPERLU
129 Map::addCreator("direct", superlu);
130#endif
131 }
132 };
133
134} // end namespace AMDiS
A CreatorMap is used to construct objects, which types depends on key words determined at run time....
Definition CreatorMap.hpp:30
Definition CreatorMap.hpp:16