1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10 #ifndef EIGEN_UMFPACKSUPPORT_H
11 #define EIGEN_UMFPACKSUPPORT_H
12
13 namespace Eigen {
14
15 /* TODO extract L, extract U, compute det, etc... */
16
17 // generic double/complex<double> wrapper functions:
18
19
umfpack_defaults(double control[UMFPACK_CONTROL],double)20 inline void umfpack_defaults(double control[UMFPACK_CONTROL], double)
21 { umfpack_di_defaults(control); }
22
umfpack_defaults(double control[UMFPACK_CONTROL],std::complex<double>)23 inline void umfpack_defaults(double control[UMFPACK_CONTROL], std::complex<double>)
24 { umfpack_zi_defaults(control); }
25
umfpack_free_numeric(void ** Numeric,double)26 inline void umfpack_free_numeric(void **Numeric, double)
27 { umfpack_di_free_numeric(Numeric); *Numeric = 0; }
28
umfpack_free_numeric(void ** Numeric,std::complex<double>)29 inline void umfpack_free_numeric(void **Numeric, std::complex<double>)
30 { umfpack_zi_free_numeric(Numeric); *Numeric = 0; }
31
umfpack_free_symbolic(void ** Symbolic,double)32 inline void umfpack_free_symbolic(void **Symbolic, double)
33 { umfpack_di_free_symbolic(Symbolic); *Symbolic = 0; }
34
umfpack_free_symbolic(void ** Symbolic,std::complex<double>)35 inline void umfpack_free_symbolic(void **Symbolic, std::complex<double>)
36 { umfpack_zi_free_symbolic(Symbolic); *Symbolic = 0; }
37
umfpack_symbolic(int n_row,int n_col,const int Ap[],const int Ai[],const double Ax[],void ** Symbolic,const double Control[UMFPACK_CONTROL],double Info[UMFPACK_INFO])38 inline int umfpack_symbolic(int n_row,int n_col,
39 const int Ap[], const int Ai[], const double Ax[], void **Symbolic,
40 const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
41 {
42 return umfpack_di_symbolic(n_row,n_col,Ap,Ai,Ax,Symbolic,Control,Info);
43 }
44
umfpack_symbolic(int n_row,int n_col,const int Ap[],const int Ai[],const std::complex<double> Ax[],void ** Symbolic,const double Control[UMFPACK_CONTROL],double Info[UMFPACK_INFO])45 inline int umfpack_symbolic(int n_row,int n_col,
46 const int Ap[], const int Ai[], const std::complex<double> Ax[], void **Symbolic,
47 const double Control [UMFPACK_CONTROL], double Info [UMFPACK_INFO])
48 {
49 return umfpack_zi_symbolic(n_row,n_col,Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Control,Info);
50 }
51
umfpack_numeric(const int Ap[],const int Ai[],const double Ax[],void * Symbolic,void ** Numeric,const double Control[UMFPACK_CONTROL],double Info[UMFPACK_INFO])52 inline int umfpack_numeric( const int Ap[], const int Ai[], const double Ax[],
53 void *Symbolic, void **Numeric,
54 const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
55 {
56 return umfpack_di_numeric(Ap,Ai,Ax,Symbolic,Numeric,Control,Info);
57 }
58
umfpack_numeric(const int Ap[],const int Ai[],const std::complex<double> Ax[],void * Symbolic,void ** Numeric,const double Control[UMFPACK_CONTROL],double Info[UMFPACK_INFO])59 inline int umfpack_numeric( const int Ap[], const int Ai[], const std::complex<double> Ax[],
60 void *Symbolic, void **Numeric,
61 const double Control[UMFPACK_CONTROL],double Info [UMFPACK_INFO])
62 {
63 return umfpack_zi_numeric(Ap,Ai,&numext::real_ref(Ax[0]),0,Symbolic,Numeric,Control,Info);
64 }
65
umfpack_solve(int sys,const int Ap[],const int Ai[],const double Ax[],double X[],const double B[],void * Numeric,const double Control[UMFPACK_CONTROL],double Info[UMFPACK_INFO])66 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const double Ax[],
67 double X[], const double B[], void *Numeric,
68 const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
69 {
70 return umfpack_di_solve(sys,Ap,Ai,Ax,X,B,Numeric,Control,Info);
71 }
72
umfpack_solve(int sys,const int Ap[],const int Ai[],const std::complex<double> Ax[],std::complex<double> X[],const std::complex<double> B[],void * Numeric,const double Control[UMFPACK_CONTROL],double Info[UMFPACK_INFO])73 inline int umfpack_solve( int sys, const int Ap[], const int Ai[], const std::complex<double> Ax[],
74 std::complex<double> X[], const std::complex<double> B[], void *Numeric,
75 const double Control[UMFPACK_CONTROL], double Info[UMFPACK_INFO])
76 {
77 return umfpack_zi_solve(sys,Ap,Ai,&numext::real_ref(Ax[0]),0,&numext::real_ref(X[0]),0,&numext::real_ref(B[0]),0,Numeric,Control,Info);
78 }
79
umfpack_get_lunz(int * lnz,int * unz,int * n_row,int * n_col,int * nz_udiag,void * Numeric,double)80 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, double)
81 {
82 return umfpack_di_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
83 }
84
umfpack_get_lunz(int * lnz,int * unz,int * n_row,int * n_col,int * nz_udiag,void * Numeric,std::complex<double>)85 inline int umfpack_get_lunz(int *lnz, int *unz, int *n_row, int *n_col, int *nz_udiag, void *Numeric, std::complex<double>)
86 {
87 return umfpack_zi_get_lunz(lnz,unz,n_row,n_col,nz_udiag,Numeric);
88 }
89
umfpack_get_numeric(int Lp[],int Lj[],double Lx[],int Up[],int Ui[],double Ux[],int P[],int Q[],double Dx[],int * do_recip,double Rs[],void * Numeric)90 inline int umfpack_get_numeric(int Lp[], int Lj[], double Lx[], int Up[], int Ui[], double Ux[],
91 int P[], int Q[], double Dx[], int *do_recip, double Rs[], void *Numeric)
92 {
93 return umfpack_di_get_numeric(Lp,Lj,Lx,Up,Ui,Ux,P,Q,Dx,do_recip,Rs,Numeric);
94 }
95
umfpack_get_numeric(int Lp[],int Lj[],std::complex<double> Lx[],int Up[],int Ui[],std::complex<double> Ux[],int P[],int Q[],std::complex<double> Dx[],int * do_recip,double Rs[],void * Numeric)96 inline int umfpack_get_numeric(int Lp[], int Lj[], std::complex<double> Lx[], int Up[], int Ui[], std::complex<double> Ux[],
97 int P[], int Q[], std::complex<double> Dx[], int *do_recip, double Rs[], void *Numeric)
98 {
99 double& lx0_real = numext::real_ref(Lx[0]);
100 double& ux0_real = numext::real_ref(Ux[0]);
101 double& dx0_real = numext::real_ref(Dx[0]);
102 return umfpack_zi_get_numeric(Lp,Lj,Lx?&lx0_real:0,0,Up,Ui,Ux?&ux0_real:0,0,P,Q,
103 Dx?&dx0_real:0,0,do_recip,Rs,Numeric);
104 }
105
umfpack_get_determinant(double * Mx,double * Ex,void * NumericHandle,double User_Info[UMFPACK_INFO])106 inline int umfpack_get_determinant(double *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
107 {
108 return umfpack_di_get_determinant(Mx,Ex,NumericHandle,User_Info);
109 }
110
umfpack_get_determinant(std::complex<double> * Mx,double * Ex,void * NumericHandle,double User_Info[UMFPACK_INFO])111 inline int umfpack_get_determinant(std::complex<double> *Mx, double *Ex, void *NumericHandle, double User_Info [UMFPACK_INFO])
112 {
113 double& mx_real = numext::real_ref(*Mx);
114 return umfpack_zi_get_determinant(&mx_real,0,Ex,NumericHandle,User_Info);
115 }
116
117
118 /** \ingroup UmfPackSupport_Module
119 * \brief A sparse LU factorization and solver based on UmfPack
120 *
121 * This class allows to solve for A.X = B sparse linear problems via a LU factorization
122 * using the UmfPack library. The sparse matrix A must be squared and full rank.
123 * The vectors or matrices X and B can be either dense or sparse.
124 *
125 * \warning The input matrix A should be in a \b compressed and \b column-major form.
126 * Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
127 * \tparam _MatrixType the type of the sparse matrix A, it must be a SparseMatrix<>
128 *
129 * \implsparsesolverconcept
130 *
131 * \sa \ref TutorialSparseSolverConcept, class SparseLU
132 */
133 template<typename _MatrixType>
134 class UmfPackLU : public SparseSolverBase<UmfPackLU<_MatrixType> >
135 {
136 protected:
137 typedef SparseSolverBase<UmfPackLU<_MatrixType> > Base;
138 using Base::m_isInitialized;
139 public:
140 using Base::_solve_impl;
141 typedef _MatrixType MatrixType;
142 typedef typename MatrixType::Scalar Scalar;
143 typedef typename MatrixType::RealScalar RealScalar;
144 typedef typename MatrixType::StorageIndex StorageIndex;
145 typedef Matrix<Scalar,Dynamic,1> Vector;
146 typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
147 typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
148 typedef SparseMatrix<Scalar> LUMatrixType;
149 typedef SparseMatrix<Scalar,ColMajor,int> UmfpackMatrixType;
150 typedef Ref<const UmfpackMatrixType, StandardCompressedFormat> UmfpackMatrixRef;
151 enum {
152 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
153 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
154 };
155
156 public:
157
158 typedef Array<double, UMFPACK_CONTROL, 1> UmfpackControl;
159
UmfPackLU()160 UmfPackLU()
161 : m_dummy(0,0), mp_matrix(m_dummy)
162 {
163 init();
164 }
165
166 template<typename InputMatrixType>
UmfPackLU(const InputMatrixType & matrix)167 explicit UmfPackLU(const InputMatrixType& matrix)
168 : mp_matrix(matrix)
169 {
170 init();
171 compute(matrix);
172 }
173
~UmfPackLU()174 ~UmfPackLU()
175 {
176 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
177 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
178 }
179
rows()180 inline Index rows() const { return mp_matrix.rows(); }
cols()181 inline Index cols() const { return mp_matrix.cols(); }
182
183 /** \brief Reports whether previous computation was successful.
184 *
185 * \returns \c Success if computation was succesful,
186 * \c NumericalIssue if the matrix.appears to be negative.
187 */
info()188 ComputationInfo info() const
189 {
190 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
191 return m_info;
192 }
193
matrixL()194 inline const LUMatrixType& matrixL() const
195 {
196 if (m_extractedDataAreDirty) extractData();
197 return m_l;
198 }
199
matrixU()200 inline const LUMatrixType& matrixU() const
201 {
202 if (m_extractedDataAreDirty) extractData();
203 return m_u;
204 }
205
permutationP()206 inline const IntColVectorType& permutationP() const
207 {
208 if (m_extractedDataAreDirty) extractData();
209 return m_p;
210 }
211
permutationQ()212 inline const IntRowVectorType& permutationQ() const
213 {
214 if (m_extractedDataAreDirty) extractData();
215 return m_q;
216 }
217
218 /** Computes the sparse Cholesky decomposition of \a matrix
219 * Note that the matrix should be column-major, and in compressed format for best performance.
220 * \sa SparseMatrix::makeCompressed().
221 */
222 template<typename InputMatrixType>
compute(const InputMatrixType & matrix)223 void compute(const InputMatrixType& matrix)
224 {
225 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
226 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
227 grab(matrix.derived());
228 analyzePattern_impl();
229 factorize_impl();
230 }
231
232 /** Performs a symbolic decomposition on the sparcity of \a matrix.
233 *
234 * This function is particularly useful when solving for several problems having the same structure.
235 *
236 * \sa factorize(), compute()
237 */
238 template<typename InputMatrixType>
analyzePattern(const InputMatrixType & matrix)239 void analyzePattern(const InputMatrixType& matrix)
240 {
241 if(m_symbolic) umfpack_free_symbolic(&m_symbolic,Scalar());
242 if(m_numeric) umfpack_free_numeric(&m_numeric,Scalar());
243
244 grab(matrix.derived());
245
246 analyzePattern_impl();
247 }
248
249 /** Provides the return status code returned by UmfPack during the numeric
250 * factorization.
251 *
252 * \sa factorize(), compute()
253 */
umfpackFactorizeReturncode()254 inline int umfpackFactorizeReturncode() const
255 {
256 eigen_assert(m_numeric && "UmfPackLU: you must first call factorize()");
257 return m_fact_errorCode;
258 }
259
260 /** Provides access to the control settings array used by UmfPack.
261 *
262 * If this array contains NaN's, the default values are used.
263 *
264 * See UMFPACK documentation for details.
265 */
umfpackControl()266 inline const UmfpackControl& umfpackControl() const
267 {
268 return m_control;
269 }
270
271 /** Provides access to the control settings array used by UmfPack.
272 *
273 * If this array contains NaN's, the default values are used.
274 *
275 * See UMFPACK documentation for details.
276 */
umfpackControl()277 inline UmfpackControl& umfpackControl()
278 {
279 return m_control;
280 }
281
282 /** Performs a numeric decomposition of \a matrix
283 *
284 * The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.
285 *
286 * \sa analyzePattern(), compute()
287 */
288 template<typename InputMatrixType>
factorize(const InputMatrixType & matrix)289 void factorize(const InputMatrixType& matrix)
290 {
291 eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
292 if(m_numeric)
293 umfpack_free_numeric(&m_numeric,Scalar());
294
295 grab(matrix.derived());
296
297 factorize_impl();
298 }
299
300 /** \internal */
301 template<typename BDerived,typename XDerived>
302 bool _solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const;
303
304 Scalar determinant() const;
305
306 void extractData() const;
307
308 protected:
309
init()310 void init()
311 {
312 m_info = InvalidInput;
313 m_isInitialized = false;
314 m_numeric = 0;
315 m_symbolic = 0;
316 m_extractedDataAreDirty = true;
317 }
318
analyzePattern_impl()319 void analyzePattern_impl()
320 {
321 umfpack_defaults(m_control.data(), Scalar());
322 int errorCode = 0;
323 errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
324 internal::convert_index<int>(mp_matrix.cols()),
325 mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
326 &m_symbolic, m_control.data(), 0);
327
328 m_isInitialized = true;
329 m_info = errorCode ? InvalidInput : Success;
330 m_analysisIsOk = true;
331 m_factorizationIsOk = false;
332 m_extractedDataAreDirty = true;
333 }
334
factorize_impl()335 void factorize_impl()
336 {
337 m_fact_errorCode = umfpack_numeric(mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
338 m_symbolic, &m_numeric, m_control.data(), 0);
339
340 m_info = m_fact_errorCode == UMFPACK_OK ? Success : NumericalIssue;
341 m_factorizationIsOk = true;
342 m_extractedDataAreDirty = true;
343 }
344
345 template<typename MatrixDerived>
grab(const EigenBase<MatrixDerived> & A)346 void grab(const EigenBase<MatrixDerived> &A)
347 {
348 mp_matrix.~UmfpackMatrixRef();
349 ::new (&mp_matrix) UmfpackMatrixRef(A.derived());
350 }
351
grab(const UmfpackMatrixRef & A)352 void grab(const UmfpackMatrixRef &A)
353 {
354 if(&(A.derived()) != &mp_matrix)
355 {
356 mp_matrix.~UmfpackMatrixRef();
357 ::new (&mp_matrix) UmfpackMatrixRef(A);
358 }
359 }
360
361 // cached data to reduce reallocation, etc.
362 mutable LUMatrixType m_l;
363 int m_fact_errorCode;
364 UmfpackControl m_control;
365
366 mutable LUMatrixType m_u;
367 mutable IntColVectorType m_p;
368 mutable IntRowVectorType m_q;
369
370 UmfpackMatrixType m_dummy;
371 UmfpackMatrixRef mp_matrix;
372
373 void* m_numeric;
374 void* m_symbolic;
375
376 mutable ComputationInfo m_info;
377 int m_factorizationIsOk;
378 int m_analysisIsOk;
379 mutable bool m_extractedDataAreDirty;
380
381 private:
UmfPackLU(const UmfPackLU &)382 UmfPackLU(const UmfPackLU& ) { }
383 };
384
385
386 template<typename MatrixType>
extractData()387 void UmfPackLU<MatrixType>::extractData() const
388 {
389 if (m_extractedDataAreDirty)
390 {
391 // get size of the data
392 int lnz, unz, rows, cols, nz_udiag;
393 umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());
394
395 // allocate data
396 m_l.resize(rows,(std::min)(rows,cols));
397 m_l.resizeNonZeros(lnz);
398
399 m_u.resize((std::min)(rows,cols),cols);
400 m_u.resizeNonZeros(unz);
401
402 m_p.resize(rows);
403 m_q.resize(cols);
404
405 // extract
406 umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
407 m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
408 m_p.data(), m_q.data(), 0, 0, 0, m_numeric);
409
410 m_extractedDataAreDirty = false;
411 }
412 }
413
414 template<typename MatrixType>
determinant()415 typename UmfPackLU<MatrixType>::Scalar UmfPackLU<MatrixType>::determinant() const
416 {
417 Scalar det;
418 umfpack_get_determinant(&det, 0, m_numeric, 0);
419 return det;
420 }
421
422 template<typename MatrixType>
423 template<typename BDerived,typename XDerived>
_solve_impl(const MatrixBase<BDerived> & b,MatrixBase<XDerived> & x)424 bool UmfPackLU<MatrixType>::_solve_impl(const MatrixBase<BDerived> &b, MatrixBase<XDerived> &x) const
425 {
426 Index rhsCols = b.cols();
427 eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
428 eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
429 eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
430
431 int errorCode;
432 Scalar* x_ptr = 0;
433 Matrix<Scalar,Dynamic,1> x_tmp;
434 if(x.innerStride()!=1)
435 {
436 x_tmp.resize(x.rows());
437 x_ptr = x_tmp.data();
438 }
439 for (int j=0; j<rhsCols; ++j)
440 {
441 if(x.innerStride()==1)
442 x_ptr = &x.col(j).coeffRef(0);
443 errorCode = umfpack_solve(UMFPACK_A,
444 mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
445 x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, m_control.data(), 0);
446 if(x.innerStride()!=1)
447 x.col(j) = x_tmp;
448 if (errorCode!=0)
449 return false;
450 }
451
452 return true;
453 }
454
455 } // end namespace Eigen
456
457 #endif // EIGEN_UMFPACKSUPPORT_H
458