00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #ifndef EIGEN_SUPERLUSUPPORT_H
00026 #define EIGEN_SUPERLUSUPPORT_H
00027
00028
00029 #define DECL_GSSVX(NAMESPACE,FNAME,FLOATTYPE,KEYTYPE) \
00030 inline float SuperLU_gssvx(superlu_options_t *options, SuperMatrix *A, \
00031 int *perm_c, int *perm_r, int *etree, char *equed, \
00032 FLOATTYPE *R, FLOATTYPE *C, SuperMatrix *L, \
00033 SuperMatrix *U, void *work, int lwork, \
00034 SuperMatrix *B, SuperMatrix *X, \
00035 FLOATTYPE *recip_pivot_growth, \
00036 FLOATTYPE *rcond, FLOATTYPE *ferr, FLOATTYPE *berr, \
00037 SuperLUStat_t *stats, int *info, KEYTYPE) { \
00038 using namespace NAMESPACE; \
00039 mem_usage_t mem_usage; \
00040 NAMESPACE::FNAME(options, A, perm_c, perm_r, etree, equed, R, C, L, \
00041 U, work, lwork, B, X, recip_pivot_growth, rcond, \
00042 ferr, berr, &mem_usage, stats, info); \
00043 return mem_usage.for_lu; \
00044 }
00045
00046 DECL_GSSVX(SuperLU_S,sgssvx,float,float)
00047 DECL_GSSVX(SuperLU_C,cgssvx,float,std::complex<float>)
00048 DECL_GSSVX(SuperLU_D,dgssvx,double,double)
00049 DECL_GSSVX(SuperLU_Z,zgssvx,double,std::complex<double>)
00050
00051 template<typename MatrixType>
00052 struct SluMatrixMapHelper;
00053
00061 struct SluMatrix : SuperMatrix
00062 {
00063 SluMatrix()
00064 {
00065 Store = &storage;
00066 }
00067
00068 SluMatrix(const SluMatrix& other)
00069 : SuperMatrix(other)
00070 {
00071 Store = &storage;
00072 storage = other.storage;
00073 }
00074
00075 SluMatrix& operator=(const SluMatrix& other)
00076 {
00077 SuperMatrix::operator=(static_cast<const SuperMatrix&>(other));
00078 Store = &storage;
00079 storage = other.storage;
00080 return *this;
00081 }
00082
00083 struct
00084 {
00085 union {int nnz;int lda;};
00086 void *values;
00087 int *innerInd;
00088 int *outerInd;
00089 } storage;
00090
00091 void setStorageType(Stype_t t)
00092 {
00093 Stype = t;
00094 if (t==SLU_NC || t==SLU_NR || t==SLU_DN)
00095 Store = &storage;
00096 else
00097 {
00098 ei_assert(false && "storage type not supported");
00099 Store = 0;
00100 }
00101 }
00102
00103 template<typename Scalar>
00104 void setScalarType()
00105 {
00106 if (ei_is_same_type<Scalar,float>::ret)
00107 Dtype = SLU_S;
00108 else if (ei_is_same_type<Scalar,double>::ret)
00109 Dtype = SLU_D;
00110 else if (ei_is_same_type<Scalar,std::complex<float> >::ret)
00111 Dtype = SLU_C;
00112 else if (ei_is_same_type<Scalar,std::complex<double> >::ret)
00113 Dtype = SLU_Z;
00114 else
00115 {
00116 ei_assert(false && "Scalar type not supported by SuperLU");
00117 }
00118 }
00119
00120 template<typename Scalar, int Rows, int Cols, int Options, int MRows, int MCols>
00121 static SluMatrix Map(Matrix<Scalar,Rows,Cols,Options,MRows,MCols>& mat)
00122 {
00123 typedef Matrix<Scalar,Rows,Cols,Options,MRows,MCols> MatrixType;
00124 ei_assert( ((Options&RowMajor)!=RowMajor) && "row-major dense matrices is not supported by SuperLU");
00125 SluMatrix res;
00126 res.setStorageType(SLU_DN);
00127 res.setScalarType<Scalar>();
00128 res.Mtype = SLU_GE;
00129
00130 res.nrow = mat.rows();
00131 res.ncol = mat.cols();
00132
00133 res.storage.lda = mat.stride();
00134 res.storage.values = mat.data();
00135 return res;
00136 }
00137
00138 template<typename MatrixType>
00139 static SluMatrix Map(SparseMatrixBase<MatrixType>& mat)
00140 {
00141 SluMatrix res;
00142 if ((MatrixType::Flags&RowMajorBit)==RowMajorBit)
00143 {
00144 res.setStorageType(SLU_NR);
00145 res.nrow = mat.cols();
00146 res.ncol = mat.rows();
00147 }
00148 else
00149 {
00150 res.setStorageType(SLU_NC);
00151 res.nrow = mat.rows();
00152 res.ncol = mat.cols();
00153 }
00154
00155 res.Mtype = SLU_GE;
00156
00157 res.storage.nnz = mat.nonZeros();
00158 res.storage.values = mat.derived()._valuePtr();
00159 res.storage.innerInd = mat.derived()._innerIndexPtr();
00160 res.storage.outerInd = mat.derived()._outerIndexPtr();
00161
00162 res.setScalarType<typename MatrixType::Scalar>();
00163
00164
00165 if (MatrixType::Flags & UpperTriangular)
00166 res.Mtype = SLU_TRU;
00167 if (MatrixType::Flags & LowerTriangular)
00168 res.Mtype = SLU_TRL;
00169 if (MatrixType::Flags & SelfAdjoint)
00170 ei_assert(false && "SelfAdjoint matrix shape not supported by SuperLU");
00171 return res;
00172 }
00173 };
00174
00175 template<typename Scalar, int Rows, int Cols, int Options, int MRows, int MCols>
00176 struct SluMatrixMapHelper<Matrix<Scalar,Rows,Cols,Options,MRows,MCols> >
00177 {
00178 typedef Matrix<Scalar,Rows,Cols,Options,MRows,MCols> MatrixType;
00179 static void run(MatrixType& mat, SluMatrix& res)
00180 {
00181 ei_assert( ((Options&RowMajor)!=RowMajor) && "row-major dense matrices is not supported by SuperLU");
00182 res.setStorageType(SLU_DN);
00183 res.setScalarType<Scalar>();
00184 res.Mtype = SLU_GE;
00185
00186 res.nrow = mat.rows();
00187 res.ncol = mat.cols();
00188
00189 res.storage.lda = mat.stride();
00190 res.storage.values = mat.data();
00191 }
00192 };
00193
00194 template<typename Derived>
00195 struct SluMatrixMapHelper<SparseMatrixBase<Derived> >
00196 {
00197 typedef Derived MatrixType;
00198 static void run(MatrixType& mat, SluMatrix& res)
00199 {
00200 if ((MatrixType::Flags&RowMajorBit)==RowMajorBit)
00201 {
00202 res.setStorageType(SLU_NR);
00203 res.nrow = mat.cols();
00204 res.ncol = mat.rows();
00205 }
00206 else
00207 {
00208 res.setStorageType(SLU_NC);
00209 res.nrow = mat.rows();
00210 res.ncol = mat.cols();
00211 }
00212
00213 res.Mtype = SLU_GE;
00214
00215 res.storage.nnz = mat.nonZeros();
00216 res.storage.values = mat._valuePtr();
00217 res.storage.innerInd = mat._innerIndexPtr();
00218 res.storage.outerInd = mat._outerIndexPtr();
00219
00220 res.setScalarType<typename MatrixType::Scalar>();
00221
00222
00223 if (MatrixType::Flags & UpperTriangular)
00224 res.Mtype = SLU_TRU;
00225 if (MatrixType::Flags & LowerTriangular)
00226 res.Mtype = SLU_TRL;
00227 if (MatrixType::Flags & SelfAdjoint)
00228 ei_assert(false && "SelfAdjoint matrix shape not supported by SuperLU");
00229 }
00230 };
00231
00232 template<typename Derived>
00233 SluMatrix SparseMatrixBase<Derived>::asSluMatrix()
00234 {
00235 return SluMatrix::Map(derived());
00236 }
00237
00239 template<typename Scalar, int Flags>
00240 MappedSparseMatrix<Scalar,Flags>::MappedSparseMatrix(SluMatrix& sluMat)
00241 {
00242 if ((Flags&RowMajorBit)==RowMajorBit)
00243 {
00244 assert(sluMat.Stype == SLU_NR);
00245 m_innerSize = sluMat.ncol;
00246 m_outerSize = sluMat.nrow;
00247 }
00248 else
00249 {
00250 assert(sluMat.Stype == SLU_NC);
00251 m_innerSize = sluMat.nrow;
00252 m_outerSize = sluMat.ncol;
00253 }
00254 m_outerIndex = sluMat.storage.outerInd;
00255 m_innerIndices = sluMat.storage.innerInd;
00256 m_values = reinterpret_cast<Scalar*>(sluMat.storage.values);
00257 m_nnz = sluMat.storage.outerInd[m_outerSize];
00258 }
00259
00260 template<typename MatrixType>
00261 class SparseLU<MatrixType,SuperLU> : public SparseLU<MatrixType>
00262 {
00263 protected:
00264 typedef SparseLU<MatrixType> Base;
00265 typedef typename Base::Scalar Scalar;
00266 typedef typename Base::RealScalar RealScalar;
00267 typedef Matrix<Scalar,Dynamic,1> Vector;
00268 typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> IntRowVectorType;
00269 typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> IntColVectorType;
00270 typedef SparseMatrix<Scalar,LowerTriangular|UnitDiagBit> LMatrixType;
00271 typedef SparseMatrix<Scalar,UpperTriangular> UMatrixType;
00272 using Base::m_flags;
00273 using Base::m_status;
00274
00275 public:
00276
00277 SparseLU(int flags = NaturalOrdering)
00278 : Base(flags)
00279 {
00280 }
00281
00282 SparseLU(const MatrixType& matrix, int flags = NaturalOrdering)
00283 : Base(flags)
00284 {
00285 compute(matrix);
00286 }
00287
00288 ~SparseLU()
00289 {
00290 }
00291
00292 inline const LMatrixType& matrixL() const
00293 {
00294 if (m_extractedDataAreDirty) extractData();
00295 return m_l;
00296 }
00297
00298 inline const UMatrixType& matrixU() const
00299 {
00300 if (m_extractedDataAreDirty) extractData();
00301 return m_u;
00302 }
00303
00304 inline const IntColVectorType& permutationP() const
00305 {
00306 if (m_extractedDataAreDirty) extractData();
00307 return m_p;
00308 }
00309
00310 inline const IntRowVectorType& permutationQ() const
00311 {
00312 if (m_extractedDataAreDirty) extractData();
00313 return m_q;
00314 }
00315
00316 Scalar determinant() const;
00317
00318 template<typename BDerived, typename XDerived>
00319 bool solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived>* x) const;
00320
00321 void compute(const MatrixType& matrix);
00322
00323 protected:
00324
00325 void extractData() const;
00326
00327 protected:
00328
00329 mutable LMatrixType m_l;
00330 mutable UMatrixType m_u;
00331 mutable IntColVectorType m_p;
00332 mutable IntRowVectorType m_q;
00333
00334 mutable SparseMatrix<Scalar> m_matrix;
00335 mutable SluMatrix m_sluA;
00336 mutable SuperMatrix m_sluL, m_sluU;
00337 mutable SluMatrix m_sluB, m_sluX;
00338 mutable SuperLUStat_t m_sluStat;
00339 mutable superlu_options_t m_sluOptions;
00340 mutable std::vector<int> m_sluEtree;
00341 mutable std::vector<RealScalar> m_sluRscale, m_sluCscale;
00342 mutable std::vector<RealScalar> m_sluFerr, m_sluBerr;
00343 mutable char m_sluEqued;
00344 mutable bool m_extractedDataAreDirty;
00345 };
00346
00347 template<typename MatrixType>
00348 void SparseLU<MatrixType,SuperLU>::compute(const MatrixType& a)
00349 {
00350 const int size = a.rows();
00351 m_matrix = a;
00352
00353 set_default_options(&m_sluOptions);
00354 m_sluOptions.ColPerm = NATURAL;
00355 m_sluOptions.PrintStat = NO;
00356 m_sluOptions.ConditionNumber = NO;
00357 m_sluOptions.Trans = NOTRANS;
00358
00359
00360 switch (Base::orderingMethod())
00361 {
00362 case NaturalOrdering : m_sluOptions.ColPerm = NATURAL; break;
00363 case MinimumDegree_AT_PLUS_A : m_sluOptions.ColPerm = MMD_AT_PLUS_A; break;
00364 case MinimumDegree_ATA : m_sluOptions.ColPerm = MMD_ATA; break;
00365 case ColApproxMinimumDegree : m_sluOptions.ColPerm = COLAMD; break;
00366 default:
00367 std::cerr << "Eigen: ordering method \"" << Base::orderingMethod() << "\" not supported by the SuperLU backend\n";
00368 m_sluOptions.ColPerm = NATURAL;
00369 };
00370
00371 m_sluA = m_matrix.asSluMatrix();
00372 memset(&m_sluL,0,sizeof m_sluL);
00373 memset(&m_sluU,0,sizeof m_sluU);
00374 m_sluEqued = 'B';
00375 int info = 0;
00376
00377 m_p.resize(size);
00378 m_q.resize(size);
00379 m_sluRscale.resize(size);
00380 m_sluCscale.resize(size);
00381 m_sluEtree.resize(size);
00382
00383 RealScalar recip_pivot_gross, rcond;
00384 RealScalar ferr, berr;
00385
00386
00387 m_sluB.setStorageType(SLU_DN);
00388 m_sluB.setScalarType<Scalar>();
00389 m_sluB.Mtype = SLU_GE;
00390 m_sluB.storage.values = 0;
00391 m_sluB.nrow = m_sluB.ncol = 0;
00392 m_sluB.storage.lda = size;
00393 m_sluX = m_sluB;
00394
00395 StatInit(&m_sluStat);
00396 SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
00397 &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0],
00398 &m_sluL, &m_sluU,
00399 NULL, 0,
00400 &m_sluB, &m_sluX,
00401 &recip_pivot_gross, &rcond,
00402 &ferr, &berr,
00403 &m_sluStat, &info, Scalar());
00404 StatFree(&m_sluStat);
00405
00406 m_extractedDataAreDirty = true;
00407
00408
00409 Base::m_succeeded = (info == 0);
00410 }
00411
00412 template<typename MatrixType>
00413 template<typename BDerived,typename XDerived>
00414 bool SparseLU<MatrixType,SuperLU>::solve(const MatrixBase<BDerived> &b, MatrixBase<XDerived> *x) const
00415 {
00416 const int size = m_matrix.rows();
00417 const int rhsCols = b.cols();
00418 ei_assert(size==b.rows());
00419
00420 m_sluOptions.Fact = FACTORED;
00421 m_sluOptions.IterRefine = NOREFINE;
00422
00423 m_sluFerr.resize(rhsCols);
00424 m_sluBerr.resize(rhsCols);
00425 m_sluB = SluMatrix::Map(b.const_cast_derived());
00426 m_sluX = SluMatrix::Map(x->derived());
00427
00428 StatInit(&m_sluStat);
00429 int info = 0;
00430 RealScalar recip_pivot_gross, rcond;
00431 SuperLU_gssvx(
00432 &m_sluOptions, &m_sluA,
00433 m_q.data(), m_p.data(),
00434 &m_sluEtree[0], &m_sluEqued,
00435 &m_sluRscale[0], &m_sluCscale[0],
00436 &m_sluL, &m_sluU,
00437 NULL, 0,
00438 &m_sluB, &m_sluX,
00439 &recip_pivot_gross, &rcond,
00440 &m_sluFerr[0], &m_sluBerr[0],
00441 &m_sluStat, &info, Scalar());
00442 StatFree(&m_sluStat);
00443
00444 return info==0;
00445 }
00446
00447
00448
00449
00450
00451
00452
00453
00454
00455 template<typename MatrixType>
00456 void SparseLU<MatrixType,SuperLU>::extractData() const
00457 {
00458 if (m_extractedDataAreDirty)
00459 {
00460 int upper;
00461 int fsupc, istart, nsupr;
00462 int lastl = 0, lastu = 0;
00463 SCformat *Lstore = static_cast<SCformat*>(m_sluL.Store);
00464 NCformat *Ustore = static_cast<NCformat*>(m_sluU.Store);
00465 Scalar *SNptr;
00466
00467 const int size = m_matrix.rows();
00468 m_l.resize(size,size);
00469 m_l.resizeNonZeros(Lstore->nnz);
00470 m_u.resize(size,size);
00471 m_u.resizeNonZeros(Ustore->nnz);
00472
00473 int* Lcol = m_l._outerIndexPtr();
00474 int* Lrow = m_l._innerIndexPtr();
00475 Scalar* Lval = m_l._valuePtr();
00476
00477 int* Ucol = m_u._outerIndexPtr();
00478 int* Urow = m_u._innerIndexPtr();
00479 Scalar* Uval = m_u._valuePtr();
00480
00481 Ucol[0] = 0;
00482 Ucol[0] = 0;
00483
00484
00485 for (int k = 0; k <= Lstore->nsuper; ++k)
00486 {
00487 fsupc = L_FST_SUPC(k);
00488 istart = L_SUB_START(fsupc);
00489 nsupr = L_SUB_START(fsupc+1) - istart;
00490 upper = 1;
00491
00492
00493 for (int j = fsupc; j < L_FST_SUPC(k+1); ++j)
00494 {
00495 SNptr = &((Scalar*)Lstore->nzval)[L_NZ_START(j)];
00496
00497
00498 for (int i = U_NZ_START(j); i < U_NZ_START(j+1); ++i)
00499 {
00500 Uval[lastu] = ((Scalar*)Ustore->nzval)[i];
00501
00502 if (Uval[lastu] != 0.0)
00503 Urow[lastu++] = U_SUB(i);
00504 }
00505 for (int i = 0; i < upper; ++i)
00506 {
00507
00508 Uval[lastu] = SNptr[i];
00509
00510 if (Uval[lastu] != 0.0)
00511 Urow[lastu++] = L_SUB(istart+i);
00512 }
00513 Ucol[j+1] = lastu;
00514
00515
00516 Lval[lastl] = 1.0;
00517 Lrow[lastl++] = L_SUB(istart + upper - 1);
00518 for (int i = upper; i < nsupr; ++i)
00519 {
00520 Lval[lastl] = SNptr[i];
00521
00522 if (Lval[lastl] != 0.0)
00523 Lrow[lastl++] = L_SUB(istart+i);
00524 }
00525 Lcol[j+1] = lastl;
00526
00527 ++upper;
00528 }
00529
00530 }
00531
00532
00533 m_l.resizeNonZeros(lastl);
00534 m_u.resizeNonZeros(lastu);
00535
00536 m_extractedDataAreDirty = false;
00537 }
00538 }
00539
00540 template<typename MatrixType>
00541 typename SparseLU<MatrixType,SuperLU>::Scalar SparseLU<MatrixType,SuperLU>::determinant() const
00542 {
00543 if (m_extractedDataAreDirty)
00544 extractData();
00545
00546
00547
00548 Scalar det = Scalar(1);
00549 for (int j=0; j<m_u.cols(); ++j)
00550 {
00551 if (m_u._outerIndexPtr()[j+1]-m_u._outerIndexPtr()[j] > 0)
00552 {
00553 int lastId = m_u._outerIndexPtr()[j+1]-1;
00554 ei_assert(m_u._innerIndexPtr()[lastId]<=j);
00555 if (m_u._innerIndexPtr()[lastId]==j)
00556 {
00557 det *= m_u._valuePtr()[lastId];
00558 }
00559 }
00560
00561 }
00562 return det;
00563 }
00564
00565 #endif // EIGEN_SUPERLUSUPPORT_H