00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020 #ifndef SELDON_FILE_SUPERLU_CXX
00021
00022 #include "SuperLU.hxx"
00023
00024
00025
00026
00027
00028 void LUextract(SuperMatrix *L, SuperMatrix *U, double *Lval, int *Lrow,
00029 int *Lcol, double *Uval, int *Urow, int *Ucol, int *snnzL,
00030 int *snnzU)
00031 {
00032 int i, j, k;
00033 int upper;
00034 int fsupc, istart, nsupr;
00035 int lastl = 0, lastu = 0;
00036 SCformat *Lstore;
00037 NCformat *Ustore;
00038 double *SNptr;
00039
00040 Lstore = static_cast<SCformat*>(L->Store);
00041 Ustore = static_cast<NCformat*>(U->Store);
00042 Lcol[0] = 0;
00043 Ucol[0] = 0;
00044
00045
00046 for (k = 0; k <= Lstore->nsuper; ++k) {
00047
00048 fsupc = L_FST_SUPC(k);
00049 istart = L_SUB_START(fsupc);
00050 nsupr = L_SUB_START(fsupc+1) - istart;
00051 upper = 1;
00052
00053
00054 for (j = fsupc; j < L_FST_SUPC(k+1); ++j) {
00055 SNptr = &(static_cast<double*>(Lstore->nzval))[L_NZ_START(j)];
00056
00057
00058 for (i = U_NZ_START(j); i < U_NZ_START(j+1); ++i) {
00059 Uval[lastu] = (static_cast<double*>(Ustore->nzval))[i];
00060 Urow[lastu++] = U_SUB(i);
00061 }
00062 for (i = 0; i < upper; ++i) {
00063 Uval[lastu] = SNptr[i];
00064 Urow[lastu++] = L_SUB(istart+i);
00065 }
00066 Ucol[j+1] = lastu;
00067
00068
00069 Lval[lastl] = 1.0;
00070 Lrow[lastl++] = L_SUB(istart + upper - 1);
00071 for (i = upper; i < nsupr; ++i) {
00072 Lval[lastl] = SNptr[i];
00073 Lrow[lastl++] = L_SUB(istart+i);
00074 }
00075 Lcol[j+1] = lastl;
00076
00077 ++upper;
00078
00079 }
00080
00081 }
00082
00083 *snnzL = lastl;
00084 *snnzU = lastu;
00085 }
00086
00087
00088 namespace Seldon
00089 {
00091 template<class T>
00092 MatrixSuperLU_Base<T>::MatrixSuperLU_Base()
00093 {
00094
00095
00096
00097
00098 n = 0;
00099 permc_spec = 2;
00100 Lstore = NULL;
00101 Ustore = NULL;
00102 StatInit(&stat);
00103 set_default_options(&options);
00104 ShowMessages();
00105 display_info = false;
00106 }
00107
00108
00110 template<class T>
00111 MatrixSuperLU_Base<T>::~MatrixSuperLU_Base()
00112 {
00113 Clear();
00114 }
00115
00116
00118
00127 template<class T>
00128 template<class Prop, class Allocator>
00129 void MatrixSuperLU_Base<T>
00130 ::GetLU(Matrix<double, Prop, ColSparse, Allocator>& Lmat,
00131 Matrix<double, Prop, ColSparse, Allocator>& Umat,
00132 bool permuted)
00133 {
00134 Lstore = static_cast<SCformat*>(L.Store);
00135 Ustore = static_cast<NCformat*>(U.Store);
00136
00137 int Lnnz = Lstore->nnz;
00138 int Unnz = Ustore->nnz;
00139
00140 int m = U.nrow;
00141 int n = U.ncol;
00142
00143 Vector<double, VectFull, Allocator> Lval(Lnnz);
00144 Vector<int, VectFull, CallocAlloc<int> > Lrow(Lnnz);
00145 Vector<int, VectFull, CallocAlloc<int> > Lcol(n + 1);
00146
00147 Vector<double, VectFull, Allocator> Uval(Unnz);
00148 Vector<int, VectFull, CallocAlloc<int> > Urow(Unnz);
00149 Vector<int, VectFull, CallocAlloc<int> > Ucol(n + 1);
00150
00151 int Lsnnz;
00152 int Usnnz;
00153 LUextract(&L, &U, Lval.GetData(), Lrow.GetData(), Lcol.GetData(),
00154 Uval.GetData(), Urow.GetData(), Ucol.GetData(), &Lsnnz, &Usnnz);
00155
00156 Lmat.SetData(m, n, Lval, Lcol, Lrow);
00157 Umat.SetData(m, n, Uval, Ucol, Urow);
00158
00159 if (!permuted)
00160 {
00161 Vector<int> row_perm_orig = perm_r;
00162 Vector<int> col_perm_orig = perm_c;
00163
00164 Vector<int> row_perm(n);
00165 Vector<int> col_perm(n);
00166 row_perm.Fill();
00167 col_perm.Fill();
00168
00169 Sort(row_perm_orig, row_perm);
00170 Sort(col_perm_orig, col_perm);
00171
00172 ApplyInversePermutation(Lmat, row_perm, col_perm);
00173 ApplyInversePermutation(Umat, row_perm, col_perm);
00174 }
00175 }
00176
00177
00179
00190 template<class T>
00191 template<class Prop, class Allocator>
00192 void MatrixSuperLU_Base<T>
00193 ::GetLU(Matrix<double, Prop, RowSparse, Allocator>& Lmat,
00194 Matrix<double, Prop, RowSparse, Allocator>& Umat,
00195 bool permuted)
00196 {
00197 Lmat.Clear();
00198 Umat.Clear();
00199
00200 Matrix<double, Prop, ColSparse, Allocator> Lmat_col;
00201 Matrix<double, Prop, ColSparse, Allocator> Umat_col;
00202 GetLU(Lmat_col, Umat_col, permuted);
00203
00204 Copy(Lmat_col, Lmat);
00205 Lmat_col.Clear();
00206 Copy(Umat_col, Umat);
00207 Umat_col.Clear();
00208 }
00209
00210
00212
00218 template<class T>
00219 const Vector<int>& MatrixSuperLU_Base<T>::GetRowPermutation() const
00220 {
00221 return perm_r;
00222 }
00223
00224
00226
00232 template<class T>
00233 const Vector<int>& MatrixSuperLU_Base<T>::GetColPermutation() const
00234 {
00235 return perm_c;
00236 }
00237
00238
00240 template<class T>
00241 void MatrixSuperLU_Base<T>::Clear()
00242 {
00243 if (n > 0)
00244 {
00245
00246 Destroy_CompCol_Matrix(&A);
00247 Destroy_SuperMatrix_Store(&B);
00248 Destroy_SuperNode_Matrix(&L);
00249 Destroy_CompCol_Matrix(&U);
00250 perm_r.Clear(); perm_c.Clear();
00251 n = 0;
00252 }
00253 }
00254
00255
00257 template<class T>
00258 void MatrixSuperLU_Base<T>::HideMessages()
00259 {
00260 display_info = false;
00261 }
00262
00263
00265 template<class T>
00266 void MatrixSuperLU_Base<T>::ShowMessages()
00267 {
00268 display_info = true;
00269 }
00270
00271
00273 template<class Prop, class Storage, class Allocator>
00274 void MatrixSuperLU<double>::
00275 FactorizeMatrix(Matrix<double, Prop, Storage, Allocator> & mat,
00276 bool keep_matrix)
00277 {
00278
00279 Clear();
00280
00281
00282 n = mat.GetN();
00283 Matrix<double, General, ColSparse> Acsr;
00284 Copy(mat, Acsr);
00285 if (!keep_matrix)
00286 mat.Clear();
00287
00288
00289 int nnz = Acsr.GetDataSize();
00290 dCreate_CompCol_Matrix(&A, n, n, nnz, Acsr.GetData(), Acsr.GetInd(),
00291 Acsr.GetPtr(), SLU_NC, SLU_D, SLU_GE);
00292
00293 perm_r.Reallocate(n);
00294 perm_c.Reallocate(n);
00295
00296
00297 int nb_rhs = 0, info;
00298 dCreate_Dense_Matrix(&B, n, nb_rhs, NULL, n, SLU_DN, SLU_D, SLU_GE);
00299
00300 dgssv(&options, &A, perm_c.GetData(), perm_r.GetData(),
00301 &L, &U, &B, &stat, &info);
00302
00303 if ((info==0)&&(display_info))
00304 {
00305 mem_usage_t mem_usage;
00306 Lstore = (SCformat *) L.Store;
00307 Ustore = (NCformat *) U.Store;
00308 cout<<"No of nonzeros in factor L = "<<Lstore->nnz<<endl;
00309 cout<<"No of nonzeros in factor U = "<<Ustore->nnz<<endl;
00310 cout<<"No of nonzeros in L+U = "<<(Lstore->nnz+Ustore->nnz)<<endl;
00311 dQuerySpace(&L, &U, &mem_usage);
00312 cout<<"Memory used for factorisation in Mo "
00313 <<mem_usage.total_needed/(1024*1024)<<endl;
00314 }
00315
00316 Acsr.Nullify();
00317 }
00318
00319
00321 template<class Allocator2>
00322 void MatrixSuperLU<double>::Solve(Vector<double, VectFull, Allocator2>& x)
00323 {
00324 trans_t trans = NOTRANS;
00325 int nb_rhs = 1, info;
00326 dCreate_Dense_Matrix(&B, x.GetM(), nb_rhs,
00327 x.GetData(), x.GetM(), SLU_DN, SLU_D, SLU_GE);
00328
00329 SuperLUStat_t stat;
00330 StatInit(&stat);
00331 dgstrs(trans, &L, &U, perm_r.GetData(),
00332 perm_c.GetData(), &B, &stat, &info);
00333 }
00334
00335
00337 template<class Allocator2>
00338 void MatrixSuperLU<double>::Solve(const SeldonTranspose& TransA,
00339 Vector<double, VectFull, Allocator2>& x)
00340 {
00341 if (TransA.NoTrans())
00342 {
00343 Solve(x);
00344 return;
00345 }
00346
00347 trans_t trans = TRANS;
00348 int nb_rhs = 1, info;
00349 dCreate_Dense_Matrix(&B, x.GetM(), nb_rhs,
00350 x.GetData(), x.GetM(), SLU_DN, SLU_D, SLU_GE);
00351
00352 SuperLUStat_t stat;
00353 StatInit(&stat);
00354 dgstrs(trans, &L, &U, perm_r.GetData(),
00355 perm_c.GetData(), &B, &stat, &info);
00356 }
00357
00358
00360 template<class Prop, class Storage, class Allocator>
00361 void MatrixSuperLU<complex<double> >::
00362 FactorizeMatrix(Matrix<complex<double>, Prop, Storage, Allocator> & mat,
00363 bool keep_matrix)
00364 {
00365
00366 Clear();
00367
00368
00369 n = mat.GetN();
00370 Matrix<complex<double>, General, ColSparse> Acsr;
00371 Copy(mat, Acsr);
00372 if (!keep_matrix)
00373 mat.Clear();
00374
00375
00376 int nnz = Acsr.GetDataSize();
00377 zCreate_CompCol_Matrix(&A, n, n, nnz,
00378 reinterpret_cast<doublecomplex*>(Acsr.GetData()),
00379 Acsr.GetInd(), Acsr.GetPtr(),
00380 SLU_NC, SLU_Z, SLU_GE);
00381
00382 perm_r.Reallocate(n);
00383 perm_c.Reallocate(n);
00384
00385 int nb_rhs = 0, info;
00386 zCreate_Dense_Matrix(&B, n, nb_rhs, NULL, n, SLU_DN, SLU_Z, SLU_GE);
00387
00388 zgssv(&options, &A, perm_c.GetData(), perm_r.GetData(),
00389 &L, &U, &B, &stat, &info);
00390
00391 if ((info==0)&&(display_info))
00392 {
00393 mem_usage_t mem_usage;
00394 Lstore = (SCformat *) L.Store;
00395 Ustore = (NCformat *) U.Store;
00396 cout<<"No of nonzeros in factor L = "<<Lstore->nnz<<endl;
00397 cout<<"No of nonzeros in factor U = "<<Ustore->nnz<<endl;
00398 cout<<"No of nonzeros in L+U = "<<(Lstore->nnz+Ustore->nnz)<<endl;
00399 zQuerySpace(&L, &U, &mem_usage);
00400 cout<<"Memory used for factorisation in Mo "
00401 <<mem_usage.total_needed/1e6<<endl;
00402 }
00403
00404 Acsr.Nullify();
00405 }
00406
00407
00409 template<class Allocator2>
00410 void MatrixSuperLU<complex<double> >::
00411 Solve(Vector<complex<double>, VectFull, Allocator2>& x)
00412 {
00413 trans_t trans = NOTRANS;
00414 int nb_rhs = 1, info;
00415 zCreate_Dense_Matrix(&B, x.GetM(), nb_rhs,
00416 reinterpret_cast<doublecomplex*>(x.GetData()),
00417 x.GetM(), SLU_DN, SLU_Z, SLU_GE);
00418
00419 zgstrs (trans, &L, &U, perm_r.GetData(),
00420 perm_c.GetData(), &B, &stat, &info);
00421 }
00422
00423
00425 template<class Allocator2>
00426 void MatrixSuperLU<complex<double> >::
00427 Solve(const SeldonTranspose& TransA,
00428 Vector<complex<double>, VectFull, Allocator2>& x)
00429 {
00430 if (TransA.NoTrans())
00431 {
00432 Solve(x);
00433 return;
00434 }
00435
00436 trans_t trans = TRANS;
00437 int nb_rhs = 1, info;
00438 zCreate_Dense_Matrix(&B, x.GetM(), nb_rhs,
00439 reinterpret_cast<doublecomplex*>(x.GetData()),
00440 x.GetM(), SLU_DN, SLU_Z, SLU_GE);
00441
00442 zgstrs (trans, &L, &U, perm_r.GetData(),
00443 perm_c.GetData(), &B, &stat, &info);
00444 }
00445
00446
00447 template<class T, class Prop, class Storage, class Allocator>
00448 void GetLU(Matrix<T, Prop, Storage, Allocator>& A, MatrixSuperLU<T>& mat_lu,
00449 bool keep_matrix = false)
00450 {
00451 mat_lu.FactorizeMatrix(A, keep_matrix);
00452 }
00453
00454
00455 template<class T, class Allocator>
00456 void SolveLU(MatrixSuperLU<T>& mat_lu, Vector<T, VectFull, Allocator>& x)
00457 {
00458 mat_lu.Solve(x);
00459 }
00460
00461
00462 template<class T, class Allocator>
00463 void SolveLU(const SeldonTranspose& TransA,
00464 MatrixSuperLU<T>& mat_lu, Vector<T, VectFull, Allocator>& x)
00465 {
00466 mat_lu.Solve(TransA, x);
00467 }
00468
00469 }
00470
00471 #define SELDON_FILE_SUPERLU_CXX
00472 #endif