10#ifndef EIGEN_KLUSUPPORT_H
11#define EIGEN_KLUSUPPORT_H
14#include "./InternalHeaderCheck.h"
37 klu_common *Common,
double) {
38 return klu_solve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), B,
42inline int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, std::complex<double> B[],
43 klu_common *Common, std::complex<double>) {
44 return klu_z_solve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs),
45 &numext::real_ref(B[0]), Common);
48inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs,
double B[],
49 klu_common *Common,
double) {
50 return klu_tsolve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), B,
54inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, std::complex<double> B[],
55 klu_common *Common, std::complex<double>) {
56 return klu_z_tsolve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs),
57 &numext::real_ref(B[0]), 0, Common);
60inline klu_numeric *klu_factor(
int Ap[],
int Ai[],
double Ax[], klu_symbolic *Symbolic, klu_common *Common,
double) {
61 return klu_factor(Ap, Ai, Ax, Symbolic, Common);
64inline klu_numeric *klu_factor(
int Ap[],
int Ai[], std::complex<double> Ax[], klu_symbolic *Symbolic,
65 klu_common *Common, std::complex<double>) {
66 return klu_z_factor(Ap, Ai, &numext::real_ref(Ax[0]), Symbolic, Common);
69template <
typename MatrixType_>
70class KLU :
public SparseSolverBase<KLU<MatrixType_> > {
73 using Base::m_isInitialized;
76 using Base::_solve_impl;
77 typedef MatrixType_ MatrixType;
78 typedef typename MatrixType::Scalar Scalar;
79 typedef typename MatrixType::RealScalar RealScalar;
80 typedef typename MatrixType::StorageIndex StorageIndex;
81 typedef Matrix<Scalar, Dynamic, 1> Vector;
82 typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
83 typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
84 typedef SparseMatrix<Scalar> LUMatrixType;
85 typedef SparseMatrix<Scalar, ColMajor, int> KLUMatrixType;
86 typedef Ref<const KLUMatrixType, StandardCompressedFormat> KLUMatrixRef;
87 enum { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime };
90 KLU() : m_dummy(0, 0), mp_matrix(m_dummy) { init(); }
92 template <
typename InputMatrixType>
93 explicit KLU(
const InputMatrixType &matrix) : mp_matrix(matrix) {
99 if (m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
100 if (m_numeric) klu_free_numeric(&m_numeric, &m_common);
103 EIGEN_CONSTEXPR
inline Index rows() const EIGEN_NOEXCEPT {
return mp_matrix.rows(); }
104 EIGEN_CONSTEXPR
inline Index cols() const EIGEN_NOEXCEPT {
return mp_matrix.cols(); }
111 ComputationInfo info()
const {
112 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
116 inline const LUMatrixType& matrixL()
const
118 if (m_extractedDataAreDirty) extractData();
122 inline const LUMatrixType& matrixU()
const
124 if (m_extractedDataAreDirty) extractData();
128 inline const IntColVectorType& permutationP()
const
130 if (m_extractedDataAreDirty) extractData();
134 inline const IntRowVectorType& permutationQ()
const
136 if (m_extractedDataAreDirty) extractData();
144 template <
typename InputMatrixType>
145 void compute(
const InputMatrixType &matrix) {
146 if (m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
147 if (m_numeric) klu_free_numeric(&m_numeric, &m_common);
148 grab(matrix.derived());
149 analyzePattern_impl();
159 template <
typename InputMatrixType>
160 void analyzePattern(
const InputMatrixType &matrix) {
161 if (m_symbolic) klu_free_symbolic(&m_symbolic, &m_common);
162 if (m_numeric) klu_free_numeric(&m_numeric, &m_common);
164 grab(matrix.derived());
166 analyzePattern_impl();
173 inline const klu_common &kluCommon()
const {
return m_common; }
181 inline klu_common &kluCommon() {
return m_common; }
189 template <
typename InputMatrixType>
190 void factorize(
const InputMatrixType &matrix) {
191 eigen_assert(m_analysisIsOk &&
"KLU: you must first call analyzePattern()");
192 if (m_numeric) klu_free_numeric(&m_numeric, &m_common);
194 grab(matrix.derived());
200 template <
typename BDerived,
typename XDerived>
201 bool _solve_impl(
const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x)
const;
204 Scalar determinant()
const;
206 void extractData()
const;
212 m_isInitialized =
false;
215 m_extractedDataAreDirty =
true;
217 klu_defaults(&m_common);
220 void analyzePattern_impl() {
222 m_analysisIsOk =
false;
223 m_factorizationIsOk =
false;
224 m_symbolic = klu_analyze(internal::convert_index<int>(mp_matrix.rows()),
225 const_cast<StorageIndex *
>(mp_matrix.outerIndexPtr()),
226 const_cast<StorageIndex *
>(mp_matrix.innerIndexPtr()), &m_common);
228 m_isInitialized =
true;
230 m_analysisIsOk =
true;
231 m_extractedDataAreDirty =
true;
235 void factorize_impl() {
236 m_numeric = klu_factor(
const_cast<StorageIndex *
>(mp_matrix.outerIndexPtr()),
237 const_cast<StorageIndex *
>(mp_matrix.innerIndexPtr()),
238 const_cast<Scalar *
>(mp_matrix.valuePtr()), m_symbolic, &m_common, Scalar());
241 m_factorizationIsOk = m_numeric ? 1 : 0;
242 m_extractedDataAreDirty =
true;
245 template <
typename MatrixDerived>
246 void grab(
const EigenBase<MatrixDerived> &A) {
247 internal::destroy_at(&mp_matrix);
248 internal::construct_at(&mp_matrix, A.derived());
251 void grab(
const KLUMatrixRef &A) {
252 if (&(A.derived()) != &mp_matrix) {
253 internal::destroy_at(&mp_matrix);
254 internal::construct_at(&mp_matrix, A);
260 mutable LUMatrixType m_l;
261 mutable LUMatrixType m_u;
262 mutable IntColVectorType m_p;
263 mutable IntRowVectorType m_q;
266 KLUMatrixType m_dummy;
267 KLUMatrixRef mp_matrix;
269 klu_numeric *m_numeric;
270 klu_symbolic *m_symbolic;
272 mutable ComputationInfo m_info;
273 int m_factorizationIsOk;
275 mutable bool m_extractedDataAreDirty;
282template<
typename MatrixType>
283void KLU<MatrixType>::extractData()
const
285 if (m_extractedDataAreDirty)
287 eigen_assert(
false &&
"KLU: extractData Not Yet Implemented");
290 int lnz, unz, rows, cols, nz_udiag;
291 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
294 m_l.resize(rows,(std::min)(rows,cols));
295 m_l.resizeNonZeros(lnz);
297 m_u.resize((std::min)(rows,cols),cols);
298 m_u.resizeNonZeros(unz);
304 umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
305 m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
306 m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
308 m_extractedDataAreDirty =
false;
312template<
typename MatrixType>
313typename KLU<MatrixType>::Scalar KLU<MatrixType>::determinant()
const
315 eigen_assert(
false &&
"KLU: extractData Not Yet Implemented");
320template <
typename MatrixType>
321template <
typename BDerived,
typename XDerived>
322bool KLU<MatrixType>::_solve_impl(
const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x)
const {
323 Index rhsCols = b.cols();
324 EIGEN_STATIC_ASSERT((XDerived::Flags & RowMajorBit) == 0, THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
325 eigen_assert(m_factorizationIsOk &&
326 "The decomposition is not in a valid state for solving, you must first call either compute() or "
327 "analyzePattern()/factorize()");
330 int info =
klu_solve(m_symbolic, m_numeric, b.rows(), rhsCols, x.const_cast_derived().data(),
331 const_cast<klu_common *
>(&m_common), Scalar());
SparseSolverBase()
Definition SparseSolverBase.h:70
int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, double B[], klu_common *Common, double)
A sparse LU factorization and solver based on KLU.
Definition KLUSupport.h:36
@ NumericalIssue
Definition Constants.h:442
@ InvalidInput
Definition Constants.h:447
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition Core:137
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:83