-
Notifications
You must be signed in to change notification settings - Fork 20
/
Copy pathLinearSolver.cpp
474 lines (443 loc) · 19.8 KB
/
LinearSolver.cpp
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
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
////////////////////////////////////////////////////////////////////////////////
#include <polysolve/LinearSolver.hpp>
#include <polysolve/LinearSolverEigen.hpp>
#include <polysolve/SaddlePointSolver.hpp>
// -----------------------------------------------------------------------------
#include <Eigen/Sparse>
#ifdef POLYSOLVE_WITH_SYMPILER
#include <polysolve/LinearSolverSympiler.hpp>
#endif
#ifdef POLYSOLVE_WITH_ACCELERATE
#include <Eigen/AccelerateSupport>
#endif
#ifdef POLYSOLVE_WITH_CHOLMOD
#include <Eigen/CholmodSupport>
#endif
#ifdef POLYSOLVE_WITH_UMFPACK
#include <Eigen/UmfPackSupport>
#endif
#ifdef POLYSOLVE_WITH_SUPERLU
#include <Eigen/SuperLUSupport>
#endif
#ifdef POLYSOLVE_WITH_MKL
#include <Eigen/PardisoSupport>
#endif
#ifdef POLYSOLVE_WITH_PARDISO
#include <polysolve/LinearSolverPardiso.hpp>
#endif
#ifdef POLYSOLVE_WITH_HYPRE
#include <polysolve/LinearSolverHypre.hpp>
#endif
#ifdef POLYSOLVE_WITH_AMGCL
#include <polysolve/LinearSolverAMGCL.hpp>
#endif
#ifdef POLYSOLVE_WITH_CUSOLVER
#include <polysolve/LinearSolverCuSolverDN.cuh>
#endif
#include <unsupported/Eigen/IterativeSolvers>
////////////////////////////////////////////////////////////////////////////////
namespace polysolve
{
////////////////////////////////////////////////////////////////////////////////
#if EIGEN_VERSION_AT_LEAST(3, 3, 0)
// Magic macro because C++ has no introspection
#define ENUMERATE_PRECOND(HelperFunctor, SolverType, DefaultPrecond, precond, name) \
do \
{ \
using namespace Eigen; \
if (precond == "Eigen::IdentityPreconditioner") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
IdentityPreconditioner>::type>(name); \
} \
else if (precond == "Eigen::DiagonalPreconditioner") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
DiagonalPreconditioner<double>>::type>(name); \
} \
else if (precond == "Eigen::IncompleteCholesky") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
IncompleteCholesky<double>>::type>(name); \
} \
else if (precond == "Eigen::LeastSquareDiagonalPreconditioner") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
LeastSquareDiagonalPreconditioner<double>>::type>(name); \
} \
else if (precond == "Eigen::IncompleteLUT") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
IncompleteLUT<double>>::type>(name); \
} \
else \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
DefaultPrecond>::type>(name); \
} \
} while (0)
#else
// Magic macro because C++ has no introspection
#define ENUMERATE_PRECOND(HelperFunctor, SolverType, DefaultPrecond, precond) \
do \
{ \
using namespace Eigen; \
if (precond == "Eigen::IdentityPreconditioner") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
IdentityPreconditioner>::type>(); \
} \
else if (precond == "Eigen::DiagonalPreconditioner") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
DiagonalPreconditioner<double>>::type>(); \
} \
else if (precond == "Eigen::IncompleteCholesky") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
IncompleteCholesky<double>>::type>(); \
} \
else if (precond == "Eigen::IncompleteLUT") \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
IncompleteLUT<double>>::type>(); \
} \
else \
{ \
return std::make_unique<typename HelperFunctor<SolverType, \
DefaultPrecond>::type>(); \
} \
} while (0)
#endif
// -----------------------------------------------------------------------------
#define RETURN_DIRECT_SOLVER_PTR(EigenSolver, Name) \
do \
{ \
return std::make_unique<LinearSolverEigenDirect<EigenSolver< \
polysolve::StiffnessMatrix>>>(Name); \
} while (0)
#define RETURN_DIRECT_DENSE_SOLVER_PTR(EigenSolver, Name) \
do \
{ \
return std::make_unique<LinearSolverEigenDense<EigenSolver< \
Eigen::MatrixXd>>>(Name); \
} while (0)
////////////////////////////////////////////////////////////////////////////////
namespace
{
template <template <class, class> class SparseSolver, typename Precond>
struct MakeSolver
{
typedef LinearSolverEigenIterative<SparseSolver<StiffnessMatrix, Precond>> type;
};
template <template <class, int, class> class SparseSolver, typename Precond>
struct MakeSolverSym
{
typedef LinearSolverEigenIterative<SparseSolver<StiffnessMatrix,
Eigen::Lower | Eigen::Upper, Precond>>
type;
};
// -----------------------------------------------------------------------------
template <
template <class, class> class SolverType,
typename DefaultPrecond = Eigen::DiagonalPreconditioner<double>>
struct PrecondHelper
{
static std::unique_ptr<LinearSolver> create(const std::string &arg, const std::string &name)
{
ENUMERATE_PRECOND(MakeSolver, SolverType, DefaultPrecond, arg, name);
}
};
template <
template <class, int, class> class SolverType,
typename DefaultPrecond = Eigen::DiagonalPreconditioner<double>>
struct PrecondHelperSym
{
static std::unique_ptr<LinearSolver> create(const std::string &arg, const std::string &name)
{
ENUMERATE_PRECOND(MakeSolverSym, SolverType, DefaultPrecond, arg, name);
}
};
} // anonymous namespace
////////////////////////////////////////////////////////////////////////////////
// Static constructor
std::unique_ptr<LinearSolver> LinearSolver::create(const std::string &solver, const std::string &precond)
{
using namespace Eigen;
if (solver.empty() || solver == "Eigen::SimplicialLDLT")
{
RETURN_DIRECT_SOLVER_PTR(SimplicialLDLT, "Eigen::SimplicialLDLT");
}
else if (solver == "Eigen::SparseLU")
{
RETURN_DIRECT_SOLVER_PTR(SparseLU, "Eigen::SparseLU");
#ifdef POLYSOLVE_WITH_SYMPILER
}
else if (solver == "Sympiler")
{
return std::make_unique<LinearSolverSympiler>();
#endif
#ifdef POLYSOLVE_WITH_ACCELERATE
}
else if (solver == "Eigen::AccelerateLLT")
{
RETURN_DIRECT_SOLVER_PTR(AccelerateLLT, "Eigen::AccelerateLLT");
}
else if (solver == "Eigen::AccelerateLDLT")
{
RETURN_DIRECT_SOLVER_PTR(AccelerateLDLT, "Eigen::AccelerateLDLT");
#endif
#ifdef POLYSOLVE_WITH_CHOLMOD
}
else if (solver == "Eigen::CholmodSupernodalLLT")
{
RETURN_DIRECT_SOLVER_PTR(CholmodSupernodalLLT, "Eigen::CholmodSupernodalLLT");
}
else if (solver == "Eigen::CholmodDecomposition")
{
RETURN_DIRECT_SOLVER_PTR(CholmodDecomposition, "Eigen::CholmodDecomposition");
}
else if (solver == "Eigen::CholmodSimplicialLLT")
{
RETURN_DIRECT_SOLVER_PTR(CholmodSimplicialLLT, "Eigen::CholmodSimplicialLLT");
}
else if (solver == "Eigen::CholmodSimplicialLDLT")
{
RETURN_DIRECT_SOLVER_PTR(CholmodSimplicialLDLT, "Eigen::CholmodSimplicialLDLT");
#endif
#ifdef POLYSOLVE_WITH_UMFPACK
#ifndef POLYSOLVE_LARGE_INDEX
}
else if (solver == "Eigen::UmfPackLU")
{
RETURN_DIRECT_SOLVER_PTR(UmfPackLU, "Eigen::UmfPackLU");
#endif
#endif
#ifdef POLYSOLVE_WITH_SUPERLU
}
else if (solver == "Eigen::SuperLU")
{
RETURN_DIRECT_SOLVER_PTR(SuperLU, "Eigen::SuperLU");
#endif
#ifdef POLYSOLVE_WITH_MKL
}
else if (solver == "Eigen::PardisoLLT")
{
RETURN_DIRECT_SOLVER_PTR(PardisoLLT, "Eigen::PardisoLLT");
}
else if (solver == "Eigen::PardisoLDLT")
{
RETURN_DIRECT_SOLVER_PTR(PardisoLDLT, "Eigen::PardisoLDLT");
}
else if (solver == "Eigen::PardisoLU")
{
RETURN_DIRECT_SOLVER_PTR(PardisoLU, "Eigen::PardisoLU");
#endif
#ifdef POLYSOLVE_WITH_PARDISO
}
else if (solver == "Pardiso")
{
return std::make_unique<LinearSolverPardiso>();
#endif
#ifdef POLYSOLVE_WITH_CUSOLVER
}
else if (solver == "cuSolverDN")
{
return std::make_unique<LinearSolverCuSolverDN<double>>();
}
else if (solver == "cuSolverDN_float")
{
return std::make_unique<LinearSolverCuSolverDN<float>>();
#endif
#ifdef POLYSOLVE_WITH_HYPRE
}
else if (solver == "Hypre")
{
return std::make_unique<LinearSolverHypre>();
#endif
#ifdef POLYSOLVE_WITH_AMGCL
}
else if (solver == "AMGCL")
{
return std::make_unique<LinearSolverAMGCL>();
#endif
#if EIGEN_VERSION_AT_LEAST(3, 3, 0)
// Available only with Eigen 3.3.0 and newer
#ifndef POLYSOLVE_LARGE_INDEX
}
else if (solver == "Eigen::LeastSquaresConjugateGradient")
{
return PrecondHelper<BiCGSTAB, LeastSquareDiagonalPreconditioner<double>>::create(precond, "Eigen::LeastSquaresConjugateGradient");
}
else if (solver == "Eigen::DGMRES")
{
return PrecondHelper<DGMRES>::create(precond, "Eigen::DGMRES");
#endif
#endif
#ifndef POLYSOLVE_LARGE_INDEX
}
else if (solver == "Eigen::ConjugateGradient")
{
return PrecondHelperSym<ConjugateGradient>::create(precond, "Eigen::ConjugateGradient");
}
else if (solver == "Eigen::BiCGSTAB")
{
return PrecondHelper<BiCGSTAB>::create(precond, "Eigen::BiCGSTAB");
}
else if (solver == "Eigen::GMRES")
{
return PrecondHelper<GMRES>::create(precond, "Eigen::GMRES");
}
else if (solver == "Eigen::MINRES")
{
return PrecondHelperSym<MINRES>::create(precond, "Eigen::MINRES");
#endif
}
else if (solver == "SaddlePointSolver")
{
return std::make_unique<SaddlePointSolver>();
}
/////DENSE Eigen
else if (solver.empty() || solver == "Eigen::PartialPivLU")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(PartialPivLU, "Eigen::PartialPivLU");
}
else if (solver.empty() || solver == "Eigen::FullPivLU")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(FullPivLU, "Eigen::FullPivLU");
}
else if (solver.empty() || solver == "Eigen::HouseholderQR")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(HouseholderQR, "Eigen::HouseholderQR");
}
else if (solver.empty() || solver == "Eigen::ColPivHouseholderQR")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(ColPivHouseholderQR, "Eigen::ColPivHouseholderQR");
}
else if (solver.empty() || solver == "Eigen::FullPivHouseholderQR")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(FullPivHouseholderQR, "Eigen::FullPivHouseholderQR");
}
else if (solver.empty() || solver == "Eigen::CompleteOrthogonalDecomposition")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(CompleteOrthogonalDecomposition, "Eigen::CompleteOrthogonalDecomposition");
}
else if (solver.empty() || solver == "Eigen::LLT")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(LLT, "Eigen::LLT");
}
else if (solver.empty() || solver == "Eigen::LDLT")
{
RETURN_DIRECT_DENSE_SOLVER_PTR(LDLT, "Eigen::LDLT");
}
// else if (solver.empty() || solver == "Eigen::BDCSVD")
// {
// RETURN_DIRECT_DENSE_SOLVER_PTR(BDCSVD, "Eigen::BDCSVD");
// }
// else if (solver.empty() || solver == "Eigen::JacobiSVD")
// {
// RETURN_DIRECT_DENSE_SOLVER_PTR(JacobiSVD, "Eigen::JacobiSVD");
// }
throw std::runtime_error("Unrecognized solver type: " + solver);
}
////////////////////////////////////////////////////////////////////////////////
// List available solvers
std::vector<std::string> LinearSolver::availableSolvers()
{
return {{
"Eigen::SimplicialLDLT",
"Eigen::SparseLU",
#ifdef POLYSOLVE_WITH_SYMPILER
"Sympiler",
#endif
#ifdef POLYSOLVE_WITH_ACCELERATE
"Eigen::AccelerateLLT",
"Eigen::AccelerateLDLT",
#endif
#ifdef POLYSOLVE_WITH_CHOLMOD
"Eigen::CholmodSupernodalLLT",
"Eigen::CholmodDecomposition",
"Eigen::CholmodSimplicialLLT",
"Eigen::CholmodSimplicialLDLT",
#endif
#ifdef POLYSOLVE_WITH_UMFPACK
"Eigen::UmfPackLU",
#endif
#ifdef POLYSOLVE_WITH_SUPERLU
"Eigen::SuperLU",
#endif
#ifdef POLYSOLVE_WITH_MKL
"Eigen::PardisoLLT",
"Eigen::PardisoLDLT",
"Eigen::PardisoLU",
#endif
#ifdef POLYSOLVE_WITH_PARDISO
"Pardiso",
#endif
#ifdef POLYSOLVE_WITH_CUSOLVER
"cuSolverDN",
"cuSolverDN_float",
#endif
#ifdef POLYSOLVE_WITH_HYPRE
"Hypre",
#endif
#ifdef POLYSOLVE_WITH_AMGCL
"AMGCL",
#endif
#if EIGEN_VERSION_AT_LEAST(3, 3, 0)
#ifndef POLYSOLVE_LARGE_INDEX
"Eigen::LeastSquaresConjugateGradient",
"Eigen::DGMRES",
#endif
#endif
"Eigen::ConjugateGradient",
"Eigen::BiCGSTAB",
"Eigen::GMRES",
"Eigen::MINRES",
"Eigen::PartialPivLU",
"Eigen::FullPivLU",
"Eigen::HouseholderQR",
"Eigen::ColPivHouseholderQR",
"Eigen::FullPivHouseholderQR",
"Eigen::CompleteOrthogonalDecomposition",
"Eigen::LLT",
"Eigen::LDLT"
// "Eigen::BDCSVD",
// "Eigen::JacobiSVD"
}};
}
std::string LinearSolver::defaultSolver()
{
// return "Eigen::BiCGSTAB";
#ifdef POLYSOLVE_WITH_PARDISO
return "Pardiso";
#else
#ifdef POLYSOLVE_WITH_HYPRE
return "Hypre";
#else
return "Eigen::BiCGSTAB";
#endif
#endif
}
// -----------------------------------------------------------------------------
// List available preconditioners
std::vector<std::string> LinearSolver::availablePrecond()
{
return {{
"Eigen::IdentityPreconditioner",
"Eigen::DiagonalPreconditioner",
"Eigen::IncompleteCholesky",
#if EIGEN_VERSION_AT_LEAST(3, 3, 0)
"Eigen::LeastSquareDiagonalPreconditioner",
#endif
#ifndef POLYSOLVE_LARGE_INDEX
"Eigen::IncompleteLUT",
#endif
}};
}
std::string LinearSolver::defaultPrecond()
{
return "Eigen::DiagonalPreconditioner";
}
////////////////////////////////////////////////////////////////////////////////
} // namespace polysolve