Skip to content

Commit d96792b

Browse files
committed
adding suitesparse::spqr
1 parent 760ff26 commit d96792b

File tree

3 files changed

+50
-0
lines changed

3 files changed

+50
-0
lines changed

CMakeLists.txt

Lines changed: 12 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -74,6 +74,7 @@ option(POLYSOLVE_WITH_ACCELERATE "Enable Apple Accelerate" ${POLYSOLVE_ON_APP
7474
option(POLYSOLVE_WITH_CHOLMOD "Enable Cholmod library" ON)
7575
option(POLYSOLVE_WITH_UMFPACK "Enable UmfPack library" ON)
7676
option(POLYSOLVE_WITH_SUPERLU "Enable SuperLU library" ON)
77+
option(POLYSOLVE_WITH_SPQR "Enable SPQR library" ON)
7778
option(POLYSOLVE_WITH_MKL "Enable MKL library" ${POLYSOLVE_NOT_ON_APPLE_SILICON})
7879
option(POLYSOLVE_WITH_CUSOLVER "Enable cuSOLVER library" OFF)
7980
option(POLYSOLVE_WITH_PARDISO "Enable Pardiso library" OFF)
@@ -269,6 +270,17 @@ if(POLYSOLVE_WITH_SUPERLU)
269270
endif()
270271
endif()
271272

273+
# SuperLU solver
274+
if(POLYSOLVE_WITH_SPQR)
275+
include(SPQR)
276+
if(TARGET SuiteSparse::SPQR)
277+
target_link_libraries(polysolve_linear PRIVATE SuiteSparse::SPQR)
278+
target_compile_definitions(polysolve_linear PUBLIC POLYSOLVE_WITH_SPQR)
279+
else()
280+
message(WARNING "SPQR Not found, solver will not be available.")
281+
endif()
282+
endif()
283+
272284
# AMGCL solver
273285
if(POLYSOLVE_WITH_AMGCL)
274286
include(amgcl)

cmake/recipes/spqr.cmake

Lines changed: 11 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,11 @@
1+
# SPQR solver
2+
3+
if(TARGET SparseSuite::SPQR)
4+
return()
5+
endif()
6+
7+
message(STATUS "Third-party: creating targets 'SuiteSparse::SPQR'")
8+
9+
# We do not have a build recipe for this, so find it as a system installed library.
10+
find_package(SPQR)
11+

src/polysolve/linear/Solver.cpp

Lines changed: 27 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -27,6 +27,24 @@ template <typename _MatrixType, typename _OrderingType = COLAMDOrdering<typename
2727
#ifdef POLYSOLVE_WITH_UMFPACK
2828
#include <Eigen/UmfPackSupport>
2929
#endif
30+
#ifdef POLYSOLVE_WITH_SPQR
31+
#include <Eigen/SPQRSupport>
32+
namespace polysolve::linear {
33+
template <>
34+
void EigenDirect<Eigen::SPQR<StiffnessMatrix>>::analyze_pattern(const StiffnessMatrix& A, const int precond_num) {
35+
m_Solver.compute(A);
36+
}
37+
template <>
38+
void EigenDirect<Eigen::SPQR<StiffnessMatrix>>::factorize(const StiffnessMatrix &A)
39+
{
40+
m_Solver.compute(A);
41+
if (m_Solver.info() == Eigen::NumericalIssue)
42+
{
43+
throw std::runtime_error("[EigenDirect] NumericalIssue encountered.");
44+
}
45+
}
46+
}
47+
#endif
3048
#ifdef POLYSOLVE_WITH_SUPERLU
3149
#include <Eigen/SuperLUSupport>
3250
#endif
@@ -345,6 +363,12 @@ namespace polysolve::linear
345363
{
346364
RETURN_DIRECT_SOLVER_PTR(SuperLU, "Eigen::SuperLU");
347365
#endif
366+
#ifdef POLYSOLVE_WITH_SPQR
367+
}
368+
else if (solver == "Eigen::SPQR")
369+
{
370+
RETURN_DIRECT_SOLVER_PTR(SPQR, "Eigen::SPQR");
371+
#endif
348372
#ifdef POLYSOLVE_WITH_MKL
349373
}
350374
else if (solver == "Eigen::PardisoLLT")
@@ -492,6 +516,9 @@ namespace polysolve::linear
492516
#ifdef POLYSOLVE_WITH_SUPERLU
493517
"Eigen::SuperLU",
494518
#endif
519+
#ifdef POLYSOLVE_WITH_SPQR
520+
"Eigen::SPQR",
521+
#endif
495522
#ifdef POLYSOLVE_WITH_MKL
496523
"Eigen::PardisoLLT",
497524
"Eigen::PardisoLDLT",

0 commit comments

Comments
 (0)