Warning: this documentation for the development version is under construction.
00001 // Copyright (C) 2011 Marc Duruflé 00002 // Copyright (C) 2010-2011 Vivien Mallet 00003 // 00004 // This file is part of the linear-algebra library Seldon, 00005 // http://seldon.sourceforge.net/. 00006 // 00007 // Seldon is free software; you can redistribute it and/or modify it under the 00008 // terms of the GNU Lesser General Public License as published by the Free 00009 // Software Foundation; either version 2.1 of the License, or (at your option) 00010 // any later version. 00011 // 00012 // Seldon is distributed in the hope that it will be useful, but WITHOUT ANY 00013 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS 00014 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for 00015 // more details. 00016 // 00017 // You should have received a copy of the GNU Lesser General Public License 00018 // along with Seldon. If not, see http://www.gnu.org/licenses/. 00019 00020 00021 #ifndef SELDON_FILE_COMPUTATION_SPARSESOLVER_CXX 00022 #define SELDON_FILE_COMPUTATION_SPARSESOLVER_CXX 00023 00024 00025 #include "Ordering.cxx" 00026 #include "SparseSolver.hxx" 00027 00028 00029 #ifdef SELDON_WITH_UMFPACK 00030 #include "camd.h" 00031 #include "colamd.h" 00032 #endif 00033 00034 namespace Seldon 00035 { 00036 00037 00038 /************************* 00039 * Default Seldon solver * 00040 *************************/ 00041 00042 00043 template<class T, class Allocator> 00044 SparseSeldonSolver<T, Allocator>::SparseSeldonSolver() 00045 { 00046 print_level = -1; 00047 symmetric_matrix = false; 00048 permtol = 0.1; 00049 } 00050 00051 00052 template<class T, class Allocator> 00053 void SparseSeldonSolver<T, Allocator>::Clear() 00054 { 00055 mat_sym.Clear(); 00056 mat_unsym.Clear(); 00057 } 00058 00059 00060 template<class T, class Allocator> 00061 void SparseSeldonSolver<T, Allocator>::HideMessages() 00062 { 00063 print_level = -1; 00064 } 00065 00066 00067 template<class T, class Allocator> 00068 void SparseSeldonSolver<T, Allocator>::ShowMessages() 00069 { 00070 print_level = 1; 00071 } 00072 00073 00074 template<class T, class Allocator> 00075 double SparseSeldonSolver<T, Allocator>::GetPivotThreshold() const 00076 { 00077 return permtol; 00078 } 00079 00080 00081 template<class T, class Allocator> 00082 void SparseSeldonSolver<T, Allocator>::SetPivotThreshold(const double& a) 00083 { 00084 permtol = a; 00085 } 00086 00087 00088 template<class T, class Allocator> 00089 template<class T0, class Storage0, class Allocator0> 00090 void SparseSeldonSolver<T, Allocator>:: 00091 FactorizeMatrix(const IVect& perm, 00092 Matrix<T0, General, Storage0, Allocator0>& mat, 00093 bool keep_matrix) 00094 { 00095 IVect inv_permutation; 00096 00097 // We convert matrix to unsymmetric format. 00098 Copy(mat, mat_unsym); 00099 00100 // Old matrix is erased if needed. 00101 if (!keep_matrix) 00102 mat.Clear(); 00103 00104 // We keep permutation array in memory, and check it. 00105 int n = mat_unsym.GetM(); 00106 if (perm.GetM() != n) 00107 throw WrongArgument("FactorizeMatrix(IVect&, Matrix&, bool)", 00108 "Numbering array is of size " 00109 + to_str(perm.GetM()) 00110 + " while the matrix is of size " 00111 + to_str(mat.GetM()) + " x " 00112 + to_str(mat.GetN()) + "."); 00113 00114 permutation_row.Reallocate(n); 00115 permutation_col.Reallocate(n); 00116 inv_permutation.Reallocate(n); 00117 inv_permutation.Fill(-1); 00118 for (int i = 0; i < n; i++) 00119 { 00120 permutation_row(i) = i; 00121 permutation_col(i) = i; 00122 inv_permutation(perm(i)) = i; 00123 } 00124 00125 for (int i = 0; i < n; i++) 00126 if (inv_permutation(i) == -1) 00127 throw WrongArgument("FactorizeMatrix(IVect&, Matrix&, bool)", 00128 "The numbering array is invalid."); 00129 00130 IVect iperm = inv_permutation; 00131 00132 // Rows of matrix are permuted. 00133 ApplyInversePermutation(mat_unsym, perm, perm); 00134 00135 // Temporary vector used for solving. 00136 xtmp.Reallocate(n); 00137 00138 // Factorization is performed. 00139 // Columns are permuted during the factorization. 00140 symmetric_matrix = false; 00141 inv_permutation.Fill(); 00142 GetLU(mat_unsym, permutation_col, inv_permutation, permtol, print_level); 00143 00144 // Combining permutations. 00145 IVect itmp = permutation_col; 00146 for (int i = 0; i < n; i++) 00147 permutation_col(i) = iperm(itmp(i)); 00148 00149 permutation_row = perm; 00150 00151 } 00152 00153 00154 template<class T, class Allocator> 00155 template<class T0, class Storage0, class Allocator0> 00156 void SparseSeldonSolver<T, Allocator>:: 00157 FactorizeMatrix(const IVect& perm, 00158 Matrix<T0, Symmetric, Storage0, Allocator0>& mat, 00159 bool keep_matrix) 00160 { 00161 IVect inv_permutation; 00162 00163 // We convert matrix to symmetric format. 00164 Copy(mat, mat_sym); 00165 00166 // Old matrix is erased if needed. 00167 if (!keep_matrix) 00168 mat.Clear(); 00169 00170 // We keep permutation array in memory, and check it. 00171 int n = mat_sym.GetM(); 00172 if (perm.GetM() != n) 00173 throw WrongArgument("FactorizeMatrix(IVect&, Matrix&, bool)", 00174 "Numbering array is of size " 00175 + to_str(perm.GetM()) 00176 + " while the matrix is of size " 00177 + to_str(mat.GetM()) + " x " 00178 + to_str(mat.GetN()) + "."); 00179 00180 permutation_row.Reallocate(n); 00181 inv_permutation.Reallocate(n); 00182 inv_permutation.Fill(-1); 00183 for (int i = 0; i < n; i++) 00184 { 00185 permutation_row(i) = perm(i); 00186 inv_permutation(perm(i)) = i; 00187 } 00188 00189 for (int i = 0; i < n; i++) 00190 if (inv_permutation(i) == -1) 00191 throw WrongArgument("FactorizeMatrix(IVect&, Matrix&, bool)", 00192 "The numbering array is invalid."); 00193 00194 // Matrix is permuted. 00195 ApplyInversePermutation(mat_sym, perm, perm); 00196 00197 // Temporary vector used for solving. 00198 xtmp.Reallocate(n); 00199 00200 // Factorization is performed. 00201 symmetric_matrix = true; 00202 GetLU(mat_sym, print_level); 00203 } 00204 00205 00206 template<class T, class Allocator> template<class Vector1> 00207 void SparseSeldonSolver<T, Allocator>::Solve(Vector1& z) 00208 { 00209 if (symmetric_matrix) 00210 { 00211 for (int i = 0; i < z.GetM(); i++) 00212 xtmp(permutation_row(i)) = z(i); 00213 00214 SolveLU(mat_sym, xtmp); 00215 00216 for (int i = 0; i < z.GetM(); i++) 00217 z(i) = xtmp(permutation_row(i)); 00218 } 00219 else 00220 { 00221 for (int i = 0; i < z.GetM(); i++) 00222 xtmp(permutation_row(i)) = z(i); 00223 00224 SolveLU(mat_unsym, xtmp); 00225 00226 for (int i = 0; i < z.GetM(); i++) 00227 z(permutation_col(i)) = xtmp(i); 00228 } 00229 } 00230 00231 00232 template<class T, class Allocator> 00233 template<class TransStatus, class Vector1> 00234 void SparseSeldonSolver<T, Allocator> 00235 ::Solve(const TransStatus& TransA, Vector1& z) 00236 { 00237 if (symmetric_matrix) 00238 Solve(z); 00239 else 00240 if (TransA.Trans()) 00241 { 00242 for (int i = 0; i < z.GetM(); i++) 00243 xtmp(i) = z(permutation_col(i)); 00244 00245 SolveLU(SeldonTrans, mat_unsym, xtmp); 00246 00247 for (int i = 0; i < z.GetM(); i++) 00248 z(i) = xtmp(permutation_row(i)); 00249 } 00250 else 00251 Solve(z); 00252 } 00253 00254 00255 /************************************************ 00256 * GetLU and SolveLU used by SeldonSparseSolver * 00257 ************************************************/ 00258 00259 00261 template<class T, class Allocator> 00262 void GetLU(Matrix<T, General, ArrayRowSparse, Allocator>& A, 00263 IVect& iperm, IVect& rperm, 00264 double permtol, int print_level) 00265 { 00266 int n = A.GetN(); 00267 00268 T fact, s, t; 00269 double tnorm, zero = 0.0; 00270 int length_lower, length_upper, jpos, jrow, i_row, j_col; 00271 int i, j, k, length, size_row, index_lu; 00272 00273 T czero, cone; 00274 SetComplexZero(czero); 00275 SetComplexOne(cone); 00276 Vector<T, VectFull, Allocator> Row_Val(n); 00277 IVect Index(n), Row_Ind(n), Index_Diag(n); 00278 Row_Val.Fill(czero); 00279 Row_Ind.Fill(-1); 00280 Index_Diag.Fill(-1); 00281 00282 Index.Fill(-1); 00283 // main loop 00284 int new_percent = 0, old_percent = 0; 00285 for (i_row = 0; i_row < n; i_row++) 00286 { 00287 // Progress bar if print level is high enough. 00288 if (print_level > 0) 00289 { 00290 new_percent = int(double(i_row) / double(n-1) * 78.); 00291 for (int percent = old_percent; percent < new_percent; percent++) 00292 { 00293 cout << "#"; 00294 cout.flush(); 00295 } 00296 old_percent = new_percent; 00297 } 00298 00299 size_row = A.GetRowSize(i_row); 00300 tnorm = zero; 00301 00302 // tnorm is the sum of absolute value of coefficients of row i_row. 00303 for (k = 0 ; k < size_row; k++) 00304 if (A.Value(i_row, k) != czero) 00305 tnorm += abs(A.Value(i_row, k)); 00306 00307 if (tnorm == zero) 00308 throw WrongArgument("GetLU(Matrix<ArrayRowSparse>&, IVect&, " 00309 "IVect&, double, int)", 00310 "The matrix is structurally singular. " 00311 "The norm of row " + to_str(i_row) 00312 + " is equal to 0."); 00313 00314 // Unpack L-part and U-part of row of A. 00315 length_upper = 1; 00316 length_lower = 0; 00317 Row_Ind(i_row) = i_row; 00318 Row_Val(i_row) = czero; 00319 Index(i_row) = i_row; 00320 00321 for (j = 0; j < size_row; j++) 00322 { 00323 k = rperm(A.Index(i_row, j)); 00324 00325 t = A.Value(i_row,j); 00326 if (k < i_row) 00327 { 00328 Row_Ind(length_lower) = k; 00329 Row_Val(length_lower) = t; 00330 Index(k) = length_lower; 00331 ++length_lower; 00332 } 00333 else if (k == i_row) 00334 { 00335 Row_Val(i_row) = t; 00336 } 00337 else 00338 { 00339 jpos = i_row + length_upper; 00340 Row_Ind(jpos) = k; 00341 Row_Val(jpos) = t; 00342 Index(k) = jpos; 00343 length_upper++; 00344 } 00345 } 00346 00347 j_col = 0; 00348 length = 0; 00349 00350 // Eliminates previous rows. 00351 while (j_col < length_lower) 00352 { 00353 // In order to do the elimination in the correct order, we must 00354 // select the smallest column index. 00355 jrow = Row_Ind(j_col); 00356 k = j_col; 00357 00358 // Determine smallest column index. 00359 for (j = j_col + 1; j < length_lower; j++) 00360 if (Row_Ind(j) < jrow) 00361 { 00362 jrow = Row_Ind(j); 00363 k = j; 00364 } 00365 00366 if (k != j_col) 00367 { 00368 // Exchanging column coefficients. 00369 j = Row_Ind(j_col); 00370 Row_Ind(j_col) = Row_Ind(k); 00371 Row_Ind(k) = j; 00372 00373 Index(jrow) = j_col; 00374 Index(j) = k; 00375 00376 s = Row_Val(j_col); 00377 Row_Val(j_col) = Row_Val(k); 00378 Row_Val(k) = s; 00379 } 00380 00381 // Zero out element in row. 00382 Index(jrow) = -1; 00383 00384 // Gets the multiplier for row to be eliminated (jrow). 00385 // first_index_upper points now on the diagonal coefficient. 00386 fact = Row_Val(j_col) * A.Value(jrow, Index_Diag(jrow)); 00387 00388 // Combines current row and row jrow. 00389 for (k = (Index_Diag(jrow)+1); k < A.GetRowSize(jrow); k++) 00390 { 00391 s = fact * A.Value(jrow,k); 00392 j = rperm(A.Index(jrow,k)); 00393 00394 jpos = Index(j); 00395 00396 if (j >= i_row) 00397 // Dealing with upper part. 00398 if (jpos == -1) 00399 { 00400 // This is a fill-in element. 00401 i = i_row + length_upper; 00402 Row_Ind(i) = j; 00403 Index(j) = i; 00404 Row_Val(i) = -s; 00405 ++length_upper; 00406 } 00407 else 00408 // This is not a fill-in element. 00409 Row_Val(jpos) -= s; 00410 else 00411 // Dealing with lower part. 00412 if (jpos == -1) 00413 { 00414 // this is a fill-in element 00415 Row_Ind(length_lower) = j; 00416 Index(j) = length_lower; 00417 Row_Val(length_lower) = -s; 00418 ++length_lower; 00419 } 00420 else 00421 // This is not a fill-in element. 00422 Row_Val(jpos) -= s; 00423 } 00424 00425 // Stores this pivot element from left to right -- no danger 00426 // of overlap with the working elements in L (pivots). 00427 Row_Val(length) = fact; 00428 Row_Ind(length) = jrow; 00429 ++length; 00430 j_col++; 00431 } 00432 00433 // Resets double-pointer to zero (U-part). 00434 for (k = 0; k < length_upper; k++) 00435 Index(Row_Ind(i_row + k)) = -1; 00436 00437 size_row = length; 00438 A.ReallocateRow(i_row,size_row); 00439 00440 // store L-part 00441 index_lu = 0; 00442 for (k = 0 ; k < length ; k++) 00443 { 00444 A.Value(i_row,index_lu) = Row_Val(k); 00445 A.Index(i_row,index_lu) = iperm(Row_Ind(k)); 00446 ++index_lu; 00447 } 00448 00449 // Saves pointer to beginning of row i_row of U. 00450 Index_Diag(i_row) = index_lu; 00451 00452 // Updates. U-matrix -- first apply dropping strategy. 00453 length = 0; 00454 for (k = 1; k <= (length_upper-1); k++) 00455 { 00456 ++length; 00457 Row_Val(i_row+length) = Row_Val(i_row+k); 00458 Row_Ind(i_row+length) = Row_Ind(i_row+k); 00459 } 00460 00461 length++; 00462 00463 // Determines next pivot. 00464 int imax = i_row; 00465 double xmax = abs(Row_Val(imax)); 00466 double xmax0 = xmax; 00467 for (k = i_row + 1; k <= i_row + length - 1; k++) 00468 { 00469 tnorm = abs(Row_Val(k)); 00470 if (tnorm > xmax && tnorm * permtol > xmax0) 00471 { 00472 imax = k; 00473 xmax = tnorm; 00474 } 00475 } 00476 00477 // Exchanges Row_Val. 00478 s = Row_Val(i_row); 00479 Row_Val(i_row) = Row_Val(imax); 00480 Row_Val(imax) = s; 00481 00482 // Updates iperm and reverses iperm. 00483 j = Row_Ind(imax); 00484 i = iperm(i_row); 00485 iperm(i_row) = iperm(j); 00486 iperm(j) = i; 00487 00488 // Reverses iperm. 00489 rperm(iperm(i_row)) = i_row; 00490 rperm(iperm(j)) = j; 00491 00492 // Copies U-part in original coordinates. 00493 int index_diag = index_lu; 00494 A.ResizeRow(i_row, size_row+length); 00495 00496 for (k = i_row ; k <= i_row + length - 1; k++) 00497 { 00498 A.Index(i_row,index_lu) = iperm(Row_Ind(k)); 00499 A.Value(i_row,index_lu) = Row_Val(k); 00500 ++index_lu; 00501 } 00502 00503 // Stores inverse of diagonal element of u. 00504 A.Value(i_row, index_diag) = cone / Row_Val(i_row); 00505 00506 } // end main loop 00507 00508 if (print_level > 0) 00509 cout << endl; 00510 00511 if (print_level > 0) 00512 cout << "The matrix takes " << 00513 int((A.GetDataSize() * (sizeof(T) + 4)) / (1024. * 1024.)) 00514 << " MB" << endl; 00515 00516 for (i = 0; i < n; i++ ) 00517 for (j = 0; j < A.GetRowSize(i); j++) 00518 A.Index(i,j) = rperm(A.Index(i,j)); 00519 } 00520 00521 00522 template<class cplx, 00523 class Allocator, class Storage2, class Allocator2> 00524 void SolveLU(const Matrix<cplx, General, ArrayRowSparse, Allocator>& A, 00525 Vector<cplx, Storage2, Allocator2>& x) 00526 { 00527 SolveLU(SeldonNoTrans, A, x); 00528 } 00529 00530 00532 00535 template<class cplx, class TransStatus, 00536 class Allocator, class Storage2, class Allocator2> 00537 void SolveLU(const TransStatus& transA, 00538 const Matrix<cplx, General, ArrayRowSparse, Allocator>& A, 00539 Vector<cplx, Storage2, Allocator2>& x) 00540 { 00541 int i, k, n, k_; 00542 cplx inv_diag; 00543 n = A.GetM(); 00544 00545 if (transA.Trans()) 00546 { 00547 // Forward solve (with U^T). 00548 for (i = 0 ; i < n ; i++) 00549 { 00550 k_ = 0; k = A.Index(i, k_); 00551 while ( k < i) 00552 { 00553 k_++; 00554 k = A.Index(i, k_); 00555 } 00556 00557 x(i) *= A.Value(i, k_); 00558 for (k = k_ + 1; k < A.GetRowSize(i); k++) 00559 x(A.Index(i, k)) -= A.Value(i, k) * x(i); 00560 } 00561 00562 // Backward solve (with L^T). 00563 for (i = n-1; i >= 0; i--) 00564 { 00565 k_ = 0; k = A.Index(i, k_); 00566 while (k < i) 00567 { 00568 x(k) -= A.Value(i, k_) * x(i); 00569 k_++; 00570 k = A.Index(i, k_); 00571 } 00572 } 00573 } 00574 else 00575 { 00576 // Forward solve. 00577 for (i = 0; i < n; i++) 00578 { 00579 k_ = 0; 00580 k = A.Index(i, k_); 00581 while (k < i) 00582 { 00583 x(i) -= A.Value(i, k_) * x(k); 00584 k_++; 00585 k = A.Index(i, k_); 00586 } 00587 } 00588 00589 // Backward solve. 00590 for (i = n-1; i >= 0; i--) 00591 { 00592 k_ = 0; 00593 k = A.Index(i, k_); 00594 while (k < i) 00595 { 00596 k_++; 00597 k = A.Index(i, k_); 00598 } 00599 00600 inv_diag = A.Value(i, k_); 00601 for (k = k_ + 1; k < A.GetRowSize(i); k++) 00602 x(i) -= A.Value(i, k) * x(A.Index(i, k)); 00603 00604 x(i) *= inv_diag; 00605 } 00606 } 00607 } 00608 00609 00611 template<class T, class Allocator> 00612 void GetLU(Matrix<T, Symmetric, ArrayRowSymSparse, Allocator>& A, 00613 int print_level) 00614 { 00615 int size_row; 00616 int n = A.GetN(); 00617 double zero = 0.0; 00618 00619 T fact, s, t; 00620 double tnorm; 00621 int length_lower, length_upper, jpos, jrow; 00622 int i_row, j_col, index_lu, length; 00623 int i, j, k; 00624 00625 Vector<T, VectFull, Allocator> Row_Val(n); 00626 IVect Index(n), Row_Ind(n); 00627 Row_Val.Zero(); 00628 Row_Ind.Fill(-1); 00629 00630 Index.Fill(-1); 00631 00632 // We convert A into an unsymmetric matrix. 00633 Matrix<T, General, ArrayRowSparse, Allocator> B; 00634 Seldon::Copy(A, B); 00635 00636 A.Clear(); 00637 A.Reallocate(n, n); 00638 00639 // Main loop. 00640 int new_percent = 0, old_percent = 0; 00641 for (i_row = 0; i_row < n; i_row++) 00642 { 00643 // Progress bar if print level is high enough. 00644 if (print_level > 0) 00645 { 00646 new_percent = int(double(i_row) / double(n-1) * 78.); 00647 for (int percent = old_percent; percent < new_percent; percent++) 00648 { 00649 cout << "#"; 00650 cout.flush(); 00651 } 00652 old_percent = new_percent; 00653 } 00654 00655 // 1-norm of the row of initial matrix. 00656 size_row = B.GetRowSize(i_row); 00657 tnorm = zero; 00658 for (k = 0 ; k < size_row; k++) 00659 tnorm += abs(B.Value(i_row, k)); 00660 00661 if (tnorm == zero) 00662 throw WrongArgument("GetLU(Matrix<ArrayRowSymSparse>&, int)", 00663 "The matrix is structurally singular. " 00664 "The norm of row " + to_str(i_row) 00665 + " is equal to 0."); 00666 00667 // Separating lower part from upper part for this row. 00668 length_upper = 1; 00669 length_lower = 0; 00670 Row_Ind(i_row) = i_row; 00671 Row_Val(i_row) = 0.0; 00672 Index(i_row) = i_row; 00673 00674 for (j = 0; j < size_row; j++) 00675 { 00676 k = B.Index(i_row,j); 00677 t = B.Value(i_row,j); 00678 if (k < i_row) 00679 { 00680 Row_Ind(length_lower) = k; 00681 Row_Val(length_lower) = t; 00682 Index(k) = length_lower; 00683 ++length_lower; 00684 } 00685 else if (k == i_row) 00686 { 00687 Row_Val(i_row) = t; 00688 } 00689 else 00690 { 00691 jpos = i_row + length_upper; 00692 Row_Ind(jpos) = k; 00693 Row_Val(jpos) = t; 00694 Index(k) = jpos; 00695 length_upper++; 00696 } 00697 } 00698 00699 // This row of B is cleared. 00700 B.ClearRow(i_row); 00701 00702 j_col = 0; 00703 length = 0; 00704 00705 // We eliminate previous rows. 00706 while (j_col <length_lower) 00707 { 00708 // In order to do the elimination in the correct order, we must 00709 // select the smallest column index. 00710 jrow = Row_Ind(j_col); 00711 k = j_col; 00712 00713 // We determine smallest column index. 00714 for (j = (j_col+1) ; j < length_lower; j++) 00715 { 00716 if (Row_Ind(j) < jrow) 00717 { 00718 jrow = Row_Ind(j); 00719 k = j; 00720 } 00721 } 00722 00723 // If needed, we exchange positions of this element in 00724 // Row_Ind/Row_Val so that it appears first. 00725 if (k != j_col) 00726 { 00727 00728 j = Row_Ind(j_col); 00729 Row_Ind(j_col) = Row_Ind(k); 00730 Row_Ind(k) = j; 00731 00732 Index(jrow) = j_col; 00733 Index(j) = k; 00734 00735 s = Row_Val(j_col); 00736 Row_Val(j_col) = Row_Val(k); 00737 Row_Val(k) = s; 00738 } 00739 00740 // Zero out element in row by setting Index to -1. 00741 Index(jrow) = -1; 00742 00743 // Gets the multiplier for row to be eliminated. 00744 fact = Row_Val(j_col) * A.Value(jrow, 0); 00745 00746 // Combines current row and row jrow. 00747 for (k = 1; k < A.GetRowSize(jrow); k++) 00748 { 00749 s = fact * A.Value(jrow, k); 00750 j = A.Index(jrow, k); 00751 00752 jpos = Index(j); 00753 if (j >= i_row) 00754 { 00755 00756 // Dealing with upper part. 00757 if (jpos == -1) 00758 { 00759 00760 // This is a fill-in element. 00761 i = i_row + length_upper; 00762 Row_Ind(i) = j; 00763 Index(j) = i; 00764 Row_Val(i) = -s; 00765 ++length_upper; 00766 } 00767 else 00768 { 00769 // This is not a fill-in element. 00770 Row_Val(jpos) -= s; 00771 } 00772 } 00773 else 00774 { 00775 // Dealing with lower part. 00776 if (jpos == -1) 00777 { 00778 // This is a fill-in element. 00779 Row_Ind(length_lower) = j; 00780 Index(j) = length_lower; 00781 Row_Val(length_lower) = -s; 00782 ++length_lower; 00783 } 00784 else 00785 { 00786 // This is not a fill-in element. 00787 Row_Val(jpos) -= s; 00788 } 00789 } 00790 } 00791 00792 // We store this pivot element (from left to right -- no 00793 // danger of overlap with the working elements in L (pivots). 00794 Row_Val(length) = fact; 00795 Row_Ind(length) = jrow; 00796 ++length; 00797 j_col++; 00798 } 00799 00800 // Resets double-pointer to zero (U-part). 00801 for (k = 0; k < length_upper; k++) 00802 Index(Row_Ind(i_row + k)) = -1; 00803 00804 // Updating U-matrix 00805 length = 0; 00806 for (k = 1; k <= length_upper - 1; k++) 00807 { 00808 ++length; 00809 Row_Val(i_row + length) = Row_Val(i_row + k); 00810 Row_Ind(i_row + length) = Row_Ind(i_row + k); 00811 } 00812 00813 length++; 00814 00815 // Copies U-part in matrix A. 00816 A.ReallocateRow(i_row, length); 00817 index_lu = 1; 00818 for (k = i_row + 1 ; k <= i_row + length - 1 ; k++) 00819 { 00820 A.Index(i_row, index_lu) = Row_Ind(k); 00821 A.Value(i_row, index_lu) = Row_Val(k); 00822 ++index_lu; 00823 } 00824 00825 // Stores the inverse of the diagonal element of u. 00826 A.Value(i_row,0) = 1.0 / Row_Val(i_row); 00827 00828 } // end main loop. 00829 00830 if (print_level > 0) 00831 cout << endl; 00832 00833 // for each row of A, we divide by diagonal value 00834 for (int i = 0; i < n; i++) 00835 for (int j = 1; j < A.GetRowSize(i); j++) 00836 A.Value(i, j) *= A.Value(i, 0); 00837 00838 if (print_level > 0) 00839 cout << "The matrix takes " << 00840 int((A.GetDataSize() * (sizeof(T) + 4)) / (1024. * 1024.)) 00841 << " MB" << endl; 00842 00843 } 00844 00845 00847 00850 template<class real, class cplx, class Allocator, 00851 class Storage2, class Allocator2> 00852 void SolveLU(const Matrix<real, Symmetric, ArrayRowSymSparse, Allocator>& A, 00853 Vector<cplx, Storage2, Allocator2>& x) 00854 { 00855 int n = A.GetM(), j_row; 00856 cplx tmp; 00857 00858 // We solve first L y = b. 00859 for (int i_col = 0; i_col < n ; i_col++) 00860 for (int k = 1; k < A.GetRowSize(i_col) ; k++) 00861 { 00862 j_row = A.Index(i_col, k); 00863 x(j_row) -= A.Value(i_col, k)*x(i_col); 00864 } 00865 00866 // Inverting by diagonal D. 00867 for (int i_col = 0; i_col < n ; i_col++) 00868 x(i_col) *= A.Value(i_col, 0); 00869 00870 // Then we solve L^t x = y. 00871 for (int i_col = n-1; i_col >=0; i_col--) 00872 { 00873 tmp = x(i_col); 00874 for (int k = 1; k < A.GetRowSize(i_col); k++) 00875 { 00876 j_row = A.Index(i_col,k); 00877 tmp -= A.Value(i_col,k) * x(j_row); 00878 } 00879 x(i_col) = tmp; 00880 } 00881 } 00882 00883 00884 /******************************************** 00885 * GetLU and SolveLU for SeldonSparseSolver * 00886 ********************************************/ 00887 00888 00889 template<class T, class Storage, class Allocator, class Alloc2> 00890 void GetLU(Matrix<T, Symmetric, Storage, Allocator>& A, 00891 SparseSeldonSolver<T, Alloc2>& mat_lu, 00892 bool keep_matrix = false) 00893 { 00894 // identity ordering 00895 IVect perm(A.GetM()); 00896 perm.Fill(); 00897 00898 // factorization 00899 mat_lu.FactorizeMatrix(perm, A, keep_matrix); 00900 } 00901 00902 00903 template<class T, class Storage, class Allocator, class Alloc2> 00904 void GetLU(Matrix<T, General, Storage, Allocator>& A, 00905 SparseSeldonSolver<T, Alloc2>& mat_lu, 00906 bool keep_matrix = false) 00907 { 00908 // identity ordering 00909 IVect perm(A.GetM()); 00910 perm.Fill(); 00911 00912 // factorization 00913 mat_lu.FactorizeMatrix(perm, A, keep_matrix); 00914 } 00915 00916 00917 template<class T, class Alloc2, class Allocator> 00918 void SolveLU(SparseSeldonSolver<T, Alloc2>& mat_lu, 00919 Vector<T, VectFull, Allocator>& x) 00920 { 00921 mat_lu.Solve(x); 00922 } 00923 00924 00925 template<class T, class Alloc2, class Allocator, class Transpose_status> 00926 void SolveLU(const Transpose_status& TransA, 00927 SparseSeldonSolver<T, Alloc2>& mat_lu, 00928 Vector<T, VectFull, Allocator>& x) 00929 { 00930 mat_lu.Solve(TransA, x); 00931 } 00932 00933 00934 /************************* 00935 * Choice of best solver * 00936 *************************/ 00937 00938 00940 template<class T> 00941 SparseDirectSolver<T>::SparseDirectSolver() 00942 { 00943 n = 0; 00944 type_ordering = SparseMatrixOrdering::AUTO; 00945 00946 type_solver = SELDON_SOLVER; 00947 00948 // We try to use an other available solver. 00949 // The order of preference is Pastix, Mumps, UmfPack and SuperLU. 00950 #ifdef SELDON_WITH_SUPERLU 00951 type_solver = SUPERLU; 00952 #endif 00953 #ifdef SELDON_WITH_UMFPACK 00954 type_solver = UMFPACK; 00955 #endif 00956 #ifdef SELDON_WITH_MUMPS 00957 type_solver = MUMPS; 00958 #endif 00959 #ifdef SELDON_WITH_PASTIX 00960 type_solver = PASTIX; 00961 #endif 00962 00963 number_threads_per_node = 1; 00964 threshold_matrix = 0; 00965 enforce_unsym_ilut = false; 00966 } 00967 00968 00970 template<class T> 00971 void SparseDirectSolver<T>::HideMessages() 00972 { 00973 mat_seldon.HideMessages(); 00974 00975 #ifdef SELDON_WITH_MUMPS 00976 mat_mumps.HideMessages(); 00977 #endif 00978 00979 #ifdef SELDON_WITH_SUPERLU 00980 mat_superlu.HideMessages(); 00981 #endif 00982 00983 #ifdef SELDON_WITH_UMFPACK 00984 mat_umf.HideMessages(); 00985 #endif 00986 00987 #ifdef SELDON_WITH_PASTIX 00988 mat_pastix.HideMessages(); 00989 #endif 00990 00991 #ifdef SELDON_WITH_PRECONDITIONING 00992 mat_ilut.SetPrintLevel(0); 00993 #endif 00994 } 00995 00996 00998 template<class T> 00999 void SparseDirectSolver<T>::ShowMessages() 01000 { 01001 mat_seldon.ShowMessages(); 01002 01003 #ifdef SELDON_WITH_MUMPS 01004 mat_mumps.ShowMessages(); 01005 #endif 01006 01007 #ifdef SELDON_WITH_SUPERLU 01008 mat_superlu.ShowMessages(); 01009 #endif 01010 01011 #ifdef SELDON_WITH_UMFPACK 01012 mat_umf.ShowMessages(); 01013 #endif 01014 01015 #ifdef SELDON_WITH_PASTIX 01016 mat_pastix.ShowMessages(); 01017 #endif 01018 01019 #ifdef SELDON_WITH_PRECONDITIONING 01020 mat_ilut.SetPrintLevel(1); 01021 #endif 01022 } 01023 01024 01026 template<class T> 01027 void SparseDirectSolver<T>::ShowFullHistory() 01028 { 01029 mat_seldon.ShowMessages(); 01030 01031 #ifdef SELDON_WITH_MUMPS 01032 mat_mumps.ShowMessages(); 01033 #endif 01034 01035 #ifdef SELDON_WITH_SUPERLU 01036 mat_superlu.ShowMessages(); 01037 #endif 01038 01039 #ifdef SELDON_WITH_UMFPACK 01040 mat_umf.ShowMessages(); 01041 #endif 01042 01043 #ifdef SELDON_WITH_PASTIX 01044 mat_pastix.ShowFullHistory(); 01045 #endif 01046 01047 #ifdef SELDON_WITH_PRECONDITIONING 01048 mat_ilut.SetPrintLevel(1); 01049 #endif 01050 } 01051 01052 01054 template<class T> 01055 void SparseDirectSolver<T>::Clear() 01056 { 01057 if (n > 0) 01058 { 01059 n = 0; 01060 mat_seldon.Clear(); 01061 01062 #ifdef SELDON_WITH_MUMPS 01063 mat_mumps.Clear(); 01064 #endif 01065 01066 #ifdef SELDON_WITH_SUPERLU 01067 mat_superlu.Clear(); 01068 #endif 01069 01070 #ifdef SELDON_WITH_UMFPACK 01071 mat_umf.Clear(); 01072 #endif 01073 01074 #ifdef SELDON_WITH_PASTIX 01075 mat_pastix.Clear(); 01076 #endif 01077 01078 #ifdef SELDON_WITH_PRECONDITIONING 01079 mat_ilut.Clear(); 01080 #endif 01081 } 01082 } 01083 01084 01086 template<class T> 01087 int SparseDirectSolver<T>::GetM() const 01088 { 01089 return n; 01090 } 01091 01092 01094 template<class T> 01095 int SparseDirectSolver<T>::GetN() const 01096 { 01097 return n; 01098 } 01099 01100 01102 template<class T> 01103 int SparseDirectSolver<T>::GetTypeOrdering() const 01104 { 01105 return type_ordering; 01106 } 01107 01108 01110 template<class T> 01111 void SparseDirectSolver<T>::SetPermutation(const IVect& num) 01112 { 01113 type_ordering = SparseMatrixOrdering::USER; 01114 permut = num; 01115 } 01116 01117 01119 template<class T> 01120 void SparseDirectSolver<T>::SelectOrdering(int type) 01121 { 01122 type_ordering = type; 01123 } 01124 01125 01127 template<class T> 01128 void SparseDirectSolver<T>::SetNumberThreadPerNode(int p) 01129 { 01130 number_threads_per_node = p; 01131 } 01132 01133 01135 template<class T> 01136 void SparseDirectSolver<T>::SelectDirectSolver(int type) 01137 { 01138 type_solver = type; 01139 } 01140 01141 01143 template<class T> 01144 void SparseDirectSolver<T>::SetNonSymmetricIlut() 01145 { 01146 enforce_unsym_ilut = true; 01147 } 01148 01149 01151 template<class T> 01152 int SparseDirectSolver<T>::GetDirectSolver() 01153 { 01154 return type_solver; 01155 } 01156 01157 01159 template<class T> 01160 double SparseDirectSolver<T>::GetThresholdMatrix() const 01161 { 01162 return threshold_matrix; 01163 } 01164 01165 01167 template<class T> template<class MatrixSparse> 01168 void SparseDirectSolver<T>::ComputeOrdering(MatrixSparse& A) 01169 { 01170 bool user_ordering = false; 01171 // We set the ordering for each direct solver interfaced. 01172 switch (type_ordering) 01173 { 01174 case SparseMatrixOrdering::AUTO : 01175 { 01176 // We choose the default strategy proposed by the direct solver that 01177 // will be called. 01178 if (type_solver == MUMPS) 01179 { 01180 #ifdef SELDON_WITH_MUMPS 01181 mat_mumps.SelectOrdering(7); 01182 #endif 01183 } 01184 else if (type_solver == PASTIX) 01185 { 01186 #ifdef SELDON_WITH_PASTIX 01187 mat_pastix.SelectOrdering(API_ORDER_SCOTCH); 01188 #endif 01189 } 01190 else if (type_solver == UMFPACK) 01191 { 01192 #ifdef SELDON_WITH_UMFPACK 01193 mat_umf.SelectOrdering(UMFPACK_ORDERING_AMD); 01194 #endif 01195 } 01196 else if (type_solver == SUPERLU) 01197 { 01198 #ifdef SELDON_WITH_SUPERLU 01199 mat_superlu.SelectOrdering(COLAMD); 01200 #endif 01201 } 01202 else 01203 { 01204 01205 type_ordering = SparseMatrixOrdering::IDENTITY; 01206 01207 #ifdef SELDON_WITH_UMFPACK 01208 type_ordering = SparseMatrixOrdering::AMD; 01209 #endif 01210 01211 #ifdef SELDON_WITH_MUMPS 01212 type_ordering = SparseMatrixOrdering::METIS; 01213 #endif 01214 01215 #ifdef SELDON_WITH_PASTIX 01216 type_ordering = SparseMatrixOrdering::SCOTCH; 01217 #endif 01218 01219 user_ordering = true; 01220 } 01221 } 01222 break; 01223 case SparseMatrixOrdering::IDENTITY : 01224 case SparseMatrixOrdering::REVERSE_CUTHILL_MCKEE : 01225 case SparseMatrixOrdering::USER : 01226 { 01227 user_ordering = true; 01228 } 01229 break; 01230 case SparseMatrixOrdering::PORD : 01231 case SparseMatrixOrdering::QAMD : 01232 { 01233 // Mumps ordering. 01234 if (type_solver == MUMPS) 01235 { 01236 #ifdef SELDON_WITH_MUMPS 01237 if (type_ordering == SparseMatrixOrdering::PORD) 01238 mat_mumps.SelectOrdering(4); 01239 else 01240 mat_mumps.SelectOrdering(6); 01241 #endif 01242 } 01243 else 01244 { 01245 user_ordering = true; 01246 } 01247 } 01248 break; 01249 case SparseMatrixOrdering::SCOTCH : 01250 { 01251 // Available for Mumps and Pastix. 01252 if (type_solver == MUMPS) 01253 { 01254 #ifdef SELDON_WITH_MUMPS 01255 mat_mumps.SelectOrdering(3); 01256 #endif 01257 } 01258 else if (type_solver == PASTIX) 01259 { 01260 #ifdef SELDON_WITH_SCOTCH 01261 mat_pastix.SelectOrdering(API_ORDER_SCOTCH); 01262 #endif 01263 } 01264 else 01265 { 01266 user_ordering = true; 01267 } 01268 } 01269 break; 01270 case SparseMatrixOrdering::METIS : 01271 { 01272 // Available for Mumps and Pastix. 01273 if (type_solver == MUMPS) 01274 { 01275 #ifdef SELDON_WITH_MUMPS 01276 mat_mumps.SelectOrdering(5); 01277 #endif 01278 } 01279 else if (type_solver == PASTIX) 01280 { 01281 #ifdef SELDON_WITH_SCOTCH 01282 mat_pastix.SelectOrdering(API_ORDER_METIS); 01283 #endif 01284 } 01285 else if (type_solver == UMFPACK) 01286 { 01287 #ifdef SELDON_WITH_SCOTCH 01288 mat_umf.SelectOrdering(UMFPACK_ORDERING_METIS); 01289 #endif 01290 } 01291 01292 else 01293 { 01294 user_ordering = true; 01295 } 01296 } 01297 break; 01298 case SparseMatrixOrdering::AMD : 01299 case SparseMatrixOrdering::COLAMD : 01300 { 01301 // Default ordering for UmfPack. 01302 if (type_solver == UMFPACK) 01303 { 01304 #ifdef SELDON_WITH_UMFPACK 01305 mat_umf.SelectOrdering(UMFPACK_ORDERING_AMD); 01306 #endif 01307 } 01308 else if (type_solver == SUPERLU) 01309 { 01310 #ifdef SELDON_WITH_SUPERLU 01311 mat_superlu.SelectOrdering(COLAMD); 01312 #endif 01313 } 01314 else 01315 { 01316 user_ordering = true; 01317 } 01318 } 01319 break; 01320 } 01321 01322 if (user_ordering) 01323 { 01324 // Case where the ordering is not natively available in the direct 01325 // solver computing separetely the ordering. 01326 FindSparseOrdering(A, permut, type_ordering); 01327 01328 if (type_solver == MUMPS) 01329 { 01330 #ifdef SELDON_WITH_MUMPS 01331 mat_mumps.SetPermutation(permut); 01332 #endif 01333 } 01334 else if (type_solver == PASTIX) 01335 { 01336 #ifdef SELDON_WITH_PASTIX 01337 mat_pastix.SetPermutation(permut); 01338 #endif 01339 } 01340 else if (type_solver == UMFPACK) 01341 { 01342 #ifdef SELDON_WITH_UMFPACK 01343 mat_umf.SetPermutation(permut); 01344 #endif 01345 } 01346 else if (type_solver == SUPERLU) 01347 { 01348 #ifdef SELDON_WITH_SUPERLU 01349 mat_superlu.SetPermutation(permut); 01350 #endif 01351 } 01352 else 01353 { 01354 } 01355 } 01356 01357 } 01358 01359 01361 01365 template<class T> template<class MatrixSparse> 01366 void SparseDirectSolver<T>::Factorize(MatrixSparse& A, bool keep_matrix) 01367 { 01368 ComputeOrdering(A); 01369 01370 n = A.GetM(); 01371 if (type_solver == UMFPACK) 01372 { 01373 #ifdef SELDON_WITH_UMFPACK 01374 mat_umf.FactorizeMatrix(A, keep_matrix); 01375 #else 01376 throw Undefined("SparseDirectSolver::Factorize(MatrixSparse&, bool)", 01377 "Seldon was not compiled with UmfPack support."); 01378 #endif 01379 } 01380 else if (type_solver == SUPERLU) 01381 { 01382 #ifdef SELDON_WITH_SUPERLU 01383 mat_superlu.FactorizeMatrix(A, keep_matrix); 01384 #else 01385 throw Undefined("SparseDirectSolver::Factorize(MatrixSparse&, bool)", 01386 "Seldon was not compiled with SuperLU support."); 01387 #endif 01388 } 01389 else if (type_solver == MUMPS) 01390 { 01391 #ifdef SELDON_WITH_MUMPS 01392 mat_mumps.FactorizeMatrix(A, keep_matrix); 01393 #else 01394 throw Undefined("SparseDirectSolver::Factorize(MatrixSparse&, bool)", 01395 "Seldon was not compiled with Mumps support."); 01396 #endif 01397 } 01398 else if (type_solver == PASTIX) 01399 { 01400 #ifdef SELDON_WITH_PASTIX 01401 mat_pastix.SetNumberThreadPerNode(number_threads_per_node); 01402 mat_pastix.FactorizeMatrix(A, keep_matrix); 01403 #else 01404 throw Undefined("SparseDirectSolver::Factorize(MatrixSparse&, bool)", 01405 "Seldon was not compiled with Pastix support."); 01406 #endif 01407 } 01408 else if (type_solver == ILUT) 01409 { 01410 #ifdef SELDON_WITH_PRECONDITIONING 01411 // Setting some parameters. 01412 if (enforce_unsym_ilut || !IsSymmetricMatrix(A)) 01413 mat_ilut.SetUnsymmetricAlgorithm(); 01414 else 01415 mat_ilut.SetSymmetricAlgorithm(); 01416 01417 // Then performing factorization. 01418 mat_ilut.FactorizeMatrix(permut, A, keep_matrix); 01419 #else 01420 throw Undefined("SparseDirectSolver::Factorize(MatrixSparse&, bool)", 01421 "Seldon was not compiled with the preconditioners."); 01422 #endif 01423 } 01424 else 01425 { 01426 mat_seldon.FactorizeMatrix(permut, A); 01427 } 01428 01429 } 01430 01431 01433 template <class T> 01434 int SparseDirectSolver<T>::GetInfoFactorization(int& ierr) const 01435 { 01436 if (type_solver == UMFPACK) 01437 { 01438 #ifdef SELDON_WITH_UMFPACK 01439 ierr = mat_umf.GetInfoFactorization(); 01440 switch (ierr) 01441 { 01442 case UMFPACK_OK : 01443 return FACTO_OK; 01444 case UMFPACK_WARNING_singular_matrix : 01445 return NUMERICALLY_SINGULAR_MATRIX; 01446 case UMFPACK_ERROR_out_of_memory : 01447 return OUT_OF_MEMORY; 01448 case UMFPACK_ERROR_invalid_Numeric_object : 01449 case UMFPACK_ERROR_invalid_Symbolic_object : 01450 case UMFPACK_ERROR_argument_missing : 01451 case UMFPACK_ERROR_different_pattern : 01452 case UMFPACK_ERROR_invalid_system : 01453 return INVALID_ARGUMENT; 01454 case UMFPACK_ERROR_n_nonpositive : 01455 return INCORRECT_NUMBER_OF_ROWS; 01456 case UMFPACK_ERROR_invalid_matrix : 01457 return MATRIX_INDICES_INCORRECT; 01458 case UMFPACK_ERROR_invalid_permutation : 01459 return INVALID_PERMUTATION; 01460 case UMFPACK_ERROR_ordering_failed : 01461 return ORDERING_FAILED; 01462 default : 01463 return INTERNAL_ERROR; 01464 } 01465 #endif 01466 } 01467 else if (type_solver == SUPERLU) 01468 { 01469 #ifdef SELDON_WITH_SUPERLU 01470 ierr = mat_superlu.GetInfoFactorization(); 01471 if (ierr > 0) 01472 return INTERNAL_ERROR; 01473 #endif 01474 } 01475 else if (type_solver == MUMPS) 01476 { 01477 #ifdef SELDON_WITH_MUMPS 01478 ierr = mat_mumps.GetInfoFactorization(); 01479 switch (ierr) 01480 { 01481 case -2 : 01482 // nz out of range 01483 return MATRIX_INDICES_INCORRECT; 01484 case -3 : 01485 // invalid job number 01486 return INVALID_ARGUMENT; 01487 case -4 : 01488 // invalid permutation 01489 return INVALID_PERMUTATION; 01490 case -5 : 01491 // problem of real workspace allocation 01492 return OUT_OF_MEMORY; 01493 case -6 : 01494 // structurally singular matrix 01495 return STRUCTURALLY_SINGULAR_MATRIX; 01496 case -7 : 01497 // problem of integer workspace allocation 01498 return OUT_OF_MEMORY; 01499 case -10 : 01500 // numerically singular matrix 01501 return NUMERICALLY_SINGULAR_MATRIX; 01502 case -13 : 01503 // allocate failed 01504 return OUT_OF_MEMORY; 01505 case -16 : 01506 // N out of range 01507 return INCORRECT_NUMBER_OF_ROWS; 01508 case -22 : 01509 // invalid pointer 01510 return INVALID_ARGUMENT; 01511 case 1 : 01512 // index out of range 01513 return MATRIX_INDICES_INCORRECT; 01514 case 0 : 01515 return FACTO_OK; 01516 default : 01517 return INTERNAL_ERROR; 01518 } 01519 #endif 01520 } 01521 01522 return FACTO_OK; 01523 } 01524 01525 01527 01530 template<class T> template<class Vector1> 01531 void SparseDirectSolver<T>::Solve(Vector1& x_solution) 01532 { 01533 if (type_solver == UMFPACK) 01534 { 01535 #ifdef SELDON_WITH_UMFPACK 01536 Seldon::SolveLU(mat_umf, x_solution); 01537 #else 01538 throw Undefined("SparseDirectSolver::Solve(Vector&)", 01539 "Seldon was not compiled with UmfPack support."); 01540 #endif 01541 } 01542 else if (type_solver == SUPERLU) 01543 { 01544 #ifdef SELDON_WITH_SUPERLU 01545 Seldon::SolveLU(mat_superlu, x_solution); 01546 #else 01547 throw Undefined("SparseDirectSolver::Solve(Vector&)", 01548 "Seldon was not compiled with SuperLU support."); 01549 #endif 01550 } 01551 else if (type_solver == MUMPS) 01552 { 01553 #ifdef SELDON_WITH_MUMPS 01554 Seldon::SolveLU(mat_mumps, x_solution); 01555 #else 01556 throw Undefined("SparseDirectSolver::Solve(Vector&)", 01557 "Seldon was not compiled with Mumps support."); 01558 #endif 01559 } 01560 else if (type_solver == PASTIX) 01561 { 01562 #ifdef SELDON_WITH_PASTIX 01563 Seldon::SolveLU(mat_pastix, x_solution); 01564 #else 01565 throw Undefined("SparseDirectSolver::Solve(Vector&)", 01566 "Seldon was not compiled with Pastix support."); 01567 #endif 01568 } 01569 else if (type_solver == ILUT) 01570 { 01571 #ifdef SELDON_WITH_PRECONDITIONING 01572 mat_ilut.Solve(x_solution); 01573 #else 01574 throw Undefined("SparseDirectSolver::Solve(Vector&)", 01575 "Seldon was not compiled with the preconditioners."); 01576 #endif 01577 } 01578 else 01579 { 01580 Seldon::SolveLU(mat_seldon, x_solution); 01581 } 01582 } 01583 01584 01586 template<class T> template<class TransStatus, class Vector1> 01587 void SparseDirectSolver<T> 01588 ::Solve(const TransStatus& TransA, Vector1& x_solution) 01589 { 01590 if (type_solver == UMFPACK) 01591 { 01592 #ifdef SELDON_WITH_UMFPACK 01593 Seldon::SolveLU(TransA, mat_umf, x_solution); 01594 #else 01595 throw Undefined("SparseDirectSolver::Solve(TransStatus, Vector&)", 01596 "Seldon was not compiled with UmpfPack support."); 01597 #endif 01598 } 01599 else if (type_solver == SUPERLU) 01600 { 01601 #ifdef SELDON_WITH_SUPERLU 01602 Seldon::SolveLU(TransA, mat_superlu, x_solution); 01603 #else 01604 throw Undefined("SparseDirectSolver::Solve(TransStatus, Vector&)", 01605 "Seldon was not compiled with SuperLU support."); 01606 #endif 01607 } 01608 else if (type_solver == MUMPS) 01609 { 01610 #ifdef SELDON_WITH_MUMPS 01611 Seldon::SolveLU(TransA, mat_mumps, x_solution); 01612 #else 01613 throw Undefined("SparseDirectSolver::Solve(TransStatus, Vector&)", 01614 "Seldon was not compiled with Mumps support."); 01615 #endif 01616 } 01617 else if (type_solver == PASTIX) 01618 { 01619 #ifdef SELDON_WITH_PASTIX 01620 Seldon::SolveLU(TransA, mat_pastix, x_solution); 01621 #else 01622 throw Undefined("SparseDirectSolver::Solve(TransStatus, Vector&)", 01623 "Seldon was not compiled with Pastix support."); 01624 #endif 01625 } 01626 else if (type_solver == ILUT) 01627 { 01628 #ifdef SELDON_WITH_PRECONDITIONING 01629 mat_ilut.Solve(TransA, x_solution); 01630 #else 01631 throw Undefined("SparseDirectSolver::Solve(TransStatus, Vector&)", 01632 "Seldon was not compiled with the preconditioners."); 01633 #endif 01634 } 01635 else 01636 { 01637 Seldon::SolveLU(TransA, mat_seldon, x_solution); 01638 } 01639 } 01640 01641 01642 #ifdef SELDON_WITH_MPI 01643 01644 01657 template<class T> template<class Tint> 01658 void SparseDirectSolver<T>:: 01659 FactorizeDistributed(MPI::Comm& comm_facto, 01660 Vector<Tint>& Ptr, Vector<Tint>& Row, Vector<T>& Val, 01661 const IVect& glob_num, bool sym, bool keep_matrix) 01662 { 01663 n = Ptr.GetM()-1; 01664 if (type_solver == MUMPS) 01665 { 01666 #ifdef SELDON_WITH_MUMPS 01667 mat_mumps.FactorizeDistributedMatrix(comm_facto, Ptr, Row, Val, 01668 glob_num, sym, keep_matrix); 01669 #else 01670 throw Undefined("SparseDirectSolver::FactorizeDistributed(MPI::Comm&," 01671 " IVect&, IVect&, Vector<T>&, IVect&, bool, bool)", 01672 "Seldon was not compiled with Mumps support."); 01673 #endif 01674 } 01675 else if (type_solver == PASTIX) 01676 { 01677 #ifdef SELDON_WITH_PASTIX 01678 mat_pastix.FactorizeDistributedMatrix(comm_facto, Ptr, Row, Val, 01679 glob_num, sym, keep_matrix); 01680 #else 01681 throw Undefined("SparseDirectSolver::FactorizeDistributed(MPI::Comm&," 01682 " IVect&, IVect&, Vector<T>&, IVect&, bool, bool)", 01683 "Seldon was not compiled with Pastix support."); 01684 #endif 01685 } 01686 else 01687 { 01688 throw Undefined("SparseDirectSolver::FactorizeDistributed(MPI::Comm&," 01689 " IVect&, IVect&, Vector<T>&, IVect&, bool, bool)", 01690 "The method is defined for Mumps and Pastix only."); 01691 } 01692 } 01693 01694 01696 01699 template<class T> template<class Vector1> 01700 void SparseDirectSolver<T>:: 01701 SolveDistributed(MPI::Comm& comm_facto, Vector1& x_solution, 01702 const IVect& glob_number) 01703 { 01704 if (type_solver == MUMPS) 01705 { 01706 #ifdef SELDON_WITH_MUMPS 01707 mat_mumps.SolveDistributed(comm_facto, x_solution, glob_number); 01708 #else 01709 throw Undefined("SparseDirectSolver::SolveDistributed(MPI::Comm&," 01710 " Vector&, IVect&)", 01711 "Seldon was not compiled with Mumps support."); 01712 #endif 01713 } 01714 else if (type_solver == PASTIX) 01715 { 01716 #ifdef SELDON_WITH_PASTIX 01717 mat_pastix.SolveDistributed(comm_facto, x_solution, glob_number); 01718 #else 01719 throw Undefined("SparseDirectSolver::SolveDistributed(MPI::Comm&," 01720 " Vector&, IVect&)", 01721 "Seldon was not compiled with Pastix support."); 01722 #endif 01723 } 01724 else 01725 { 01726 throw Undefined("SparseDirectSolver::SolveDistributed(MPI::Comm&," 01727 " Vector&, IVect&)", 01728 "The method is defined for Mumps and Pastix only."); 01729 } 01730 01731 } 01732 01733 01739 template<class T> template<class TransStatus, class Vector1> 01740 void SparseDirectSolver<T>:: 01741 SolveDistributed(MPI::Comm& comm_facto, const TransStatus& TransA, 01742 Vector1& x_solution, const IVect& glob_number) 01743 { 01744 if (type_solver == MUMPS) 01745 { 01746 #ifdef SELDON_WITH_MUMPS 01747 mat_mumps.SolveDistributed(comm_facto, TransA, x_solution, 01748 glob_number); 01749 #else 01750 throw Undefined("SparseDirectSolver::SolveDistributed(TransStatus, " 01751 "MPI::Comm&, Vector&, IVect&)", 01752 "Seldon was not compiled with Mumps support."); 01753 #endif 01754 } 01755 else if (type_solver == PASTIX) 01756 { 01757 #ifdef SELDON_WITH_PASTIX 01758 mat_pastix.SolveDistributed(comm_facto, TransA, x_solution, glob_number); 01759 #else 01760 throw Undefined("SparseDirectSolver::SolveDistributed(TransStatus, " 01761 "MPI::Comm&, Vector&, IVect&)", 01762 "Seldon was not compiled with Pastix support."); 01763 #endif 01764 } 01765 else 01766 { 01767 throw Undefined("SparseDirectSolver::SolveDistributed(TransStatus, " 01768 "MPI::Comm&, Vector&, IVect&)", 01769 "The method is defined for Mumps and Pastix only."); 01770 } 01771 } 01772 #endif 01773 01774 01775 /************************* 01776 * Solve and SparseSolve * 01777 *************************/ 01778 01779 01781 01788 template <class T0, class Prop0, class Storage0, class Allocator0, 01789 class T1, class Storage1, class Allocator1> 01790 void SparseSolve(Matrix<T0, Prop0, Storage0, Allocator0>& M, 01791 Vector<T1, Storage1, Allocator1>& Y) 01792 { 01793 #ifdef SELDON_CHECK_DIMENSIONS 01794 CheckDim(M, Y, "SparseSolve(M, Y)"); 01795 #endif 01796 01797 SparseDirectSolver<T0> matrix_lu; 01798 01799 matrix_lu.Factorize(M); 01800 matrix_lu.Solve(Y); 01801 } 01802 01803 01805 template <class T, class Prop0, class Allocator0, class Allocator1> 01806 void Solve(Matrix<T, Prop0, ColSparse, Allocator0>& M, 01807 Vector<T, VectFull, Allocator1>& Y) 01808 { 01809 SparseSolve(M, Y); 01810 } 01811 01812 01814 template <class T, class Prop0, class Allocator0, class Allocator1> 01815 void Solve(Matrix<T, Prop0, RowSparse, Allocator0>& M, 01816 Vector<T, VectFull, Allocator1>& Y) 01817 { 01818 SparseSolve(M, Y); 01819 } 01820 01821 01822 } // namespace Seldon. 01823 01824 01825 #endif