aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--cmake/FindKLU.cmake51
-rw-r--r--test/CMakeLists.txt16
-rw-r--r--test/klu_support.cpp32
-rw-r--r--unsupported/Eigen/KLUSupport41
-rw-r--r--unsupported/Eigen/src/KLUSupport/KLUSupport.h364
5 files changed, 504 insertions, 0 deletions
diff --git a/cmake/FindKLU.cmake b/cmake/FindKLU.cmake
new file mode 100644
index 000000000..2783b63d2
--- /dev/null
+++ b/cmake/FindKLU.cmake
@@ -0,0 +1,51 @@
+# KLU lib usually requires linking to a blas library.
+# It is up to the user of this module to find a BLAS and link to it.
+
+if (KLU_INCLUDES AND KLU_LIBRARIES)
+ set(KLU_FIND_QUIETLY TRUE)
+endif (KLU_INCLUDES AND KLU_LIBRARIES)
+
+find_path(KLU_INCLUDES
+ NAMES
+ klu.h
+ PATHS
+ $ENV{KLUDIR}
+ ${INCLUDE_INSTALL_DIR}
+ PATH_SUFFIXES
+ suitesparse
+ ufsparse
+)
+
+if(KLU_LIBRARIES)
+
+ if(NOT KLU_LIBDIR)
+ get_filename_component(KLU_LIBDIR ${KLU_LIBRARIES} PATH)
+ endif(NOT KLU_LIBDIR)
+
+ find_library(COLAMD_LIBRARY colamd PATHS ${KLU_LIBDIR} $ENV{KLUDIR} ${LIB_INSTALL_DIR})
+ if(COLAMD_LIBRARY)
+ set(KLU_LIBRARIES ${KLU_LIBRARIES} ${COLAMD_LIBRARY})
+ endif ()
+
+ find_library(AMD_LIBRARY amd PATHS ${KLU_LIBDIR} $ENV{KLUDIR} ${LIB_INSTALL_DIR})
+ if(AMD_LIBRARY)
+ set(KLU_LIBRARIES ${KLU_LIBRARIES} ${AMD_LIBRARY})
+ endif ()
+
+ find_library(SUITESPARSE_LIBRARY SuiteSparse PATHS ${KLU_LIBDIR} $ENV{KLUDIR} ${LIB_INSTALL_DIR})
+ if(SUITESPARSE_LIBRARY)
+ set(KLU_LIBRARIES ${KLU_LIBRARIES} ${SUITESPARSE_LIBRARY})
+ endif ()
+
+ find_library(CHOLMOD_LIBRARY cholmod PATHS $ENV{KLU_LIBDIR} $ENV{KLUDIR} ${LIB_INSTALL_DIR})
+ if(CHOLMOD_LIBRARY)
+ set(KLU_LIBRARIES ${KLU_LIBRARIES} ${CHOLMOD_LIBRARY})
+ endif()
+
+endif(KLU_LIBRARIES)
+
+include(FindPackageHandleStandardArgs)
+find_package_handle_standard_args(KLU DEFAULT_MSG
+ KLU_INCLUDES KLU_LIBRARIES)
+
+mark_as_advanced(KLU_INCLUDES KLU_LIBRARIES AMD_LIBRARY COLAMD_LIBRARY CHOLMOD_LIBRARY SUITESPARSE_LIBRARY)
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index e73ab92b4..8bd086ce3 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -68,6 +68,17 @@ else()
ei_add_property(EIGEN_MISSING_BACKENDS "UmfPack, ")
endif()
+find_package(KLU)
+if(KLU_FOUND)
+ add_definitions("-DEIGEN_KLU_SUPPORT")
+ include_directories(${KLU_INCLUDES})
+ set(SPARSE_LIBS ${SPARSE_LIBS} ${KLU_LIBRARIES} ${EIGEN_BLAS_LIBRARIES})
+ set(KLU_ALL_LIBS ${KLU_LIBRARIES} ${EIGEN_BLAS_LIBRARIES})
+ ei_add_property(EIGEN_TESTED_BACKENDS "KLU, ")
+else()
+ ei_add_property(EIGEN_MISSING_BACKENDS "KLU, ")
+endif()
+
find_package(SuperLU 4.0)
if(SUPERLU_FOUND)
add_definitions("-DEIGEN_SUPERLU_SUPPORT")
@@ -297,6 +308,11 @@ if(UMFPACK_FOUND)
ei_add_test(umfpack_support "" "${UMFPACK_ALL_LIBS}")
endif()
+if(KLU_FOUND OR SuiteSparse_FOUND)
+ message("ADDING KLU TEST")
+ ei_add_test(klu_support "" "${KLU_ALL_LIBS}")
+endif()
+
if(SUPERLU_FOUND)
ei_add_test(superlu_support "" "${SUPERLU_ALL_LIBS}")
endif()
diff --git a/test/klu_support.cpp b/test/klu_support.cpp
new file mode 100644
index 000000000..8b1fdeb41
--- /dev/null
+++ b/test/klu_support.cpp
@@ -0,0 +1,32 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+
+#define EIGEN_NO_DEBUG_SMALL_PRODUCT_BLOCKS
+#include "sparse_solver.h"
+
+#include <unsupported/Eigen/KLUSupport>
+
+template<typename T> void test_klu_support_T()
+{
+ KLU<SparseMatrix<T, ColMajor> > klu_colmajor;
+ KLU<SparseMatrix<T, RowMajor> > klu_rowmajor;
+
+ check_sparse_square_solving(klu_colmajor);
+ check_sparse_square_solving(klu_rowmajor);
+
+ //check_sparse_square_determinant(umfpack_colmajor);
+ //check_sparse_square_determinant(umfpack_rowmajor);
+}
+
+void test_klu_support()
+{
+ CALL_SUBTEST_1(test_klu_support_T<double>());
+ CALL_SUBTEST_2(test_klu_support_T<std::complex<double> >());
+}
+
diff --git a/unsupported/Eigen/KLUSupport b/unsupported/Eigen/KLUSupport
new file mode 100644
index 000000000..b23d90535
--- /dev/null
+++ b/unsupported/Eigen/KLUSupport
@@ -0,0 +1,41 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+
+#ifndef EIGEN_KLUSUPPORT_MODULE_H
+#define EIGEN_KLUSUPPORT_MODULE_H
+
+#include <Eigen/SparseCore>
+
+#include <Eigen/src/Core/util/DisableStupidWarnings.h>
+
+extern "C" {
+#include <btf.h>
+#include <klu.h>
+ }
+
+/** \ingroup Support_modules
+ * \defgroup KLUSupport_Module KLUSupport module
+ *
+ * This module provides an interface to the KLU library which is part of the <a href="http://www.suitesparse.com">suitesparse</a> package.
+ * It provides the following factorization class:
+ * - class KLU: a sparse LU factorization, well-suited for circuit simulation.
+ *
+ * \code
+ * #include <Eigen/KLUSupport>
+ * \endcode
+ *
+ * In order to use this module, the klu and btf headers must be accessible from the include paths, and your binary must be linked to the klu library and its dependencies.
+ * The dependencies depend on how umfpack has been compiled.
+ * For a cmake based project, you can use our FindKLU.cmake module to help you in this task.
+ *
+ */
+
+#include "src/KLUSupport/KLUSupport.h"
+
+#include <Eigen/src/Core/util/ReenableStupidWarnings.h>
+
+#endif // EIGEN_KLUSUPPORT_MODULE_H
diff --git a/unsupported/Eigen/src/KLUSupport/KLUSupport.h b/unsupported/Eigen/src/KLUSupport/KLUSupport.h
new file mode 100644
index 000000000..d2781202e
--- /dev/null
+++ b/unsupported/Eigen/src/KLUSupport/KLUSupport.h
@@ -0,0 +1,364 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2017 Kyle Macfarlan <kyle.macfarlan@gmail.com>
+//
+// This Source Code Form is subject to the terms of the Mozilla
+// Public License v. 2.0. If a copy of the MPL was not distributed
+// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
+
+#ifndef EIGEN_KLUSUPPORT_H
+#define EIGEN_KLUSUPPORT_H
+
+namespace Eigen {
+
+/* TODO extract L, extract U, compute det, etc... */
+
+/** \ingroup KLUSupport_Module
+ * \brief A sparse LU factorization and solver based on KLU
+ *
+ * This class allows to solve for A.X = B sparse linear problems via a LU factorization
+ * using the KLU library. The sparse matrix A must be squared and full rank.
+ * The vectors or matrices X and B can be either dense or sparse.
+ *
+ * \warning The input matrix A should be in a \b compressed and \b column-major form.
+ * Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
+ * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
+ *
+ * \implsparsesolverconcept
+ *
+ * \sa \ref TutorialSparseSolverConcept, class SparseLU
+ */
+
+
+inline int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, int ldim, int nrhs, double B [ ], klu_common *Common, double) {
+ return klu_solve(Symbolic, Numeric, ldim, nrhs, B, Common);
+}
+
+inline int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, int ldim, int nrhs, std::complex<double>B[], klu_common *Common, std::complex<double>) {
+ return klu_z_solve(Symbolic, Numeric, ldim, nrhs, &numext::real_ref(B[0]), Common);
+}
+
+inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, int ldim, int nrhs, double B[], klu_common *Common, double) {
+ return klu_tsolve(Symbolic, Numeric, ldim, nrhs, B, Common);
+}
+
+inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, int ldim, int nrhs, std::complex<double>B[], klu_common *Common, std::complex<double>) {
+ return klu_z_tsolve(Symbolic, Numeric, ldim, nrhs, &numext::real_ref(B[0]), 0, Common);
+}
+
+inline klu_numeric* klu_factor(int Ap [ ], int Ai [ ], double Ax [ ], klu_symbolic *Symbolic, klu_common *Common, double) {
+ return klu_factor(Ap, Ai, Ax, Symbolic, Common);
+}
+
+inline klu_numeric* klu_factor(int Ap[], int Ai[], std::complex<double> Ax[], klu_symbolic *Symbolic, klu_common *Common, std::complex<double>) {
+ return klu_z_factor(Ap, Ai, &numext::real_ref(Ax[0]), Symbolic, Common);
+}
+
+
+template<typename _MatrixType>
+class KLU : public SparseSolverBase<KLU<_MatrixType> >
+{
+ protected:
+ typedef SparseSolverBase<KLU<_MatrixType> > Base;
+ using Base::m_isInitialized;
+ public:
+ using Base::_solve_impl;
+ typedef _MatrixType MatrixType;
+ typedef typename MatrixType::Scalar Scalar;
+ typedef typename MatrixType::RealScalar RealScalar;
+ typedef typename MatrixType::StorageIndex StorageIndex;
+ typedef Matrix<Scalar,Dynamic,1> Vector;
+ typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
+ typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
+ typedef SparseMatrix<Scalar> LUMatrixType;
+ typedef SparseMatrix<Scalar,ColMajor,int> KLUMatrixType;
+ typedef Ref<const KLUMatrixType, StandardCompressedFormat> KLUMatrixRef;
+ enum {
+ ColsAtCompileTime = MatrixType::ColsAtCompileTime,
+ MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
+ };
+
+ public:
+
+ KLU()
+ : m_dummy(0,0), mp_matrix(m_dummy)
+ {
+ init();
+ }
+
+ template<typename InputMatrixType>
+ explicit KLU(const InputMatrixType& matrix)
+ : mp_matrix(matrix)
+ {
+ init();
+ compute(matrix);
+ }
+
+ ~KLU()
+ {
+ if(m_symbolic) klu_free_symbolic(&m_symbolic,&m_common);
+ if(m_numeric) klu_free_numeric(&m_numeric,&m_common);
+ }
+
+ inline Index rows() const { return mp_matrix.rows(); }
+ inline Index cols() const { return mp_matrix.cols(); }
+
+ /** \brief Reports whether previous computation was successful.
+ *
+ * \returns \c Success if computation was succesful,
+ * \c NumericalIssue if the matrix.appears to be negative.
+ */
+ ComputationInfo info() const
+ {
+ eigen_assert(m_isInitialized && "Decomposition is not initialized.");
+ return m_info;
+ }
+
+ inline const LUMatrixType& matrixL() const
+ {
+ if (m_extractedDataAreDirty) extractData();
+ return m_l;
+ }
+
+ inline const LUMatrixType& matrixU() const
+ {
+ if (m_extractedDataAreDirty) extractData();
+ return m_u;
+ }
+
+ inline const IntColVectorType& permutationP() const
+ {
+ if (m_extractedDataAreDirty) extractData();
+ return m_p;
+ }
+
+ inline const IntRowVectorType& permutationQ() const
+ {
+ if (m_extractedDataAreDirty) extractData();
+ return m_q;
+ }
+
+ /** Computes the sparse Cholesky decomposition of \a matrix
+ * Note that the matrix should be column-major, and in compressed format for best performance.
+ * \sa SparseMatrix::makeCompressed().
+ */
+ template<typename InputMatrixType>
+ void compute(const InputMatrixType& matrix)
+ {
+ if(m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
+ if(m_numeric) klu_free_numeric(&m_numeric, &m_common);
+ grab(matrix.derived());
+ analyzePattern_impl();
+ factorize_impl();
+ }
+
+ /** Performs a symbolic decomposition on the sparcity of \a matrix.
+ *
+ * This function is particularly useful when solving for several problems having the same structure.
+ *
+ * \sa factorize(), compute()
+ */
+ template<typename InputMatrixType>
+ void analyzePattern(const InputMatrixType& matrix)
+ {
+ if(m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
+ if(m_numeric) klu_free_numeric(&m_numeric, &m_common);
+
+ grab(matrix.derived());
+
+ analyzePattern_impl();
+ }
+
+
+ /** Provides access to the control settings array used by KLU.
+ *
+ * See KLU documentation for details.
+ */
+ inline const klu_common& kluCommon() const
+ {
+ return m_common;
+ }
+
+ /** Provides access to the control settings array used by UmfPack.
+ *
+ * If this array contains NaN's, the default values are used.
+ *
+ * See KLU documentation for details.
+ */
+ inline klu_common& kluCommon()
+ {
+ return m_common;
+ }
+
+ /** Performs a numeric decomposition of \a matrix
+ *
+ * The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.
+ *
+ * \sa analyzePattern(), compute()
+ */
+ template<typename InputMatrixType>
+ void factorize(const InputMatrixType& matrix)
+ {
+ eigen_assert(m_analysisIsOk && "KLU: you must first call analyzePattern()");
+ if(m_numeric)
+ klu_free_numeric(&m_numeric,&m_common);
+
+ grab(matrix.derived());
+
+ factorize_impl();
+ }
+
+ /** \internal */
+ template<typename BDerived,typename XDerived>
+ bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
+
+ Scalar determinant() const;
+
+ void extractData() const;
+
+ protected:
+
+ void init()
+ {
+ m_info = InvalidInput;
+ m_isInitialized = false;
+ m_numeric = 0;
+ m_symbolic = 0;
+ m_extractedDataAreDirty = true;
+
+ klu_defaults(&m_common);
+ }
+
+ void analyzePattern_impl()
+ {
+ m_info = InvalidInput;
+ m_analysisIsOk = false;
+ m_factorizationIsOk = false;
+ m_symbolic = klu_analyze(internal::convert_index<int>(mp_matrix.rows()),
+ const_cast<StorageIndex*>(mp_matrix.outerIndexPtr()), const_cast<StorageIndex*>(mp_matrix.innerIndexPtr()),
+ &m_common);
+ if (m_symbolic) {
+ m_isInitialized = true;
+ m_info = Success;
+ m_analysisIsOk = true;
+ m_extractedDataAreDirty = true;
+ }
+ }
+
+ void factorize_impl()
+ {
+
+ m_numeric = klu_factor(const_cast<StorageIndex*>(mp_matrix.outerIndexPtr()), const_cast<StorageIndex*>(mp_matrix.innerIndexPtr()), const_cast<Scalar*>(mp_matrix.valuePtr()),
+ m_symbolic, &m_common, Scalar());
+
+
+ m_info = m_numeric ? Success : NumericalIssue;
+ m_factorizationIsOk = m_numeric ? 1 : 0;
+ m_extractedDataAreDirty = true;
+ }
+
+ template<typename MatrixDerived>
+ void grab(const EigenBase<MatrixDerived> &A)
+ {
+ mp_matrix.~KLUMatrixRef();
+ ::new (&mp_matrix) KLUMatrixRef(A.derived());
+ }
+
+ void grab(const KLUMatrixRef &A)
+ {
+ if(&(A.derived()) != &mp_matrix)
+ {
+ mp_matrix.~KLUMatrixRef();
+ ::new (&mp_matrix) KLUMatrixRef(A);
+ }
+ }
+
+ // cached data to reduce reallocation, etc.
+ mutable LUMatrixType m_l;
+
+ mutable LUMatrixType m_u;
+ mutable IntColVectorType m_p;
+ mutable IntRowVectorType m_q;
+
+ KLUMatrixType m_dummy;
+ KLUMatrixRef mp_matrix;
+
+ klu_numeric* m_numeric;
+ klu_symbolic* m_symbolic;
+ klu_common m_common;
+ mutable ComputationInfo m_info;
+ int m_factorizationIsOk;
+ int m_analysisIsOk;
+ mutable bool m_extractedDataAreDirty;
+
+ private:
+ KLU(const KLU& ) { }
+};
+
+
+template<typename MatrixType>
+void KLU<MatrixType>::extractData() const
+{
+ if (m_extractedDataAreDirty)
+ {
+ eigen_assert(false && "KLU: extractData Not Yet Implemented");
+
+// // get size of the data
+// int lnz, unz, rows, cols, nz_udiag;
+// umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
+//
+// // allocate data
+// m_l.resize(rows,(std::min)(rows,cols));
+// m_l.resizeNonZeros(lnz);
+//
+// m_u.resize((std::min)(rows,cols),cols);
+// m_u.resizeNonZeros(unz);
+//
+// m_p.resize(rows);
+// m_q.resize(cols);
+//
+// // extract
+// umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
+// m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
+// m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
+//
+// m_extractedDataAreDirty = false;
+ }
+}
+
+template<typename MatrixType>
+typename KLU<MatrixType>::Scalar KLU<MatrixType>::determinant() const
+{
+ eigen_assert(false && "KLU: extractData Not Yet Implemented");
+ return Scalar();
+}
+
+template<typename MatrixType>
+template<typename BDerived,typename XDerived>
+bool KLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
+{
+ Index rhsCols = b.cols();
+ eigen_assert((BDerived::Flags&RowMajorBit)==0 && "KLU backend does not support non col-major rhs yet");
+ eigen_assert((XDerived::Flags&RowMajorBit)==0 && "KLU backend does not support non col-major result yet");
+ eigen_assert(b.derived().data() != x.derived().data() && " KLU does not support inplace solve");
+ eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
+
+ x = b;
+ int info = 0;
+ if (true/*(MatrixType::Flags&RowMajorBit) == 0*/)
+ {
+ info = klu_solve(m_symbolic, m_numeric, b.rows(), rhsCols, x.const_cast_derived().data(), const_cast<klu_common*>(&m_common), Scalar());
+ }
+ else
+ {
+ info = klu_tsolve(m_symbolic, m_numeric, b.rows(), rhsCols, x.const_cast_derived().data(), const_cast<klu_common*>(&m_common), Scalar());
+ }
+
+ m_info = info!=0 ? Success : NumericalIssue;
+ return true;
+}
+
+} // end namespace Eigen
+
+#endif // EIGEN_KLUSUPPORT_H