-
Notifications
You must be signed in to change notification settings - Fork 286
Expand file tree
/
Copy pathlinear_solvers_application.cpp
More file actions
181 lines (148 loc) · 7.1 KB
/
linear_solvers_application.cpp
File metadata and controls
181 lines (148 loc) · 7.1 KB
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
// KRATOS _ _ ____ _
// | | (_)_ __ ___ __ _ _ __/ ___| ___ | |_ _____ _ __ ___
// | | | | '_ \ / _ \/ _` | '__\___ \ / _ \| \ \ / / _ \ '__/ __|
// | |___| | | | | __/ (_| | | ___) | (_) | |\ V / __/ | \__ |
// |_____|_|_| |_|\___|\__,_|_| |____/ \___/|_| \_/ \___|_| |___/ Application
//
// Author: Thomas Oberbichler
// System includes
// External includes
// Project includes
#include "includes/define.h"
#include "includes/registry.h"
#include "linear_solvers_application.h"
#include "custom_factories/dense_linear_solver_factory.h"
#include "custom_solvers/eigen_sparse_cg_solver.h"
#include "custom_solvers/eigen_sparse_lu_solver.h"
#include "custom_solvers/eigen_sparse_qr_solver.h"
#include "custom_solvers/eigen_direct_solver.h"
#include "custom_solvers/eigen_cholmod_solver.hpp" // EigenCholmodSolver
#include "custom_solvers/eigen_spqr_solver.hpp" // EigenSPQRSolver
#include "custom_solvers/eigen_umfpack_solver.hpp" // EigenUmfPackSolver
#if defined USE_EIGEN_MKL
#include "custom_solvers/eigen_pardiso_lu_solver.h"
#include "custom_solvers/eigen_pardiso_llt_solver.h"
#include "custom_solvers/eigen_pardiso_ldlt_solver.h"
#include "custom_solvers/mkl_ilu.hpp" // MKLILU0Smoother, MKLILUTSmoother
#include "mkl_service.h"
#endif
Kratos::KratosApplication* CreateApplication()
{
return new Kratos::KratosLinearSolversApplication();
}
namespace Kratos
{
void KratosLinearSolversApplication::Register()
{
KRATOS_INFO("") << " Kratos _ _ ____ _\n"
<< " | | (_)_ __ ___ __ _ _ __/ ___| ___ | |_ _____ _ __ ___\n"
<< " | | | | '_ \\ / _ \\/ _` | '__\\___ \\ / _ \\| \\ \\ / / _ \\ '__/ __|\n"
<< " | |___| | | | | __/ (_| | | ___) | (_) | |\\ V / __/ | \\__ \\\n"
<< " |_____|_|_| |_|\\___|\\__,_|_| |____/ \\___/|_| \\_/ \\___|_| |___/\n"
<< "Initializing KratosLinearSolversApplication..." << std::endl;
// Adding the eigen library to the list of registered libraries
if (!Registry::HasItem("libraries.eigen")) {
Registry::AddItem<std::string>("libraries.eigen");
}
RegisterDenseLinearSolvers();
using complex = std::complex<double>;
// Sparse LU Solver
using SparseLUType = EigenDirectSolver<EigenSparseLUSolver<double>>;
static auto SparseLUFactory = SparseLUType::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("sparse_lu", SparseLUFactory);
// Complex Sparse LU Solver
using ComplexSparseLUType = EigenDirectSolver<EigenSparseLUSolver<complex>>;
static auto ComplexSparseLUFactory = ComplexSparseLUType::Factory();
KRATOS_REGISTER_COMPLEX_LINEAR_SOLVER("sparse_lu_complex", ComplexSparseLUFactory);
// Sparse QR Solver
using SparseQRType = EigenDirectSolver<EigenSparseQRSolver<double>>;
static auto SparseQRFactory = SparseQRType::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("sparse_qr", SparseQRFactory);
// Sparse CG Solver
using SparseCGType = EigenDirectSolver<EigenSparseCGSolver<double>>;
static auto SparseCGFactory = SparseCGType::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("sparse_cg", SparseCGFactory);
#ifdef USE_EIGEN_MKL
{
MKLVersion mkl_version;
mkl_get_version(&mkl_version);
KRATOS_INFO("") << "Using Intel MKL version "
<< mkl_version.MajorVersion << "." << mkl_version.MinorVersion << "." << mkl_version.UpdateVersion
<< " build " << mkl_version.Build << std::endl
<< "MKL processor: " << mkl_version.Processor << std::endl;
}
if (!Registry::HasItem("libraries.mkl")) {
Registry::AddItem<std::string>("libraries.mkl");
}
// Pardiso LU Solver
using PardisoLUType = EigenDirectSolver<EigenPardisoLUSolver<double>>;
static auto PardisoLUFactory = PardisoLUType::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("pardiso_lu", PardisoLUFactory);
// Pardiso LDLT Solver
using PardisoLDLTType = EigenDirectSolver<EigenPardisoLDLTSolver<double>>;
static auto PardisoLDLTFactory = PardisoLDLTType::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("pardiso_ldlt", PardisoLDLTFactory);
// Pardiso LLT Solver
using PardisoLLTType = EigenDirectSolver<EigenPardisoLLTSolver<double>>;
static auto PardisoLLTFactory = PardisoLLTType::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("pardiso_llt", PardisoLLTFactory);
// Complex Pardiso LU Solver
using ComplexPardisoLUType = EigenDirectSolver<EigenPardisoLUSolver<complex>>;
static auto ComplexPardisoLUFactory = ComplexPardisoLUType::Factory();
KRATOS_REGISTER_COMPLEX_LINEAR_SOLVER("pardiso_lu_complex", ComplexPardisoLUFactory);
// Complex Pardiso LDLT Solver
using ComplexPardisoLDLTType = EigenDirectSolver<EigenPardisoLDLTSolver<complex>>;
static auto ComplexPardisoLDLTFactory = ComplexPardisoLDLTType::Factory();
KRATOS_REGISTER_COMPLEX_LINEAR_SOLVER("pardiso_ldlt_complex", ComplexPardisoLDLTFactory);
// Complex Pardiso LLT Solver
using ComplexPardisoLLTType = EigenDirectSolver<EigenPardisoLLTSolver<complex>>;
static auto ComplexPardisoLLTFactory = ComplexPardisoLLTType::Factory();
KRATOS_REGISTER_COMPLEX_LINEAR_SOLVER("pardiso_llt_complex", ComplexPardisoLLTFactory);
// ILU0 smoother.
static auto mkl_ilu0_factory = StandardLinearSolverFactory<
TUblasSparseSpace<double>,
TUblasDenseSpace<double>,
MKLILU0Smoother<TUblasSparseSpace<double>,TUblasDenseSpace<double>>
>();
KratosComponents<LinearSolverFactory<
TUblasSparseSpace<double>,
TUblasDenseSpace<double>>>::Add("mkl_ilu0", mkl_ilu0_factory);
// ILUT smoother.
static auto mkl_ilut_factory = StandardLinearSolverFactory<
TUblasSparseSpace<double>,
TUblasDenseSpace<double>,
MKLILUTSmoother<TUblasSparseSpace<double>,TUblasDenseSpace<double>>
>();
KratosComponents<LinearSolverFactory<
TUblasSparseSpace<double>,
TUblasDenseSpace<double>>>::Add("mkl_ilut", mkl_ilut_factory);
#endif // defined USE_EIGEN_MKL
#ifdef KRATOS_USE_EIGEN_SUITESPARSE
// Real CHOLMOD.
{
static auto factory = EigenDirectSolver<EigenCholmodSolver<double>>::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("cholmod", factory);
}
// Real SPQR.
{
static auto factory = EigenDirectSolver<EigenSPQRSolver<double>>::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("spqr", factory);
}
// Complex SPQR.
{
static auto factory = EigenDirectSolver<EigenSPQRSolver<complex>>::Factory();
KRATOS_REGISTER_COMPLEX_LINEAR_SOLVER("spqr_complex", factory);
}
// Real UMFPACK.
{
static auto factory = EigenDirectSolver<EigenUmfPackSolver<double>>::Factory();
KRATOS_REGISTER_LINEAR_SOLVER("umfpack", factory);
}
// Complex UMFPACK.
{
static auto factory = EigenDirectSolver<EigenUmfPackSolver<complex>>::Factory();
KRATOS_REGISTER_COMPLEX_LINEAR_SOLVER("umfpack_complex", factory);
}
#endif // KRATOS_USE_EIGEN_SUITESPARSE
}
} // namespace Kratos