10#ifndef EIGEN_KLUSUPPORT_H
11#define EIGEN_KLUSUPPORT_H
38inline int klu_solve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, std::complex<double>B[], klu_common *Common, std::complex<double>) {
39 return klu_z_solve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), &numext::real_ref(B[0]), Common);
42inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs,
double B[], klu_common *Common,
double) {
43 return klu_tsolve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), B, Common);
46inline int klu_tsolve(klu_symbolic *Symbolic, klu_numeric *Numeric, Index ldim, Index nrhs, std::complex<double>B[], klu_common *Common, std::complex<double>) {
47 return klu_z_tsolve(Symbolic, Numeric, internal::convert_index<int>(ldim), internal::convert_index<int>(nrhs), &numext::real_ref(B[0]), 0, Common);
50inline klu_numeric* klu_factor(
int Ap [ ],
int Ai [ ],
double Ax [ ], klu_symbolic *Symbolic, klu_common *Common,
double) {
51 return klu_factor(Ap, Ai, Ax, Symbolic, Common);
54inline klu_numeric* klu_factor(
int Ap[],
int Ai[], std::complex<double> Ax[], klu_symbolic *Symbolic, klu_common *Common, std::complex<double>) {
55 return klu_z_factor(Ap, Ai, &numext::real_ref(Ax[0]), Symbolic, Common);
59template<
typename _MatrixType>
64 using Base::m_isInitialized;
66 using Base::_solve_impl;
67 typedef _MatrixType MatrixType;
68 typedef typename MatrixType::Scalar Scalar;
69 typedef typename MatrixType::RealScalar RealScalar;
70 typedef typename MatrixType::StorageIndex StorageIndex;
78 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
79 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
85 : m_dummy(0,0), mp_matrix(m_dummy)
90 template<
typename InputMatrixType>
104 EIGEN_CONSTEXPR
inline Index rows()
const EIGEN_NOEXCEPT {
return mp_matrix.rows(); }
105 EIGEN_CONSTEXPR
inline Index cols()
const EIGEN_NOEXCEPT {
return mp_matrix.cols(); }
114 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
118 inline const LUMatrixType& matrixL()
const
120 if (m_extractedDataAreDirty) extractData();
124 inline const LUMatrixType& matrixU()
const
126 if (m_extractedDataAreDirty) extractData();
130 inline const IntColVectorType& permutationP()
const
132 if (m_extractedDataAreDirty) extractData();
136 inline const IntRowVectorType& permutationQ()
const
138 if (m_extractedDataAreDirty) extractData();
146 template<
typename InputMatrixType>
151 grab(matrix.derived());
152 analyzePattern_impl();
162 template<
typename InputMatrixType>
168 grab(matrix.derived());
170 analyzePattern_impl();
200 template<
typename InputMatrixType>
203 eigen_assert(m_analysisIsOk &&
"KLU: you must first call analyzePattern()");
207 grab(matrix.derived());
213 template<
typename BDerived,
typename XDerived>
217 Scalar determinant()
const;
219 void extractData()
const;
227 m_isInitialized =
false;
230 m_extractedDataAreDirty =
true;
235 void analyzePattern_impl()
238 m_analysisIsOk =
false;
239 m_factorizationIsOk =
false;
240 m_symbolic = klu_analyze(internal::convert_index<int>(mp_matrix.rows()),
241 const_cast<StorageIndex*
>(mp_matrix.outerIndexPtr()),
const_cast<StorageIndex*
>(mp_matrix.innerIndexPtr()),
244 m_isInitialized =
true;
246 m_analysisIsOk =
true;
247 m_extractedDataAreDirty =
true;
251 void factorize_impl()
254 m_numeric = klu_factor(
const_cast<StorageIndex*
>(mp_matrix.outerIndexPtr()),
const_cast<StorageIndex*
>(mp_matrix.innerIndexPtr()),
const_cast<Scalar*
>(mp_matrix.valuePtr()),
255 m_symbolic, &m_common, Scalar());
259 m_factorizationIsOk = m_numeric ? 1 : 0;
260 m_extractedDataAreDirty =
true;
263 template<
typename MatrixDerived>
264 void grab(
const EigenBase<MatrixDerived> &A)
266 mp_matrix.~KLUMatrixRef();
267 ::new (&mp_matrix) KLUMatrixRef(A.derived());
270 void grab(
const KLUMatrixRef &A)
272 if(&(A.derived()) != &mp_matrix)
274 mp_matrix.~KLUMatrixRef();
275 ::new (&mp_matrix) KLUMatrixRef(A);
281 mutable LUMatrixType m_l;
282 mutable LUMatrixType m_u;
283 mutable IntColVectorType m_p;
284 mutable IntRowVectorType m_q;
287 KLUMatrixType m_dummy;
288 KLUMatrixRef mp_matrix;
290 klu_numeric* m_numeric;
291 klu_symbolic* m_symbolic;
294 int m_factorizationIsOk;
296 mutable bool m_extractedDataAreDirty;
303template<
typename MatrixType>
304void KLU<MatrixType>::extractData()
const
306 if (m_extractedDataAreDirty)
308 eigen_assert(
false &&
"KLU: extractData Not Yet Implemented");
311 int lnz, unz, rows, cols, nz_udiag;
312 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
315 m_l.resize(rows,(std::min)(rows,cols));
316 m_l.resizeNonZeros(lnz);
318 m_u.resize((std::min)(rows,cols),cols);
319 m_u.resizeNonZeros(unz);
325 umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
326 m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
327 m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
329 m_extractedDataAreDirty =
false;
333template<
typename MatrixType>
334typename KLU<MatrixType>::Scalar KLU<MatrixType>::determinant()
const
336 eigen_assert(
false &&
"KLU: extractData Not Yet Implemented");
341template<
typename MatrixType>
342template<
typename BDerived,
typename XDerived>
343bool KLU<MatrixType>::_solve_impl(
const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x)
const
345 Index rhsCols = b.cols();
346 EIGEN_STATIC_ASSERT((XDerived::Flags&
RowMajorBit)==0, THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
347 eigen_assert(m_factorizationIsOk &&
"The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");
350 int info =
klu_solve(m_symbolic, m_numeric, b.rows(), rhsCols, x.const_cast_derived().data(),
const_cast<klu_common*
>(&m_common), Scalar());
Definition KLUSupport.h:61
void compute(const InputMatrixType &matrix)
Computes the sparse Cholesky decomposition of matrix Note that the matrix should be column-major,...
Definition KLUSupport.h:147
klu_common & kluCommon()
Provides access to the control settings array used by UmfPack.
Definition KLUSupport.h:189
void factorize(const InputMatrixType &matrix)
Performs a numeric decomposition of matrix.
Definition KLUSupport.h:201
void analyzePattern(const InputMatrixType &matrix)
Performs a symbolic decomposition on the sparcity of matrix.
Definition KLUSupport.h:163
const klu_common & kluCommon() const
Provides access to the control settings array used by KLU.
Definition KLUSupport.h:178
ComputationInfo info() const
Reports whether previous computation was successful.
Definition KLUSupport.h:112
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:50
A base class for sparse solvers.
Definition SparseSolverBase.h:68
ComputationInfo
Enum for reporting the status of a computation.
Definition Constants.h:440
@ NumericalIssue
The provided data did not satisfy the prerequisites.
Definition Constants.h:444
@ InvalidInput
The inputs are invalid, or the algorithm has been improperly called.
Definition Constants.h:449
@ Success
Computation was successful.
Definition Constants.h:442
const unsigned int RowMajorBit
for a matrix, this means that the storage order is row-major.
Definition Constants.h:66
Namespace containing all symbols from the Eigen library.
Definition LDLT.h:16
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
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:34